load_tem_BC_prop Subroutine

public subroutine load_tem_BC_prop(me, offset, nElems, basename, myPart, comm)

load bc property header from lua file, boundaryID from bnd.lsb

Arguments

Type IntentOptional Attributes Name
type(tem_BC_prop_type), intent(inout) :: me

Boundary condition construct to load the data into

integer(kind=long_k), intent(in) :: offset

Offset of the local set of elements in the global list

integer, intent(in) :: nElems

Local number of elements with this property

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

Name of the file, the data is stored in, will be appended with ".ascii" for the header information and ".lsb" or ".msb" for the binary data.

integer, intent(in) :: myPart

Partition to load

integer, intent(in) :: comm

Communicator to use


Calls

proc~~load_tem_bc_prop~~CallsGraph proc~load_tem_bc_prop load_tem_BC_prop aot_get_val aot_get_val proc~load_tem_bc_prop->aot_get_val aot_table_close aot_table_close proc~load_tem_bc_prop->aot_table_close aot_table_open aot_table_open proc~load_tem_bc_prop->aot_table_open close_config close_config proc~load_tem_bc_prop->close_config mpi_bcast mpi_bcast proc~load_tem_bc_prop->mpi_bcast mpi_comm_rank mpi_comm_rank proc~load_tem_bc_prop->mpi_comm_rank mpi_comm_split mpi_comm_split proc~load_tem_bc_prop->mpi_comm_split mpi_file_close mpi_file_close proc~load_tem_bc_prop->mpi_file_close mpi_file_open mpi_file_open proc~load_tem_bc_prop->mpi_file_open mpi_file_read_all mpi_file_read_all proc~load_tem_bc_prop->mpi_file_read_all mpi_file_set_view mpi_file_set_view proc~load_tem_bc_prop->mpi_file_set_view mpi_type_commit mpi_type_commit proc~load_tem_bc_prop->mpi_type_commit mpi_type_contiguous mpi_type_contiguous proc~load_tem_bc_prop->mpi_type_contiguous mpi_type_free mpi_type_free proc~load_tem_bc_prop->mpi_type_free mpi_type_size mpi_type_size proc~load_tem_bc_prop->mpi_type_size open_config_file open_config_file proc~load_tem_bc_prop->open_config_file proc~check_mpi_error check_mpi_error proc~load_tem_bc_prop->proc~check_mpi_error proc~tem_create_endiansuffix tem_create_EndianSuffix proc~load_tem_bc_prop->proc~tem_create_endiansuffix proc~tem_open tem_open proc~load_tem_bc_prop->proc~tem_open mpi_error_string mpi_error_string proc~check_mpi_error->mpi_error_string proc~tem_abort tem_abort proc~check_mpi_error->proc~tem_abort proc~newunit newunit proc~tem_open->proc~newunit proc~tem_open->proc~tem_abort proc~upper_to_lower upper_to_lower proc~tem_open->proc~upper_to_lower mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~load_tem_bc_prop~~CalledByGraph proc~load_tem_bc_prop load_tem_BC_prop proc~init_tem_bc_prop init_tem_bc_prop proc~init_tem_bc_prop->proc~load_tem_bc_prop

