tem_alloc_levelDesc Subroutine

private subroutine tem_alloc_levelDesc(me, minLevel, maxLevel, initlen, nStencils)

Allocate level descriptor and initilize its variables

Arguments

Type IntentOptional Attributes Name
type(tem_levelDesc_type), intent(out), allocatable :: me(:)
integer, intent(in) :: minLevel
integer, intent(in) :: maxLevel
integer, intent(in) :: initlen
integer, intent(in) :: nStencils

Calls

proc~~tem_alloc_leveldesc~~CallsGraph proc~tem_alloc_leveldesc tem_alloc_levelDesc interface~init~22 init proc~tem_alloc_leveldesc->interface~init~22 proc~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real

Called by

proc~~tem_alloc_leveldesc~~CalledByGraph proc~tem_alloc_leveldesc tem_alloc_levelDesc proc~tem_init_elemlevels tem_init_elemLevels proc~tem_init_elemlevels->proc~tem_alloc_leveldesc proc~tem_create_leveldesc tem_create_levelDesc proc~tem_create_leveldesc->proc~tem_init_elemlevels proc~tem_dimbydim_construction tem_dimByDim_construction proc~tem_dimbydim_construction->proc~tem_create_leveldesc proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_dimbydim_construction

Contents

Source Code


Source Code

  subroutine tem_alloc_levelDesc( me, minLevel, maxLevel, initlen, nStencils )
    ! ---------------------------------------------------------------------------
    type(tem_levelDesc_type), allocatable, intent(out)  :: me(:)
    integer, intent(in) :: minLevel, maxLevel, initlen, nStencils
    ! ---------------------------------------------------------------------------
    integer :: iLevel
    ! ---------------------------------------------------------------------------

    write(logUnit(5),*) 'Allocating level Descriptor ...'
    allocate( me(minLevel:maxLevel) )
    me(minLevel:maxLevel)%nElems = 0

    ! Initial length for the element array
    ! We don't know how many elements there will be per level, just using
    ! the average here as a starting point. This should at least for single
    ! level meshes reduce the number of expands, we need to do.
    ! Giving a small space for additional elements (ghosts and halos, might
    ! avoid immidiate doubling, for elements beyond the fluids).

    ! initialize dynamic require array for requireNeighNeigh
    ! and element type
    write(logUnit(7),*) 'Initializing elemList, require and neigh array'
    do iLevel = minLevel, maxLevel
      call init( me = me( iLevel )%require )
      call init( me = me( iLevel )%elem,    length = initlen )
      ! For each of the neighbor lists create the horizontal (neighbor)
      ! relations of the size of the total number of stencils used in this
      ! scheme
      allocate( me( iLevel )%neigh( nStencils ) )
    end do

  end subroutine tem_alloc_levelDesc