dump_tem_BC_realArray Subroutine

private subroutine dump_tem_BC_realArray(offset, arraylen, nElems, propdat, basename, comm)

dump bc properties

Mpi IO

Arguments

Type IntentOptional Attributes Name
integer(kind=long_k), intent(in) :: offset

Offset of the local set of elements in the global list

integer, intent(in) :: arraylen

Length of real array for each element

integer, intent(in) :: nElems

Local number of elements with this property

real(kind=rk), intent(in) :: propdat(arraylen,nElems)

Real-valued property data for each element to write out

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

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

integer, intent(in) :: comm

Communicator to use


Calls

proc~~dump_tem_bc_realarray~~CallsGraph proc~dump_tem_bc_realarray dump_tem_BC_realArray mpi_file_close mpi_file_close proc~dump_tem_bc_realarray->mpi_file_close mpi_file_open mpi_file_open proc~dump_tem_bc_realarray->mpi_file_open mpi_file_set_view mpi_file_set_view proc~dump_tem_bc_realarray->mpi_file_set_view mpi_file_write_all mpi_file_write_all proc~dump_tem_bc_realarray->mpi_file_write_all mpi_type_commit mpi_type_commit proc~dump_tem_bc_realarray->mpi_type_commit mpi_type_contiguous mpi_type_contiguous proc~dump_tem_bc_realarray->mpi_type_contiguous mpi_type_free mpi_type_free proc~dump_tem_bc_realarray->mpi_type_free mpi_type_size mpi_type_size proc~dump_tem_bc_realarray->mpi_type_size proc~check_mpi_error check_mpi_error proc~dump_tem_bc_realarray->proc~check_mpi_error proc~tem_create_endiansuffix tem_create_EndianSuffix proc~dump_tem_bc_realarray->proc~tem_create_endiansuffix 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 mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~dump_tem_bc_realarray~~CalledByGraph proc~dump_tem_bc_realarray dump_tem_BC_realArray proc~dump_tem_bc_normal dump_tem_BC_normal proc~dump_tem_bc_normal->proc~dump_tem_bc_realarray proc~dump_tem_bc_qval dump_tem_BC_qVal proc~dump_tem_bc_qval->proc~dump_tem_bc_realarray

Source Code

  subroutine dump_tem_BC_realArray( offset, arraylen, nElems, propdat, basename, &
    &                               comm                                         )
    ! ---------------------------------------------------------------------------
    !> Offset of the local set of elements in the global list
    integer(kind=long_k), intent(in) :: offset
    !> Length of real array for each element
    integer, intent(in) :: arraylen
    !> Local number of elements with this property
    integer, intent(in) :: nElems
    !> Real-valued property data for each element to write out
    real(kind=rk), intent(in) :: propdat(arraylen, nElems)
    !> Name of the file, the data is stored in, will be appended with
    !! ".lua" for the header information and ".lsb" or ".msb" for the
    !! binary data.
    character(len=*), intent(in) :: basename
    !> Communicator to use
    integer, intent(in) :: comm
    ! ---------------------------------------------------------------------------
    integer, parameter :: root = 0
    character(len=256) :: datafile
    character(len=4) :: EndianSuffix
    ! ---------------------------------------------------------------------------
    integer(kind=MPI_OFFSET_KIND)     :: displacement
    integer :: fh, etype, ftype, ierror, iostatus( MPI_STATUS_SIZE ), typesize
    ! ---------------------------------------------------------------------------

    EndianSuffix = tem_create_EndianSuffix()
    datafile = trim(basename)//trim(EndianSuffix)

    if (nElems > 0) then
      !> Mpi IO
      write(logUnit(1),*) 'Write qVal to file: ' // trim(datafile)

      ! Open the binary file for MPI I/O (Write)
      call MPI_FILE_OPEN( comm, trim(datafile),                  &
        &                 ior(MPI_MODE_WRONLY,MPI_MODE_CREATE),  &
        &                 MPI_INFO_NULL, fh, iError              )
      call check_mpi_error( iError,'File open in dump_tem_BC_realArray')

      ! Create a contiguous type to describe the vector per element
      call MPI_TYPE_CONTIGUOUS( arrayLen, rk_mpi, etype, iError )
      call check_mpi_error( iError,'contiguous etype in dump_tem_BC_realArray')
      call MPI_TYPE_COMMIT(     etype, iError )
      call check_mpi_error( iError,'commit etype in dump_tem_BC_realArray')
      call MPI_TYPE_SIZE(etype, typesize, iError )
      call check_mpi_error(iError,'typesize in dump_tem_BC_realArray')

      ! 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 dump_tem_BC_realArray')
      call MPI_TYPE_COMMIT( ftype, iError )
      call check_mpi_error( iError,'commit ftype in dump_tem_BC_realArray')

      ! 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 dump_tem_BC_realArray')

      ! Read data from the file
      call MPI_FILE_WRITE_ALL( fh, propdat, nElems, etype, iostatus, iError )
      call check_mpi_error( iError,'File write all in dump_tem_BC_realArray')

      !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 dump_tem_BC_realArray')
      call MPI_TYPE_FREE (ftype, iError)
      call check_mpi_error( iError,'free ftype in dump_tem_BC_realArray')
      call MPI_FILE_CLOSE(fh,    iError)
      call check_mpi_error( iError,'close file in dump_tem_BC_realArray')
      ! END IO-part
    end if

  end subroutine dump_tem_BC_realArray