Building neighor array
is neighID and stencil in element_type still used after this?
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | iStencil |
Index of your neighbor list?? |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc(tree%global%minLevel:) |
Level descriptor for each level of your mesh (starting from minimum level). |
||
type(treelmesh_type), | intent(in) | :: | tree |
Tree representation of your mesh. |
||
type(tem_stencilHeader_type), | intent(in) | :: | computeStencil |
The stencil you build the horizontal dependencies for. |
subroutine tem_build_horizontalDependencies( iStencil, levelDesc, tree, & & computeStencil ) ! --------------------------------------------------------------------------- !> Tree representation of your mesh. type(treelmesh_type), intent(in) :: tree !> Level descriptor for each level of your mesh (starting from minimum !! level). type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:) !> Index of your neighbor list?? integer, intent(in) :: iStencil !> The stencil you build the horizontal dependencies for. type(tem_stencilHeader_type), intent(in) :: computeStencil ! --------------------------------------------------------------------------- integer :: iLevel, iIndex integer :: nElemsWithNeigh ! --------------------------------------------------------------------------- write(logUnit(4),"(A,I0)") 'Build neighbor array for stencil: ', iStencil do iLevel = tree%global%minLevel,tree%global%maxLevel if( computeStencil%useAll .and. iStencil == 1 ) then ! JQ: is it possible the 1st stencil does not use all? ! all fluid and ghost and halo elements nElemsWithNeigh = levelDesc(iLevel)%nElems else if( computeStencil%useAll ) then ! @todo: the logic of this IF condition is not clear ! all fluid and ghost elements no halos nElemsWithNeigh = levelDesc(iLevel)%nElems & & - levelDesc(iLevel)%elem%nElems( eT_halo ) else ! only the fluids and ghosts??? that have the stencil iStencil nElemsWithNeigh = computeStencil%elemLvl( iLevel )%nVals end if write(logUnit(4),"(A,I2,A,I0)") ' level: ', iLevel, & & ', nElems to treat: ', nElemsWithNeigh ! allocate the nghElems array with the number of elements this stencil is ! used by allocate( levelDesc( iLevel )%neigh( iStencil )% & & nghElems( computeStencil%QQN, nElemsWithNeigh ) ) iIndex = 0 if( computeStencil%QQN > 0 ) then ! we init everything by default with 0 to find a bug as fast as possible levelDesc( iLevel )%neigh( iStencil )%nghElems(:,:) = 0 if (computeStencil%useAll) then ! write(logUnit(5),*) ' fluid elements ', iIndex call tem_build_listHorizontalDep( & & levelDesc = levelDesc( iLevel ), & & iStencil = iStencil, & & posInSortElem = levelDesc( iLevel )%totalPnt, & & nElems = nElemsWithNeigh, & & iIndex = iIndex ) else ! stencil not using all elements ! write(logUnit(5),*) ' not all elems ', iIndex call tem_build_treeHorizontalDep( & & levelDesc = levelDesc( iLevel ), & & computeStencil = computeStencil, & & iStencil = iStencil, & & list = computeStencil%elemLvl( iLevel )%val, & & nElems = nElemsWithNeigh, & & tree = tree ) end if end if end do ! iLevel end subroutine tem_build_horizontalDependencies