This subroutine identifies the parent treelm elements of the surface data points.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_surfData_type), | intent(inout) | :: | me |
datatype to store the surface information |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc |
the level descriptor incl. ghost and halo elements as well as the communicator information on the level iLevel |
||
type(treelmesh_type), | intent(in) | :: | globTree |
global Tree information |
||
integer, | intent(in) | :: | iLevel |
the current level |
subroutine tem_init_surfData( me, levelDesc, globTree, iLevel ) ! --------------------------------------------------------------------------- !> datatype to store the surface information type( tem_surfData_type ), intent(inout) :: me !> the level descriptor incl. ghost and halo elements as well as the !! communicator information on the level iLevel type( tem_levelDesc_type ), intent(inout) :: levelDesc !> global Tree information type( treelmesh_type ), intent(in) :: globTree !> the current level integer, intent(in) :: iLevel ! --------------------------------------------------------------------------- ! counters integer :: iCoord, iVal integer :: pos integer( kind=long_k ) :: tmpTreeID type( dyn_longArray_type ) :: tmpTreeIDs integer, allocatable :: posInTotal(:) integer, allocatable :: tmpPos(:) real( kind=rk ) :: huge_real type(tem_path_type) :: tmpPath type(tem_path_type) :: minHaloPath type(tem_path_type) :: maxHaloPath integer :: nFluids integer :: firstHalo integer :: lastHalo logical :: wasAdded ! --------------------------------------------------------------------------- huge_real = huge(huge_real) if( .not. allocated(me%parentIDs(iLevel)%ptrs))then ! allocate the array of parent positions in the local tree id list allocate( me%parentIDs(iLevel)%ptrs( me%nUniquePoints_total )) end if nFluids = globTree%nElems ! correct the halo position (in case of serial runs no halos available) if( globTree%global%nParts > 1 )then firstHalo = levelDesc%elem%nElems( eT_fluid ) & & + levelDesc%elem%nElems( eT_ghostFromCoarser ) & & + levelDesc%elem%nElems( eT_ghostFromFiner ) & & + 1 lastHalo = levelDesc%elem%nElems( eT_fluid ) & & + levelDesc%elem%nElems( eT_ghostFromCoarser ) & & + levelDesc%elem%nElems( eT_ghostFromFiner ) & & + levelDesc%elem%nElems( eT_halo ) ! get the path of the min and max halo minHaloPath = tem_PathOf( levelDesc%total( firstHalo )) maxHaloPath = tem_PathOf( levelDesc%total( lastHalo )) end if allocate( tmpPos( me%nUniquePoints_total )) ! for all coordinates do iCoord = 1, me%nUniquePoints_total pos = 0 ! if( any( me%pointCoords( (iCoord-1)*3+1:(iCoord-1)*3+3) < & ! & huge_real ))then ! ... get the parent treeID (first integer from real coordinates then ! treeID) tmpTreeID = tem_IdOfCoord( & & tem_CoordOfReal( globTree, & & me%pointCoords((iCoord-1)*3+1:(iCoord-1)*3+3 ), & & iLevel )) call append( me = tmpTreeIDs, & & val = tmpTreeID, & & pos = tmpPos( iCoord ), & & wasAdded = wasAdded ) end do allocate( posInTotal( tmpTreeIDs%nVals )) !$omp parallel !$omp do private( tmpPath, pos ) do iVal=1, tmpTreeIDs%nVals pos = 0 tmpPath = tem_PathOf( tmpTreeIDs%val(iVal) ) ! compare the paths to check wether the surface point is located ! on this proc if( tem_PathComparison( tmpPath, globTree%pathList( 1 )) > -1 .and. & & tem_PathComparison( tmpPath, globTree%pathList( nFluids )) < 1) then ! search in the fluid elements pos = tem_treeIDinTotal( tmpTreeIDs%val(iVal), levelDesc, eT_fluid ) else if( tem_PathComparison( tmpPath, minHaloPath ) > -1 .and. & & tem_PathComparison( tmpPath, maxHaloPath ) < 1 )then ! surf point is not located on this proc but might be a halo pos = tem_treeIDinTotal( tmpTreeIDs%val(iVal), levelDesc, eT_halo) end if posInTotal( iVal ) = pos end do !$omp end do ! end if !$omp do do iCoord = 1, me%nUniquePoints_total ! ... store the new coordinate me%parentIDs(iLevel)%ptrs(iCoord) = posInTotal( tmpPos(iCoord) ) if( me%parentIDs(iLevel)%ptrs(iCoord) > 0 )then ! ... and set the corresponding property bits levelDesc%property( me%ParentIDs(iLevel)%ptrs( iCoord )) = & & ibset( levelDesc%property( me%ParentIDs(iLevel)%ptrs( iCoord )), & & prp_hasIBM ) end if end do !$omp end do !$omp end parallel end subroutine tem_init_surfData