tem_calc_vrtx_coord Subroutine

public subroutine tem_calc_vrtx_coord(tree, vrtx, subTree, boundary, useQVal)

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

type(tem_BC_prop_type), intent(in), optional :: boundary

boundary information incl. q-Values

logical, intent(in), optional :: useQVal

use the qValue information?


Calls

proc~~tem_calc_vrtx_coord~~CallsGraph proc~tem_calc_vrtx_coord tem_calc_vrtx_coord proc~tem_treeidfrom_subtree tem_treeIDfrom_subTree proc~tem_calc_vrtx_coord->proc~tem_treeidfrom_subtree proc~tem_calc_vrtx_coord_noqval tem_calc_vrtx_coord_noqval proc~tem_calc_vrtx_coord->proc~tem_calc_vrtx_coord_noqval proc~tem_idofcoord tem_IdOfCoord proc~tem_calc_vrtx_coord->proc~tem_idofcoord proc~tem_init_vrtx_prop tem_init_vrtx_prop proc~tem_calc_vrtx_coord->proc~tem_init_vrtx_prop interface~append~11 append proc~tem_calc_vrtx_coord->interface~append~11 proc~tem_invertrealrkarray tem_invertRealRkArray proc~tem_calc_vrtx_coord->proc~tem_invertrealrkarray proc~tem_unify_vrtx tem_unify_vrtx proc~tem_calc_vrtx_coord->proc~tem_unify_vrtx proc~tem_coordofid tem_CoordOfId proc~tem_calc_vrtx_coord->proc~tem_coordofid proc~tem_calc_vrtxof_qval tem_calc_vrtxOf_qVal proc~tem_calc_vrtx_coord->proc~tem_calc_vrtxof_qval proc~qsort_vrtx qsort_vrtx proc~tem_calc_vrtx_coord->proc~qsort_vrtx proc~tem_calc_vrtx_coord_noqval->proc~tem_treeidfrom_subtree proc~tem_calc_vrtx_coord_noqval->proc~tem_idofcoord proc~tem_calc_vrtx_coord_noqval->proc~tem_init_vrtx_prop proc~tem_calc_vrtx_coord_noqval->proc~tem_coordofid interface~init~22 init proc~tem_init_vrtx_prop->interface~init~22 proc~tem_appendint2darray tem_appendInt2dArray interface~append~11->proc~tem_appendint2darray proc~tem_appendintlist tem_appendIntList interface~append~11->proc~tem_appendintlist proc~tem_appenddp2darray tem_appendDp2dArray interface~append~11->proc~tem_appenddp2darray proc~tem_appendintlong2darray tem_appendIntLong2dArray interface~append~11->proc~tem_appendintlong2darray proc~tem_appendintlong1darray tem_appendIntLong1dArray interface~append~11->proc~tem_appendintlong1darray proc~tem_appendsp1darray tem_appendSp1dArray interface~append~11->proc~tem_appendsp1darray proc~tem_appendsp2darray tem_appendSp2dArray interface~append~11->proc~tem_appendsp2darray proc~tem_appendlonglist tem_appendLongList interface~append~11->proc~tem_appendlonglist proc~tem_appenddp1darray tem_appendDp1dArray interface~append~11->proc~tem_appenddp1darray proc~tem_appendintlongarrayto1darray tem_appendIntLongArrayTo1dArray interface~append~11->proc~tem_appendintlongarrayto1darray proc~tem_appendint1darray tem_appendInt1dArray interface~append~11->proc~tem_appendint1darray proc~tem_unify_vrtx->interface~append~11 proc~tem_originofid tem_originOfId proc~tem_unify_vrtx->proc~tem_originofid proc~tem_posofid tem_PosOfId proc~tem_unify_vrtx->proc~tem_posofid proc~tem_posoflong tem_posOfLong proc~tem_unify_vrtx->proc~tem_posoflong proc~tem_unify_vrtx->interface~init~22 interface~destroy~22 destroy proc~tem_unify_vrtx->interface~destroy~22 proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof proc~tem_elemsize tem_ElemSize proc~tem_calc_vrtxof_qval->proc~tem_elemsize proc~tem_baryofid tem_BaryOfId proc~tem_calc_vrtxof_qval->proc~tem_baryofid proc~qsort_vrtx->proc~qsort_vrtx proc~partition partition proc~qsort_vrtx->proc~partition

Called by

proc~~tem_calc_vrtx_coord~~CalledByGraph proc~tem_calc_vrtx_coord tem_calc_vrtx_coord 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


