tem_calc_vrtx_coord_noqval Subroutine

private subroutine tem_calc_vrtx_coord_noqval(tree, vrtx, subTree, vrtxTreeID)

Run over all 8 vertices for each element in the treeID list, calculate its coordinates and add its position to the map.

Arguments

Type IntentOptional Attributes Name
type(treelmesh_type), intent(in) :: tree

fluid mesh

type(tem_vrtx_type), intent(inout) :: vrtx

Vertex data

type(tem_subTree_type), intent(in), optional :: subTree

optional subTree information

integer(kind=long_k), intent(out), allocatable :: vrtxTreeID(:)

boundary information incl. q-Values


Calls

proc~~tem_calc_vrtx_coord_noqval~~CallsGraph proc~tem_calc_vrtx_coord_noqval tem_calc_vrtx_coord_noqval proc~tem_treeidfrom_subtree tem_treeIDfrom_subTree proc~tem_calc_vrtx_coord_noqval->proc~tem_treeidfrom_subtree proc~tem_idofcoord tem_IdOfCoord proc~tem_calc_vrtx_coord_noqval->proc~tem_idofcoord proc~tem_coordofid tem_CoordOfId proc~tem_calc_vrtx_coord_noqval->proc~tem_coordofid proc~tem_init_vrtx_prop tem_init_vrtx_prop proc~tem_calc_vrtx_coord_noqval->proc~tem_init_vrtx_prop proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof interface~init~22 init proc~tem_init_vrtx_prop->interface~init~22 proc~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real

Called by

proc~~tem_calc_vrtx_coord_noqval~~CalledByGraph proc~tem_calc_vrtx_coord_noqval tem_calc_vrtx_coord_noqval proc~tem_calc_vrtx_coord tem_calc_vrtx_coord proc~tem_calc_vrtx_coord->proc~tem_calc_vrtx_coord_noqval proc~hvs_output_init hvs_output_init proc~hvs_output_init->proc~tem_calc_vrtx_coord proc~tem_init_tracker tem_init_tracker proc~tem_init_tracker->proc~hvs_output_init

Contents


Source Code

  subroutine tem_calc_vrtx_coord_noqval( tree, vrtx, subTree, vrtxTreeID )
    ! ---------------------------------------------------------------------------
    !> fluid mesh
    type(treelmesh_type),     intent(in) :: tree
    !> Vertex data
    type(tem_vrtx_type), intent(inout) :: vrtx
    !> optional subTree information
    type(tem_subTree_type), optional, intent(in) :: subTree
    !> boundary information incl. q-Values
    integer(kind=long_k), allocatable, intent(out) :: vrtxTreeID(:)
    ! ---------------------------------------------------------------------------
    ! counters
    integer :: iVrtx, iElem
    integer :: local_nElems
    integer(kind=long_k) :: vrtxID
    integer :: elemCoord(4)
    integer :: locVrtx(4)
    integer :: vrtxAnchor(4)
    integer :: iLevel
    integer(kind=long_k), allocatable :: treeID(:)
    ! ---------------------------------------------------------------------------

    vrtx%maxVertices = 0

    if (allocated(treeID)) deallocate(treeID)
    if (allocated(vrtx%refine)) deallocate(vrtx%refine)

    if( present( subTree ))then
      local_nElems = subTree%nElems
      allocate( treeID( local_nElems ))
      call tem_treeIDfrom_subTree( subTree, tree, treeID, (/1,local_nElems/) )
    else
      local_nElems = tree%nElems
      allocate( treeID( local_nElems ))
      treeID = tree%treeID
    end if
    vrtx%maxVertices = 8*local_nElems
    allocate( vrtx%refine( local_nElems ))
    vrtx%refine = .false.

    ! initialize the vertex type
    call tem_init_vrtx_prop( vrtx = vrtx)

    ! allocate the list of all vrtxTreeIDs including dublicates
    allocate( vrtxTreeID( vrtx%maxVertices ))

    write(logUnit(6),*) 'DEBUG: Filling the global vrtxTreeID ...'

    ! map the treeIDs to those of a tree with a bounding cube length
    ! twice as big -> treeIDs correspond to those 1 refinement level
    ! higher in the bigger tree
    !                       -------------------------------
    !                       |              |              |
    !                       |              |              |
    !                       |              |              |
    !                       |              |              |
    !                       |              |              |
    ! ----------------      -------------------------------
    ! |      |       |      |      |       |              |
    ! |      |       |      |      |       |              |
    ! ----------------  --> ----------------              |
    ! |      |       |      |      |       |              |
    ! |      |       |      |      |       |              |
    ! ----------------      -------------------------------
    !

    do iVrtx=1,8
      locVrtx = tem_coordOfID(int(iVrtx, kind=long_k))
      do iElem = 1, local_nElems
        elemCoord = tem_coordOfID(treeID(iElem))

        ! since the coordinates of the individual vertices are on level 1
        ! the level for the vrtxAnchor is increased by 1 matching the
        ! requirements of the new tree (bounding cube size)
        vrtxAnchor = elemCoord + locVrtx
        ! retransforming the coords to the treeID on the 'new tree'
        vrtxID = tem_IDofCoord(vrtxAnchor)
        ! get the treeID on the highest refinement level possible as
        ! a unique identifier
        !NEC$ NOVECTOR
        do iLevel=vrtxAnchor(4)+1,globalMaxLevels
          vrtxID = vrtxID*8_long_k + 1_long_k
        end do
        vrtxTreeID( (iElem-1)*8 + iVrtx ) = vrtxID
      end do
    end do

    write(logUnit(6),*) 'DEBUG: Filled it.'


  end subroutine tem_calc_vrtx_coord_noqval