tem_cano_checkNeigh Subroutine

public subroutine tem_cano_checkNeigh(me, inTree, countElems, map2global)

Check for neighboring elements to a given point.

This routine can be used to find close by elements for points that just lie outside the domain and for which values can be extrapolated, though they do not strictly reside in the computational domain described by inTree. The closest neighboring element is identified globally and added to the map2global if any is found. This only works for points, and none of the points in the given list should have been found inside the domain.

Arguments

Type IntentOptional Attributes Name
type(tem_canonicalND_type), intent(in) :: me(:)

canonicalND objects to check, only points will be considered

type(treelmesh_type), intent(in) :: inTree

Global tree

integer, intent(inout) :: countElems(globalMaxLevels)

How many elements there will be for each level in the track

type(dyn_intarray_type), intent(inout) :: map2global

growing array for the map2global


Calls

proc~~tem_cano_checkneigh~~CallsGraph proc~tem_cano_checkneigh tem_cano_checkNeigh proc~tem_idofcoord tem_IdOfCoord proc~tem_cano_checkneigh->proc~tem_idofcoord proc~tem_firstidatlevel tem_FirstIdAtLevel proc~tem_cano_checkneigh->proc~tem_firstidatlevel interface~append~9 append proc~tem_cano_checkneigh->interface~append~9 proc~tem_coordofid tem_CoordOfId proc~tem_cano_checkneigh->proc~tem_coordofid proc~tem_posofid tem_PosOfId proc~tem_cano_checkneigh->proc~tem_posofid proc~tem_coordofreal tem_CoordOfReal proc~tem_cano_checkneigh->proc~tem_coordofreal proc~tem_baryofid tem_BaryOfId proc~tem_cano_checkneigh->proc~tem_baryofid mpi_allreduce mpi_allreduce proc~tem_cano_checkneigh->mpi_allreduce proc~tem_levelof tem_LevelOf proc~tem_cano_checkneigh->proc~tem_levelof proc~append_da_label append_da_label interface~append~9->proc~append_da_label proc~append_da_veclabel append_da_veclabel interface~append~9->proc~append_da_veclabel proc~tem_coordofid->proc~tem_levelof proc~tem_pathof tem_PathOf proc~tem_posofid->proc~tem_pathof proc~tem_pathcomparison tem_PathComparison proc~tem_posofid->proc~tem_pathcomparison proc~tem_baryofid->proc~tem_coordofid proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_baryofid->proc~tem_elemsizelevel interface~sortedposofval~5 sortedposofval proc~append_da_label->interface~sortedposofval~5 interface~expand~9 expand proc~append_da_label->interface~expand~9 proc~append_da_veclabel->interface~expand~9 proc~sortposofval_label sortposofval_label interface~sortedposofval~5->proc~sortposofval_label proc~expand_da_label expand_da_label interface~expand~9->proc~expand_da_label

Called by

proc~~tem_cano_checkneigh~~CalledByGraph proc~tem_cano_checkneigh tem_cano_checkNeigh proc~tem_shape_subtreefromgeominters tem_shape_subTreeFromGeomInters proc~tem_shape_subtreefromgeominters->proc~tem_cano_checkneigh proc~tem_shape2subtree tem_shape2subTree proc~tem_shape2subtree->proc~tem_shape_subtreefromgeominters proc~tem_create_subtree_of tem_create_subTree_of proc~tem_create_subtree_of->proc~tem_shape2subtree proc~tem_init_convergence tem_init_convergence proc~tem_init_convergence->proc~tem_create_subtree_of proc~tem_init_tracker_subtree tem_init_tracker_subTree proc~tem_init_tracker_subtree->proc~tem_create_subtree_of proc~tem_write_debugmesh tem_write_debugMesh proc~tem_write_debugmesh->proc~tem_create_subtree_of proc~tem_create_subtree_of_st_funlist tem_create_subTree_of_st_funList proc~tem_create_subtree_of_st_funlist->proc~tem_create_subtree_of

Contents

Source Code


