load_tem_global Subroutine

public subroutine load_tem_global(me, dirname, myPart, nParts, comm)

A routine to load global informations from the header file in the given directory.

Read the header only on the root process, broadcast to all others

Broadcast the header informations to all processes.

Arguments

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

Structure to store header in

character(len=*), intent(in) :: dirname

Directory containing the mesh informations

integer, intent(in) :: myPart

The process local part (= 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~~load_tem_global~~CallsGraph proc~load_tem_global load_tem_global proc~open_config_file open_config_file proc~load_tem_global->proc~open_config_file mpi_bcast mpi_bcast proc~load_tem_global->mpi_bcast proc~aot_table_open aot_table_open proc~load_tem_global->proc~aot_table_open proc~load_tem_prophead load_tem_prophead proc~load_tem_global->proc~load_tem_prophead proc~aot_table_close aot_table_close proc~load_tem_global->proc~aot_table_close proc~aot_get_val~2 aot_get_val proc~load_tem_global->proc~aot_get_val~2 proc~close_config close_config proc~load_tem_global->proc~close_config proc~load_tem_prophead->mpi_bcast proc~load_tem_prophead->proc~aot_table_open proc~load_tem_prophead->proc~aot_table_close proc~aot_get_val aot_get_val proc~load_tem_prophead->proc~aot_get_val

Called by

proc~~load_tem_global~~CalledByGraph proc~load_tem_global load_tem_global proc~tem_global_mesh_read tem_global_mesh_read proc~tem_global_mesh_read->proc~load_tem_global proc~load_tem load_tem proc~load_tem->proc~load_tem_global 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 load_tem_global( me, dirname, myPart, nParts, comm )
    ! -------------------------------------------------------------------- !
    !> Structure to store header in
    type(tem_global_type), intent(out) :: me
    !> Directory containing the mesh informations
    character(len=*), intent(in) :: dirname
    !> The process local part (= 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=300) :: headname
    integer :: iError
    integer :: root
    integer :: i
    logical :: ex
    integer :: thandle, sub_handle
    type( flu_State ) :: conf ! lua flu state to read lua file
    ! -------------------------------------------------------------------- !

    root = 0

    me%comm = comm

    me%myPart = myPart
    me%nParts = nParts
    me%dirname = trim(adjustl(dirname))
    headname = trim(me%dirname)//'header.lua'
    write(logUnit(1), *) 'Load mesh header from file: '//trim(headname)

    if (myPart == root) then
      inquire(file=trim(headname), exist=ex)
      if (.not. ex) then
        write(*,*) 'File ',trim(headname),' not found. Aborting.'
        stop
      endif
      !! Read the header only on the root process, broadcast to all others
      ! open mesh header file
      call open_config_file(L = conf, filename = trim(headname))
      ! load label
      call aot_get_val( L       = conf,     &
        &               key     = 'label',  &
        &               val     = me%label, &
        &               ErrCode = iError    )
      call aot_get_val( L       = conf,       &
        &               key     = 'comment',  &
        &               val     = me%comment, &
        &               ErrCode = iError      )

      ! Open boundingbox table
      call aot_table_open( L = conf, thandle = thandle, key='boundingbox' )

      ! Read the origin
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'origin'    )
      do i = 1,3
        call aot_get_val( L       = conf,         &
          &               thandle = sub_handle,   &
          &               pos     = i,            &
          &               val     = me%origin(i), &
          &               ErrCode = iError        )
      end do
      call aot_table_close( L = conf, thandle = sub_handle )

      ! Read the bounding cube length
      call aot_get_val( L       = conf,                  &
        &               thandle = thandle,               &
        &               key     = 'length',              &
        &               val     = me%BoundingCubeLength, &
        &               ErrCode = iError                 )
      ! Close boundingbox table again
      call aot_table_close( L = conf, thandle = thandle )

      call aot_get_val( L       = conf,      &
        &               key     = 'nElems',  &
        &               val     = me%nElems, &
        &               ErrCode = iError     )

      call aot_get_val( L       = conf,        &
        &               key     = 'minLevel',  &
        &               val     = me%minLevel, &
        &               ErrCode = iError       )

      call aot_get_val( L       = conf,        &
        &               key     = 'maxLevel',  &
        &               val     = me%maxLevel, &
        &               ErrCode = iError       )

      call aot_get_val( L       = conf,           &
        &               key     = 'nProperties',  &
        &               val     = me%nProperties, &
        &               ErrCode = iError          )

      ! Read the effective bounding cube parameters
      ! Open the effective bounding box table
      call aot_table_open( L = conf, thandle = thandle, key='effBoundingbox' )
      ! Read the origin
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'origin'    )
      do i = 1,3
        call aot_get_val( L       = conf,                    &
          &               thandle = sub_handle,              &
          &               pos     = i,                       &
          &               val     = me%effboundingcube(i,1), &
          &               ErrCode = iError                   )
      end do
      call aot_table_close( L = conf, thandle = sub_handle )
      me%effOrigin = me%effboundingcube(:,1)
      ! Read the effective length (min and max)
      call aot_table_open( L       = conf,       &
        &                  parent  = thandle,    &
        &                  thandle = sub_handle, &
        &                  key     = 'effLength' )
      do i = 1,3
        call aot_get_val( L       = conf,            &
          &               thandle = sub_handle,      &
          &               pos     = i,               &
          &               val     = me%effLength(i), &
          &               ErrCode = iError           )
        me%effboundingcube(i,2) = me%effboundingcube(i,1) + me%effLength(i)
      end do
      call aot_table_close( L = conf, thandle = sub_handle )

    end if
    write(logUnit(1),*) 'The real bounding cube is...'
    write(logUnit(1),*) '  min: ',me%effBoundingCube(:,1)
    write(logUnit(1),*) '  max: ',me%effBoundingCube(:,2)

    !! Broadcast the header informations to all processes.
    call MPI_Bcast(me%nElems, 1, long_k_mpi, root, me%comm, iError)
    call MPI_Bcast(me%label, LabelLen, MPI_CHARACTER, root, me%comm, iError)
    call MPI_Bcast(me%comment, LabelLen, MPI_CHARACTER, root, me%comm, iError)
    call MPI_Bcast(me%BoundingCubeLength, 1, rk_mpi, root, me%comm, iError)
    call MPI_Bcast(me%Origin, 3, rk_mpi, root, me%comm, iError)
    call MPI_Bcast(me%minLevel, 1, MPI_INTEGER, root, me%comm, iError)
    call MPI_Bcast(me%maxLevel, 1, MPI_INTEGER, root, me%comm, iError)
    call MPI_Bcast(me%nProperties, 1, MPI_INTEGER, root, me%comm, iError)
    call MPI_Bcast(me%effBoundingCube, 6, rk_mpi, root, me%comm, iError)

    if (associated(me%Property)) deallocate(me%property)
    allocate(me%Property(me%nProperties))

    call load_tem_prophead( me     = me%Property, &
      &                     myPart = myPart,      &
      &                     comm   = me%comm,     &
      &                     conf   = conf,        &
      &                     root   = root         )

    if (myPart == root) then
      call close_config(conf)
    end if

  end subroutine load_tem_global