Source Code

  subroutine load_tem_BC_prop( me, offset, nElems, basename, myPart, comm )
    ! ---------------------------------------------------------------------------
    !> Boundary condition construct to load the data into
    type(tem_BC_prop_type), intent(inout) :: me
    !> Offset of the local set of elements in the global list
    integer(kind=long_k), intent(in) :: offset
    !> Local number of elements with this property
    integer, intent(in) :: nElems
    !> Name of the file, the data is stored in, will be appended with
    !! ".ascii" for the header information and ".lsb" or ".msb" for the
    !! binary data.
    character(len=*), intent(in) :: basename
    !> Partition to load
    integer, intent(in) :: myPart
    !> Communicator to use
    integer, intent(in) :: comm
    ! ---------------------------------------------------------------------------
    type( flu_State ) :: conf ! lua flu state to read lua file
    integer :: i
    integer, parameter :: root = 0
    integer :: BCcomm
    integer :: color, iError
    logical :: participant !< If the local rank is a participant in BC
    character(len=4) :: EndianSuffix
    character(len=256) :: headerfile
    character(len=256) :: datafile
    integer :: thandle, typesize
    integer(kind=long_k), allocatable :: buffer(:)
    integer(kind=long_k), allocatable :: globbuffer(:)
    integer(kind=MPI_OFFSET_KIND)     :: displacement
    integer :: fh, etype, ftype, iostatus( MPI_STATUS_SIZE )
    integer :: file_rec_len
    integer :: myBCrank

    ! ---------------------------------------------------------------------------

    headerfile = trim(basename)//'.lua'
    EndianSuffix = tem_create_EndianSuffix()
    datafile = trim(basename)//trim(EndianSuffix)

    if (me%header%nElems > 0) then
      write(logUnit(1), *) 'Load boundary ID from file: '//trim(datafile)
    end if

    if (myPart == root) then
      ! Read the header only on the root process, broadcast to all others
      ! open mesh header file
      call open_config_file( L = conf, filename = headerfile )
      call aot_get_val( L       = conf,      &
        &               key     = 'nSides',  &
        &               val     = me%nSides, &
        &               ErrCode = iError     )
      call aot_get_val( L       = conf,        &
        &               key     = 'nBCtypes',  &
        &               val     = me%nBCtypes, &
        &               ErrCode = iError       )
    end if

    call MPI_Bcast(me%nSides, 1, MPI_INTEGER, root, comm, iError)
    call MPI_Bcast(me%nBCtypes, 1, MPI_INTEGER, root, comm, iError)

    allocate(me%BC_label(me%nBCtypes))
    allocate(me%hasQVal(me%nBCtypes))
    me%hasQVal(:) = .false.

    if (myPart == root) then
      ! Now read the list of boundary labels
      call aot_table_open( L = conf, thandle = thandle, key = 'bclabel' )
      do i=1,me%nBCtypes
        call aot_get_val( L       = conf,           &
          &               thandle = thandle,        &
          &               pos     = i,              &
          &               val     = me%BC_label(i), &
          &               ErrCode = iError          )
      end do
      call aot_table_close( L = conf, thandle = thandle )
      call close_config( conf )
    end if

    call MPI_Bcast( me%BC_label, LabelLen*me%nBCtypes, MPI_CHARACTER, &
      &             root, comm, iError                                )

    allocate(me%boundary_ID(me%nSides, nElems))

    participant = ( nElems > 0 )

    if (participant) then
      color = 1
    else
      color = MPI_UNDEFINED
    end if

    ! Split the communicator
    call MPI_COMM_SPLIT(comm, color, myPart, BCcomm, iError)

    if (nElems > 0) then

      allocate( buffer(me%nSides * nElems) )

      if (me%header%nElems*me%nSides > io_buffer_size) then
        write(logUnit(5), *) 'read with MPI'


        ! Create a contiguous type to describe the vector per element
        call MPI_TYPE_CONTIGUOUS( me%nSides, long_k_mpi, etype, iError )
        call check_mpi_error(iError,'contiguous etype in load_tem_BC_prop')
        call MPI_TYPE_COMMIT(     etype, iError )
        call check_mpi_error(iError,'commit etype in load_tem_BC_prop')
        call MPI_TYPE_SIZE(etype, typesize, iError )
        call check_mpi_error(iError,'typesize in load_tem_BC_prop')

        ! Calculate displacement for file view
        displacement = offset * typesize * 1_MPI_OFFSET_KIND

        ! Create a MPI CONTIGUOUS as ftype for file view
        call MPI_TYPE_CONTIGUOUS(nElems, etype, ftype, iError)
        call check_mpi_error(iError,'contiguous ftype in load_tem_BC_prop')
        call MPI_TYPE_COMMIT( ftype, iError )
        call check_mpi_error( iError, 'commit ftype in load_tem_BC_prop')

        ! Open the binary file for MPI I/O (read)
        call MPI_FILE_OPEN( BCcomm, trim(datafile), MPI_MODE_RDONLY,   &
          &                 MPI_INFO_NULL, fh, iError                  )
        call check_mpi_error( iError, 'Open File in load_tem_BC_prop')

        ! Set the view for each process on the file above
        call MPI_FILE_SET_VIEW( fh, displacement, etype, ftype, "native",  &
          &                     MPI_INFO_NULL, iError                      )
        call check_mpi_error( iError,'Set File view in load_tem_BC_prop')

        ! Read data from the file
        call MPI_FILE_READ_ALL( fh, buffer, nElems, etype, iostatus, iError )
        call check_mpi_error( iError,'Read All in load_tem_BC_prop')

        !Free the MPI_Datatypes which were created and close the file
        call MPI_TYPE_FREE (etype, iError)
        call check_mpi_error( iError,'free etype in load_tem_BC_prop')
        call MPI_TYPE_FREE (ftype, iError)
        call check_mpi_error( iError,'free ftype in load_tem_BC_prop')
        call MPI_FILE_CLOSE(fh,    iError)
        call check_mpi_error( iError,'close file in load_tem_BC_prop')
        ! END IO-part

      else

        ! File is so small, it probably is faster to read it on a single process
        ! and broadcast the data.
        call MPI_Comm_rank(BCcomm, myBCrank, iError)


        allocate(globbuffer(me%header%nElems*me%nSides))

        if (myBCrank == 0) then
          write(logUnit(5), *) 'read on a single process'
          inquire(iolength = file_rec_len) globbuffer
          call tem_open( newunit = fh,             &
            &            file    = trim(datafile), &
            &            recl    = file_rec_len,   &
            &            action  = 'read',         &
            &            access  = 'direct',       &
            &            form    = 'unformatted'   )

          read(fh, rec=1) globbuffer

          close(fh)
         end if

        call MPI_Bcast( globbuffer, int(me%nSides*me%header%nElems), &
          &             long_k_mpi, 0, BCcomm, iError                )

        buffer = globbuffer(offset*me%nSides+1:(offset+nElems)*me%nSides)
        deallocate(globbuffer)

      end if

      do i=1,nElems
        me%boundary_ID(:,i) = buffer( ((i-1)*me%nSides+1) : (i*me%nSides) )
      end do

      deallocate(buffer)

   end if

  end subroutine load_tem_BC_prop