Source Code

  subroutine tem_cano_checkNeigh(me, inTree, countElems, map2global)
    ! --------------------------------------------------------------------------
    !> canonicalND objects to check, only points will be considered
    type(tem_canonicalND_type ),intent(in) :: me(:)
    !> Global tree
    type(treelmesh_type), intent(in) :: inTree
    !> How many elements there will be for each level in the track
    integer, intent( inout ) :: countElems( globalMaxLevels )
    !> growing array for the map2global
    type(dyn_intArray_type), intent(inout) :: map2global
    ! --------------------------------------------------------------------------
    integer :: iCano, tLevel, elemPos, dPos, maxLevel, iQQN
    integer :: closest_neigh(size(me))
    logical :: wasAdded
    integer :: coordOfId(4)
    integer :: nCanoNDs
    integer :: iError
    integer(kind=long_k) :: treeID, tOffset, neighID
    integer :: minproc(size(me))
    real(kind=rk) :: mindist(size(me)), dist
    real(kind=rk) :: minglobdist(size(me))
    real(kind=rk) :: elemcenter(3)
    ! --------------------------------------------------------------------------
    maxLevel = inTree%global%maxLevel
    mindist = huge(dist)
    closest_neigh = 0
    nCanoNDs = size(me)

    ! Create subTree by intersecting canoND objects with the inTree
    do iCano=1,nCanoNDs
      if (trim(me(iCano)%kind) == 'point') then
        treeID = tem_IdOfCoord( tem_CoordOfReal(inTree,                     &
          &                                     me(iCano)%point%coord(1:3), &
          &                                     maxLevel)                   )
        ! add a neighbor element of this point if any of the neighbor exist to
        ! interpolate to a point
        coordOfId = tem_CoordOfId( treeID )
        tOffset = tem_FirstIdAtLevel( coordOfId(4) )
        directionLoop: do iQQN = 1, qQQQ
          neighID = tem_IdOfCoord(                          &
            &         [ coordOfId(1) + qOffset( iQQN, 1 ),  &
            &           coordOfId(2) + qOffset( iQQN, 2 ),  &
            &           coordOfId(3) + qOffset( iQQN, 3 ),  &
            &           coordOfId(4) ], tOffset)
          elemPos = tem_PosOfId(neighID, inTree%treeID)
          if (elemPos > 0) then
            elemcenter = tem_BaryOfId(inTree, neighID)
            dist = sqrt(  (me(iCano)%point%coord(1) - elemcenter(1))**2 &
              &         + (me(iCano)%point%coord(2) - elemcenter(2))**2 &
              &         + (me(iCano)%point%coord(3) - elemcenter(3))**2 )
            if (dist < mindist(iCano)) then
              mindist(iCano) = dist
              closest_neigh(iCano) = elemPos
            else if (dist <= mindist(iCano)) then
              ! Make sure to pick the first element in the SFC ordering
              ! as the nearest neighbor, if there are multiple with the
              ! same distance (<= comparison to avoid == comparison for
              ! reals)
              closest_neigh(iCano) = min(elemPos, closest_neigh(iCano))
            end if
          end if
        end do directionLoop
      end if
    end do

    call MPI_Allreduce( mindist, minglobdist, nCanoNDs,             &
      &                 rk_mpi, MPI_MIN, inTree%global%comm, iError )

    minproc = inTree%global%nParts
    do iCano=1,nCanoNDs
      if ( minglobdist(iCano) < huge(dist) ) then
        if ( mindist(iCano) <= minglobdist(iCano) ) then
          minproc(iCano) = inTree%global%myPart
        end if
      end if
    end do

    call MPI_Allreduce( MPI_IN_PLACE, minproc, nCanoNDs,                 &
      &                 MPI_INTEGER, MPI_MIN, inTree%global%comm, iError )

    do iCano=1,nCanoNDs
      if (trim(me(iCano)%kind) == 'point') then
        found_elem: if (minproc(iCano) == inTree%global%myPart) then
          ! append the position in inTree to the map (note that already existing
          ! ones are omitted)
          call append( me       = map2global,           &
            &          pos      = dpos,                 &
            &          val      = closest_neigh(iCano), &
            &          wasAdded = wasAdded              )

          ! Count up if it was added
          if (wasAdded) then
            tLevel   = tem_levelOf( inTree%treeID(closest_neigh(iCano)) )
            countElems( tLevel ) = countElems( tLevel ) + 1
          end if ! wasAdded
        end if found_elem
      end if
    end do

  end subroutine tem_cano_checkNeigh