tem_load_internal Subroutine

public subroutine tem_load_internal(me, conf, thandle, myPart, nParts, comm)

Load an internally generated mesh

The predefined meshes are illustrated in the following configuration snippet:

 -- use the predefined full cube
 mesh = { predefined = 'cube',
          -- for a single element in Z-direction use slice:
          -- predefined = 'slice',

          -- for just a line (single element in Y and Z) use line:
          -- predefined = 'line',
          -- If you set the refinementLevel to 0, a single element will
          -- be created with treeID=1 on level 1, with the bounding cube
          -- redefined such, that the element fills the original
          -- specification for the bounding cube.
          -- The line can also be generated with an arbitrary number of
          -- elements instead of a level. In this case the length refers
          -- to the overall length of the given elements not the extent of
          -- the cube. Periodicity is assumed in X direction as well.
          -- To create the line with a certain number of elements set the
          -- element_count instead of the refinementLevel. If both are
          -- given, the element_count overrules the refinementLevel.

          -- for a line with boundary conditions line_bounded:
          -- predefined = 'line_bounded'
          -- Like line, but needs to have `element_count` given, not
          -- `refinementLevel`. This mesh defines boundary conditions in the
          -- X direction as west (on the left) and east (on the right), that
          -- then need to be defined in the solver configuration.
 -- origin of the cube
 origin = {0.,0.,0.},
 -- length of the cube
 length = 10.,
 -- refinement level to resolve the cube
 refinementLevel = 4 }

Arguments

Type IntentOptional Attributes Name
type(treelmesh_type), intent(out) :: me

Structure to load the mesh to

type(flu_State) :: conf

Directory containing the mesh informations

integer, intent(in) :: thandle

Handle for the table to read the description of the mesh from.

integer, intent(in) :: myPart

Partition to use on the calling process (= MPI Rank in comm)

integer, intent(in) :: nParts

Number of partitions, the mesh is partitioned into (= Number of MPI processes in comm).

integer, intent(in) :: comm

MPI Communicator to use


Calls

proc~~tem_load_internal~~CallsGraph proc~tem_load_internal tem_load_internal proc~generate_treelm_cube generate_treelm_cube proc~tem_load_internal->proc~generate_treelm_cube proc~generate_treelm_slice generate_treelm_slice proc~tem_load_internal->proc~generate_treelm_slice proc~aot_table_open aot_table_open proc~tem_load_internal->proc~aot_table_open proc~aot_table_close aot_table_close proc~tem_load_internal->proc~aot_table_close proc~aot_get_val~2 aot_get_val proc~tem_load_internal->proc~aot_get_val~2 proc~aot_table_length aot_table_length proc~tem_load_internal->proc~aot_table_length proc~generate_treelm_single generate_treelm_single proc~tem_load_internal->proc~generate_treelm_single proc~tem_abort tem_abort proc~tem_load_internal->proc~tem_abort proc~generate_treelm_elements generate_treelm_elements proc~tem_load_internal->proc~generate_treelm_elements proc~generate_treelm_line generate_treelm_line proc~tem_load_internal->proc~generate_treelm_line proc~tem_firstidatlevel tem_FirstIdAtLevel proc~generate_treelm_cube->proc~tem_firstidatlevel proc~tem_lastidatlevel tem_LastIdAtLevel proc~generate_treelm_cube->proc~tem_lastidatlevel proc~generate_treelm_slice->proc~tem_firstidatlevel proc~tem_coordof_2d_id tem_CoordOf_2d_Id proc~generate_treelm_slice->proc~tem_coordof_2d_id proc~tem_idofcoord tem_IdOfCoord proc~generate_treelm_slice->proc~tem_idofcoord mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~generate_treelm_elements->proc~tem_firstidatlevel proc~generate_treelm_elements->proc~tem_idofcoord proc~generate_treelm_line->proc~tem_firstidatlevel proc~generate_treelm_line->proc~tem_idofcoord

Called by

proc~~tem_load_internal~~CalledByGraph proc~tem_load_internal tem_load_internal proc~load_tem load_tem proc~load_tem->proc~tem_load_internal proc~tem_restart_readheader tem_restart_readHeader proc~tem_restart_readheader->proc~load_tem proc~tem_load_restart tem_load_restart proc~tem_load_restart->proc~tem_restart_readheader