Source Code

  subroutine tem_calc_vrtx_coord( tree, vrtx, subTree, boundary, useQVal )
    ! ---------------------------------------------------------------------------
    !> 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
    type(tem_BC_prop_type), optional, intent(in) :: boundary
    !> use the qValue information?
    logical, optional, intent(in) :: useQVal
    ! ---------------------------------------------------------------------------
    ! counters
    integer :: iVrtx, iElem
    integer :: local_nElems
    ! store all treeIDs for the vertices of each element in vrtxTreeID
    ! in the order
    ! -----------------------------------------------------------------
    ! | vrtx1 vrtx2 vrtx3 vrtx4 vrtx5 vrtx6 vrtx7 vrtx8 |     ...     |
    ! -----------------------------------------------------------------
    !                 iElem = 1                            iElem =...
    integer(kind=long_k), allocatable :: vrtxTreeID(:)
    integer(kind=long_k), allocatable :: sortedVrtxTreeID(:)
    integer(kind=long_k) :: vrtxID
    integer :: elemCoord(4)
    integer :: locVrtx(4)
    integer :: vrtxAnchor(4)
    integer :: iLevel
    ! tree with bounding cube length twice as big as in tree (treeID array will
    ! not be filled!!!!)
    type(treelmesh_type) :: bigTree
    type(tem_property_type), pointer :: tree_property(:) => NULL()
    integer(kind=long_k), allocatable :: treeID(:)
    ! type(tem_global_type) :: global
    ! counters
    integer :: globCounter
    integer :: uniqueCounter
    integer :: nElemsQVal
    integer :: iBCElem
    real(kind=rk) :: coord(3)
    logical :: local_useQVal
    ! ---------------------------------------------------------------------------

    if( present( useQVal ))then
      local_useQVal = useQVal
    else
      local_useQVal = .false.
    end if

    if( present( subTree ))then
      local_nElems = subTree%nElems
    else
      local_nElems = tree%nElems
    end if

    uniqueCounter = 0

    if (local_useQval) then
      vrtx%maxVertices = 0
      if (associated(tree_property)) deallocate(tree_property)
      if (allocated(treeID)) deallocate(treeID)
      if (allocated(vrtx%refine)) deallocate(vrtx%refine)

      if( present( subTree ))then
        ! global = subTree%global
        allocate( tree_property( subTree%global%nProperties ))
        tree_property = subTree%Property
        allocate( treeID( local_nElems ))
        allocate( vrtx%refine( local_nElems ))
        call tem_treeIDfrom_subTree( subTree, tree, treeID, (/1,local_nElems/) )
        ! possible q-Values attached calculate max number of vertices (8 per
        ! 'normal' element, 20 per element with q-Values) and allocate the array
        do iElem = 1, local_nElems
          if( btest(subTree%elemPropertyBits( iElem ), prp_hasQVal )) then
            vrtx%maxVertices = vrtx%maxVertices + 20
            vrtx%refine( iElem ) = .true.
          else
            vrtx%maxVertices = vrtx%maxVertices + 8
            vrtx%refine( iElem ) = .false.
          end if
        end do
      else
        ! global = tree%global
        allocate( tree_property( tree%global%nProperties ))
        tree_property = tree%Property
        allocate( treeID( local_nElems ))
        allocate( vrtx%refine( local_nElems ))
        treeID = tree%treeID
        ! possible q-Values attached calculate max number of vertices (8 per
        ! 'normal' element, 20 per element with q-Values) and allocate the array
        do iElem = 1, local_nElems
          if (btest(tree%elemPropertyBits( iElem ), prp_hasQVal) ) then
            vrtx%maxVertices = vrtx%maxVertices + 20
            vrtx%refine( iElem ) = .true.
          else
            vrtx%maxVertices = vrtx%maxVertices + 8
            vrtx%refine( iElem ) = .false.
          end if
        end do
      end if

      ! 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 ...'

      ! initialize counters
      globCounter = 0

      ! 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
      !                       -------------------------------
      !                       |              |              |
      !                       |              |              |
      !                       |              |              |
      !                       |              |              |
      !                       |              |              |
      ! ----------------      -------------------------------
      ! |      |       |      |      |       |              |
      ! |      |       |      |      |       |              |
      ! ----------------  --> ----------------              |
      ! |      |       |      |      |       |              |
      ! |      |       |      |      |       |              |
      ! ----------------      -------------------------------
      !

      nElemsQVal = 0
      do iElem = 1, local_nElems
        ! if element has q-Values it has to be refined once
        if( vrtx%refine( iElem ))then
           if( present( subTree ) .and. .not. subTree%useGlobalMesh )then
             do iBCElem = 1, size( tree%property(2)%elemID )
               if( tree%property(2)%elemID( iBCElem ) .eq.                       &
                 &                              subTree%map2global( iElem )) then
                 nElemsQVal = iBCElem
                 exit
               end if
             end do
           else
             nElemsQVal = nElemsQVal + 1
           end if
          ! calculate the vertices for the element incl. q-Values
          do iVrtx = 1, 20
            ! check if the q-Value for vertex iVrtx is greater than 0.5 (this
            ! means the point might be shared between elements) or the q-Value
            ! is -1.0 (this means that no q-Value is set in the corresponding
            ! direction)
            if( ( abs(boundary%qVal( vrtxMap(5, iVrtx), nElemsQVal) -0.5_rk) .le.&
              &                                            eps)              .or.&
!              & ( abs(boundary%qVal( vrtxMap(5, iVrtx), nElemsQVal) +0.5_rk) .gt.&
!              &                                       1.0_rk + eps)          .or.&
              & ( abs(boundary%qVal( vrtxMap(5, iVrtx), nElemsQVal) +1.0_rk) .le.&
              &                                            eps) )then

              ! for the 8 corners get them from the treeID
              if( iVrtx .le. 8)then
                elemCoord = tem_coordOfID(treeID(iElem))
              else ! for the 12 intermediate get them from the refined treeID
                ! refine the element by 1 level
                elemCoord = tem_coordOfID( treeID(iElem)*8_long_k + 1_long_k )
              end if
              ! 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 + vrtxMap( 1:4, iVrtx )
              ! 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
              do iLevel=vrtxAnchor(4)+1,globalMaxLevels
                vrtxID = vrtxID*8_long_k + 1_long_k
              end do
            else ! q-Value .ne. 0.5 (assume vertex is unique)
              uniqueCounter = uniqueCounter - 1
              vrtxID = uniqueCounter
              ! calculate the vertex based on the q-Value and append it to the
              ! growing array of vertices
              coord = tem_calc_vrtxOf_qVal(                                      &
                &         treeID = treeID(iElem),                                &
                &         tree   = tree,                                         &
                &         qVal   = boundary%qVal( vrtxMap(5, iVrtx), nElemsQVal),&
                &         iVrtx  = iVrtx )
              call append( me = vrtx%coord, val = coord )
            end if
            vrtxTreeID( globCounter+iVrtx ) = vrtxID
          end do
          globCounter = globCounter + 20
        else ! no q-Values
          elemCoord = tem_coordOfID(treeID(iElem))
          do iVrtx=1,8
            locVrtx = tem_coordOfID(int(iVrtx, kind=long_k))
            ! 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
            do iLevel=vrtxAnchor(4)+1,globalMaxLevels
              vrtxID = vrtxID*8_long_k + 1_long_k
            end do
            vrtxTreeID( globCounter+iVrtx ) = vrtxID
          end do
          globCounter = globCounter + 8
        end if
      end do

      write(logUnit(6),*) 'DEBUG: Filled it.'
    else
      call tem_calc_vrtx_coord_noqval( tree, vrtx, subTree, vrtxTreeID )
    end if

    allocate( sortedVrtxTreeID( vrtx%maxVertices ))
    sortedVrtxTreeID = vrtxTreeID

    write(logUnit(6),*) 'DEBUG: Start sorting ...'

    ! sort the treeID array
    call qsort_vrtx( sortedVrtxTreeID )

    write(logUnit(6),*) 'DEBUG: Done sorting, start inverting coords ...'

    ! in case q-Values are present reorganize the growing array of coords
    ! such that it is in the correct order
    if( (-1)*uniqueCounter .gt. 1 )then
      call tem_invertRealRkArray( me = vrtx%coord%val, nElems = vrtx%coord%nVals )
    end if

    write(logUnit(6),*) 'DEBUG: Done inverting, start to unify ...'

    ! redefine the tree bounding cube size
    bigTree%global%origin = tree%global%origin
    bigTree%global%BoundingCubeLength = 2.0_rk * tree%global%BoundingCubeLength

    ! make sorted array vrtxTreeID unique and map the elements to the right
    ! vertex real coordinates
    call tem_unify_vrtx( inList   = sortedVrtxTreeID,                          &
      &                  origList = vrtxTreeID,                                &
      &                  coord    = vrtx%coord,                                &
      &                  map      = vrtx%map2global,                           &
      &                  tree     = bigTree,                                   &
      &                  nElems   = local_nElems,                              &
      &                  nUnique  = (-1)*uniqueCounter,                        &
      &                  refine = vrtx%refine )

    write(logUnit(6),*) 'DEBUG: Done unifying.'

    ! update the number of calculated vertices
    vrtx%nVertices = vrtx%coord%nVals

  end subroutine tem_calc_vrtx_coord