create all the neighbors of an element's parent
Create all elements required up to the actual existing fluid element these include the neighbors of the parents. In a level jump >1, these intermediate levels have to provide valid quantities over two of their computation updates to account for the recursive algorithm.
Here the fromCoarser interpolation should be handed in.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=long_k), | intent(in) | :: | targetID |
requested element position (child element) in LevelDesc elem list |
||
integer, | intent(in) | :: | level |
requested element level |
||
type(tem_stencilHeader_type), | intent(in) | :: | stencil |
current stencil definition |
||
type(treelmesh_type), | intent(in) | :: | tree |
tree information |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc(tree%global%minLevel:) |
the level descriptor to be filled |
||
type(tem_path_type), | intent(in) | :: | pathFirst(:) |
first treeID path in every process |
||
type(tem_path_type), | intent(in) | :: | pathLast(:) |
last treeID path in every process |
||
type(tem_comm_env_type), | intent(in) | :: | proc |
process |
recursive subroutine create_allParentNeighbors( & & targetID, level, stencil, tree, levelDesc, & & pathFirst, pathLast, proc ) ! -------------------------------------------------------------------- ! !> requested element position (child element) in LevelDesc elem list integer(kind=long_k), intent(in) :: targetID !> requested element level integer, intent(in) :: level !> first treeID path in every process type(tem_path_type), intent(in) :: pathFirst(:) !> last treeID path in every process type(tem_path_type), intent(in) :: pathLast(:) !> tree information type(treelmesh_type), intent(in) :: tree !> the level descriptor to be filled type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:) !> process type(tem_comm_env_type), intent(in) :: proc !> current stencil definition type(tem_stencilHeader_type), intent(in) :: stencil ! -------------------------------------------------------------------- ! integer(kind=long_k) :: parentID, neighID ! current tree ID integer :: coarserLevel, cPos, parentNesting, addedPos, iStencilElem integer :: neighIDpos ! -------------------------------------------------------------------- ! ! exit if we have reached the minimal level if ( level == tree%global%minlevel ) return ! Get the parent of the current treeID parentID = tem_parentOf( targetID ) ! ... and identify the parent parentNesting = -1 call identify_elements( TreeID = parentID, & & tree = tree, & & pathFirst = pathFirst, & & pathLast = pathLast, & & levelDesc = levelDesc, & & elemPos = cPos, & & proc = proc, & & stencil = stencil, & & nesting = -1 ) if ( cPos <= 0 ) then write(dbgUnit(3),*) ' Element not found: ', parentID write(dbgUnit(3),*) ' cPos: ', cPos end if ! identify the stencil neighbors of the parent. ! Here we should identify the fromCoarser interpolation stencil neighbors ! instead of the compute stencil neighbors coarserLevel = level - 1 call identify_stencilNeigh( iElem = cPos, & & iLevel = coarserLevel, & & tree = tree, & & iStencil = 1, & & pathFirst = pathFirst, & & pathLast = pathLast, & & levelDesc = levelDesc, & & proc = proc, & & stencil = stencil, & & nesting = parentNesting ) ! adding neighs of neighs for interpolation at ML ! this is for the stencil of the fine ghost, composed by coarse fluid do iStencilElem = 1, stencil%QQN neighIDpos = levelDesc(coarserLevel)%elem%stencil%val(cPos) & & %val(1)%tIDpos(iStencilElem) if( neighIDpos > 0 ) then neighID = & & levelDesc( coarserLevel )%elem%neighID%val(cPos)%val(neighIDpos) ! This call might add new halo elements if ( neighID > 0_long_k ) then call append( me = levelDesc( coarserLevel )%require, & & val = neighID, & & pos = addedPos ) end if end if ! neighIDpos > 0 enddo end subroutine create_allParentNeighbors