Contents

Source Code


Source Code

  subroutine tem_load_internal( me, conf, thandle, myPart, nParts, comm )
    ! -------------------------------------------------------------------- !
    !> Structure to load the mesh to
    type(treelmesh_type), intent(out) :: me
    !> Directory containing the mesh informations
    type(flu_State) :: conf
    !> Handle for the table to read the description
    !! of the mesh from.
    integer, intent(in) :: thandle
    !> Partition to use on the calling process (= MPI Rank in comm)
    integer, intent(in) :: myPart
    !> Number of partitions, the mesh is partitioned into (= Number of MPI
    !! processes in comm).
    integer, intent(in) :: nParts
    !> MPI Communicator to use
    integer, intent(in) :: comm
    ! -------------------------------------------------------------------- !
    character(len=40) :: meshtype
    real(kind=rk) :: origin(3)
    real(kind=rk) :: length
    integer :: elementcount
    integer :: level
    integer :: sub_handle
    integer :: iError, i
    integer :: orig_err(3)
    ! -------------------------------------------------------------------- !

    call aot_get_val( L       = conf,         &
      &               thandle = thandle,      &
      &               val     = meshtype,     &
      &               ErrCode = iError,       &
      &               key     = 'predefined', &
      &               default = 'cube'        )

    select case(trim(meshtype))
    case('cube')
      write(logUnit(1),*)'Creating a full cubical mesh without boundaries'

      ! Get the origin of the cube:
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'origin'    )
      if (aot_table_length(L = conf, thandle = sub_handle) == 3) then
        do i=1,3
          call aot_get_val( L       = conf,       &
            &               thandle = sub_handle, &
            &               val     = origin(i),  &
            &               ErrCode = iError,     &
            &               pos     = i           )
        end do
      else
        write(logUnit(1),*)'Origin specified for the mesh to generate has to ' &
          &            //'have 3 entries!'
        write(logUnit(1),*)'Default to {0., 0., 0.} '
        origin = 0.0_rk
      end if
      call aot_table_close( L=conf, thandle=sub_handle )

      ! Get the length of the cube:
      call aot_get_val( L       = conf,     &
        &               thandle = thandle,  &
        &               val     = length,   &
        &               ErrCode = iError,   &
        &               key     = 'length', &
        &               default = 1.0_rk    )

      ! Get the refinement level:
      call aot_get_val( L       = conf,             &
        &               thandle = thandle,          &
        &               val     = level,            &
        &               ErrCode = iError,           &
        &               key     = 'refinementLevel' )

      ! Now generate a mesh with the gathered configuration parameters
      if (level > 0) then
        ! Only if the level is greater than 0, actually create mesh.
        call generate_treelm_cube( me     = me,     &
          &                        origin = origin, &
          &                        length = length, &
          &                        level  = level,  &
          &                        myPart = myPart, &
          &                        nParts = nParts, &
          &                        comm   = comm    )
      else
        ! Meshes on level 0 have only one element, but we need treeids > 0,
        ! thus we create single elements on level 1 with treeid=1, and
        ! redefine the bounding cube accordingly.
        call generate_treelm_single( me     = me,       &
          &                          origin = origin,   &
          &                          length = 2*length, &
          &                          myPart = myPart,   &
          &                          nParts = nParts,   &
          &                          comm   = comm      )
      end if


    case('slice')
      write(logUnit(1),*)'Creating a slice (single element in Z) mesh' &
        &          //' without boundaries'

      ! Get the origin of the cube:
      origin = 0.0_rk
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'origin'    )
      if (aot_table_length(L = conf, thandle = sub_handle) >= 2) then
        do i=1,2
          call aot_get_val( L       = conf,       &
            &               thandle = sub_handle, &
            &               val     = origin(i),  &
            &               ErrCode = iError,     &
            &               pos     = i           )
        end do
      else
        write(logUnit(1),*)'Origin specified for the mesh to generate has to ' &
          &            // 'have 2 entries!'
        write(logUnit(1),*)'Default to {0., 0.} '
      end if
      call aot_get_val( L       = conf,       &
        &               thandle = sub_handle, &
        &               val     = origin(3),  &
        &               ErrCode = iError,     &
        &               Default = 0.0_rk,     &
        &               pos     = 3           )
      call aot_table_close( L=conf, thandle=sub_handle )

      ! Get the length of the cube:
      call aot_get_val( L       = conf,     &
        &               thandle = thandle,  &
        &               val     = length,   &
        &               ErrCode = iError,   &
        &               key     = 'length', &
        &               default = 1.0_rk    )

      ! Get the refinement level:
      call aot_get_val( L       = conf,             &
        &               thandle = thandle,          &
        &               val     = level,            &
        &               ErrCode = iError,           &
        &               key     = 'refinementLevel' )

      ! Now generate a mesh with the gathered configuration parameters
      if (level > 0) then
        ! Only if the level is greater than 0, actually create mesh.
        call generate_treelm_slice( me     = me,     &
          &                         origin = origin, &
          &                         length = length, &
          &                         level  = level,  &
          &                         myPart = myPart, &
          &                         nParts = nParts, &
          &                         comm   = comm    )
      else
        ! Meshes on level 0 have only one element, but we need treeids > 0,
        ! thus we create single elements on level 1 with treeid=1, and
        ! redefine the bounding cube accordingly.
        call generate_treelm_single( me     = me,       &
          &                          origin = origin,   &
          &                          length = 2*length, &
          &                          myPart = myPart,   &
          &                          nParts = nParts,   &
          &                          comm   = comm      )
      end if


    case('line')
      write(logUnit(1),*) 'Creating a line (single element in Y and Z) mesh' &
        &                 // ' without boundaries'

      ! Get the origin of the cube:
      origin = 0.0_rk
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'origin'    )
      if (aot_table_length(L = conf, thandle = sub_handle) >= 1) then
        call aot_get_val( L       = conf,       &
          &               thandle = sub_handle, &
          &               val     = origin(1),  &
          &               ErrCode = iError,     &
          &               pos     = 1           )
      else
        write(logUnit(1),*) 'Origin specified for the mesh to generate has' &
          &                 // ' to have 1 entry!'
        write(logUnit(1),*)'Default to {0.} '
      end if
      call aot_get_val( L       = conf,       &
        &               thandle = sub_handle, &
        &               val     = origin(2),  &
        &               ErrCode = iError,     &
        &               Default = 0.0_rk,     &
        &               pos     = 2           )
      call aot_get_val( L       = conf,       &
        &               thandle = sub_handle, &
        &               val     = origin(3),  &
        &               ErrCode = iError,     &
        &               Default = 0.0_rk,     &
        &               pos     = 3           )
      call aot_table_close( L=conf, thandle=sub_handle )

      ! Get the length of the cube:
      call aot_get_val( L       = conf,     &
        &               thandle = thandle,  &
        &               val     = length,   &
        &               ErrCode = iError,   &
        &               key     = 'length', &
        &               default = 1.0_rk    )

      ! Get the refinement level:
      call aot_get_val( L       = conf,              &
        &               thandle = thandle,           &
        &               val     = level,             &
        &               ErrCode = iError,            &
        &               key     = 'refinementLevel', &
        &               default = -1                 )

      ! Get the element count:
      call aot_get_val( L       = conf,            &
        &               thandle = thandle,         &
        &               val     = elementcount,    &
        &               ErrCode = iError,          &
        &               key     = 'element_count', &
        &               default = -1               )

      if ( (elementcount < 1) .and. (level < 0) ) then
        write(logUnit(1),*) 'For the line mesh you need to specify either'
        write(logUnit(1),*) 'the element_count or the refinementLevel!'
        write(logUnit(1),*) 'STOPPING...'
        call tem_abort()
      end if

      ! Now generate a mesh with the gathered configuration parameters
      if (elementcount > 1) then
        call generate_treelm_elements( me           = me,              &
          &                            origin       = origin,          &
          &                            length       = length,          &
          &                            elementcount = elementcount,    &
          &                            myPart       = myPart,          &
          &                            nParts       = nParts,          &
          &                            comm         = comm,            &
          &                            predefined   = trim(meshtype),  &
          &                            bclabel      = 'internal 1D BC' )
      else
        if ((level > 0) .and. (elementcount < 1)) then
          ! Only if the level is greater than 0, actually create mesh.
          call generate_treelm_line( me     = me,     &
            &                        origin = origin, &
            &                        length = length, &
            &                        level  = level,  &
            &                        myPart = myPart, &
            &                        nParts = nParts, &
            &                        comm   = comm    )
        else
          ! Meshes on level 0 have only one element, but we need treeids > 0,
          ! thus we create single elements on level 1 with treeid=1, and
          ! redefine the bounding cube accordingly.
          call generate_treelm_single( me     = me,       &
            &                          origin = origin,   &
            &                          length = 2*length, &
            &                          myPart = myPart,   &
            &                          nParts = nParts,   &
            &                          comm   = comm      )
        end if
      end if

    case('line_bounded')
      write(logUnit(1),*) 'Creating a line (single element in Y and Z) mesh' &
        &                 // ' WITH boundaries'
      write(logUnit(1),*) 'Boundaries west and east have to be provided!'

      ! Get the origin of the cube:
      origin = 0.0_rk
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'origin'    )
      if (aot_table_length(L = conf, thandle = sub_handle) >= 1) then
        call aot_get_val( L       = conf,       &
          &               thandle = sub_handle, &
          &               val     = origin(1),  &
          &               ErrCode = iError,     &
          &               pos     = 1           )
      else
        write(logUnit(1),*) 'Origin specified for the mesh to generate has' &
          &                 // ' to have 1 entry!'
        write(logUnit(1),*)'Default to {0.} '
      end if
      call aot_get_val( L       = conf,       &
        &               thandle = sub_handle, &
        &               val     = origin(2),  &
        &               ErrCode = iError,     &
        &               Default = 0.0_rk,     &
        &               pos     = 2           )
      call aot_get_val( L       = conf,       &
        &               thandle = sub_handle, &
        &               val     = origin(3),  &
        &               ErrCode = iError,     &
        &               Default = 0.0_rk,     &
        &               pos     = 3           )
      call aot_table_close( L=conf, thandle=sub_handle )

      ! Get the length of the cube:
      call aot_get_val( L       = conf,     &
        &               thandle = thandle,  &
        &               val     = length,   &
        &               ErrCode = iError,   &
        &               key     = 'length', &
        &               default = 1.0_rk    )

      ! Get the element count:
      call aot_get_val( L       = conf,            &
        &               thandle = thandle,         &
        &               val     = elementcount,    &
        &               ErrCode = iError,          &
        &               key     = 'element_count', &
        &               default = -1               )

      if ( (elementcount < 1) ) then
        write(logUnit(1),*) 'For the bounded line mesh you need to specify'
        write(logUnit(1),*) 'the element_count!'
        write(logUnit(1),*) 'STOPPING...'
        call tem_abort()
      end if

      ! Now generate a mesh with the gathered configuration parameters
      call generate_treelm_elements( me           = me,                      &
        &                            origin       = origin,                  &
        &                            length       = length,                  &
        &                            elementcount = elementcount,            &
        &                            myPart       = myPart,                  &
        &                            nParts       = nParts,                  &
        &                            comm         = comm,                    &
        &                            predefined   = trim(meshtype),          &
        &                            bclabel      = 'bounded internal 1D BC' )

    case('single')
      write(logUnit(1),*) 'Creating a single element mesh without boundaries'

      ! Get the origin of the cube:
      call aot_get_val( L       = conf,                    &
        &               thandle = thandle,                 &
        &               val     = origin,                  &
        &               ErrCode = orig_err,                &
        &               key     = 'origin',                &
        &               default = [0.0_rk, 0.0_rk, 0.0_rk] )

      ! Get the length of the cube:
      call aot_get_val( L       = conf,     &
        &               thandle = thandle,  &
        &               val     = length,   &
        &               ErrCode = iError,   &
        &               key     = 'length', &
        &               default = 1.0_rk    )

      ! Now generate a mesh with the gathered configuration parameters
      call generate_treelm_single( me     = me,     &
        &                          origin = origin, &
        &                          length = length, &
        &                          myPart = myPart, &
        &                          nParts = nParts, &
        &                          comm   = comm    )


    case('default')
      write(logUnit(1),*) 'Do not know how to generate mesh ' &
        &                 // trim(meshtype)

    end select

  end subroutine tem_load_internal