tem_unify_vrtx Subroutine

private subroutine tem_unify_vrtx(inList, origList, coord, map, tree, nElems, nUnique, refine)

This subroutine takes the sorted list as an input and unifies its entries the result is used to create a unique array of vertex coordinates and a map for the 8 vertices of each element.

Arguments

Type IntentOptional Attributes Name
integer(kind=long_k), intent(inout), allocatable :: inList(:)
integer(kind=long_k), intent(inout), allocatable :: origList(:)
type(grw_real2darray_type) :: coord
integer, intent(inout), allocatable :: map(:,:)
type(treelmesh_type), intent(in) :: tree
integer, intent(in) :: nElems
integer, intent(in) :: nUnique

number of unique vertices (from q-Values)

logical, intent(in) :: refine(:)

Calls

proc~~tem_unify_vrtx~~CallsGraph proc~tem_unify_vrtx tem_unify_vrtx interface~append~2 append proc~tem_unify_vrtx->interface~append~2 interface~destroy~25 destroy proc~tem_unify_vrtx->interface~destroy~25 interface~init~24 init proc~tem_unify_vrtx->interface~init~24 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_appenddp1darray tem_appendDp1dArray interface~append~2->proc~tem_appenddp1darray proc~tem_appenddp2darray tem_appendDp2dArray interface~append~2->proc~tem_appenddp2darray proc~tem_appendint1darray tem_appendInt1dArray interface~append~2->proc~tem_appendint1darray proc~tem_appendint2darray tem_appendInt2dArray interface~append~2->proc~tem_appendint2darray proc~tem_appendintlist tem_appendIntList interface~append~2->proc~tem_appendintlist proc~tem_appendintlong1darray tem_appendIntLong1dArray interface~append~2->proc~tem_appendintlong1darray proc~tem_appendintlong2darray tem_appendIntLong2dArray interface~append~2->proc~tem_appendintlong2darray proc~tem_appendintlongarrayto1darray tem_appendIntLongArrayTo1dArray interface~append~2->proc~tem_appendintlongarrayto1darray proc~tem_appendlonglist tem_appendLongList interface~append~2->proc~tem_appendlonglist proc~tem_appendsp1darray tem_appendSp1dArray interface~append~2->proc~tem_appendsp1darray proc~tem_appendsp2darray tem_appendSp2dArray interface~append~2->proc~tem_appendsp2darray proc~destroy_ga2d_real destroy_ga2d_real interface~destroy~25->proc~destroy_ga2d_real proc~init_ga2d_real init_ga2d_real interface~init~24->proc~init_ga2d_real proc~tem_coordofid tem_CoordOfId proc~tem_originofid->proc~tem_coordofid proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_originofid->proc~tem_elemsizelevel proc~tem_pathcomparison tem_PathComparison proc~tem_posofid->proc~tem_pathcomparison proc~tem_pathof tem_PathOf proc~tem_posofid->proc~tem_pathof proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof

Called by

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

Source Code

  subroutine tem_unify_vrtx( inList, origList, coord, map, tree, nElems,       &
    &                        nUnique, refine )
    ! ---------------------------------------------------------------------------
    integer(kind=long_k), allocatable, intent(inout) :: inList(:)
    integer(kind=long_k), allocatable, intent(inout) :: origList(:)
    type(grw_real2dArray_type) :: coord
    integer, allocatable, intent(inout) :: map(:,:)
    type(treelmesh_type), intent(in) :: tree
    integer, intent(in) :: nElems
    !> number of unique vertices (from q-Values)
    integer, intent(in) :: nUnique
    logical, intent(in) :: refine(:)
    ! ---------------------------------------------------------------------------
    ! counters
    integer :: count1, count2
    integer :: iElem, iVrtx, pos, counter
    type( grw_longArray_type ) :: unique
    real(kind=rk) :: tmp_vrtx(3)
    ! ---------------------------------------------------------------------------
    count1 = nUnique+1
    count2 = nUnique+1

    call init( me = unique, length = 10 )

    ! append the unique treeIDs to the unique array
    do iElem = 1, nUnique
      call append( me  = unique,                                               &
        &          val = inList( iElem ))
    end do

    ! at first make the list of vertex treeIDs and their coordinates unique
    do while( count1 .lt. size( inList))
      count1 = count2
      call append( me  = unique,                                               &
        &          val = inList( count1 ))

      ! get the real coordinates of the unique treeID ...
      tmp_vrtx = tem_originOfId( tree, inList( count1 ))
      ! ... and store them in by definition unique array
      call append( me  = coord,                                                &
        &          val = tmp_vrtx )

      do while(( inList( count1 ) .eq. inList( count2 ))                       &
        &        .and. count2 .lt. size( inList ))
        count2 = count2+1
      end do
    end do

    deallocate( inList )

    counter = 0
    ! map the original treeID list to the unique and by this to the unique coord
    ! array
    do iElem = 1, nElems
      if( refine( iElem ))then
        do iVrtx = 1, 20
          ! if the element has a valid treeID search for it using the regular
          ! posOfID
          if( origList( counter + iVrtx ) .gt. 0 )then
            pos = tem_posOfID( origList( counter + iVrtx ),                    &
              &                unique%val( nUnique+1:unique%nVals ))
            pos = pos + nUnique
          ! if not do a simple search in the first part of the unique array
          ! (the negative entries)
          else
            pos = tem_posOfLong( origList( counter + iVrtx ),                  &
              &                  unique%val( 1:nUnique ))
          end if
          ! map the vertex
          call append( array     = map,                                        &
            &          position1 = iElem,                                      &
            &          position2 = iVrtx,                                      &
            &          value     = pos )
        end do
        counter = counter + 20
      else
        do iVrtx = 1, 8
          pos = tem_posOfID( origList( counter+iVrtx ),                        &
            &                unique%val( nUnique+1:unique%nVals ))
          pos = pos + nUnique
          call append( array     = map,                                        &
            &          position1 = iElem,                                      &
            &          position2 = iVrtx,                                      &
            &          value     = pos )
        end do
        counter = counter + 8
      end if
    end do

    deallocate(origList)
    call destroy( me = unique )

  end subroutine tem_unify_vrtx