tem_load_weights Subroutine

public subroutine tem_load_weights(me, weights, success)

Arguments

Type IntentOptional Attributes Name
type(treelmesh_type), intent(in) :: me
real(kind=rk), intent(out) :: weights(me%nElems)
logical, intent(out) :: success

Calls

proc~~tem_load_weights~~CallsGraph proc~tem_load_weights tem_load_weights proc~tem_create_endiansuffix tem_create_EndianSuffix proc~tem_load_weights->proc~tem_create_endiansuffix mpi_type_contiguous mpi_type_contiguous proc~tem_load_weights->mpi_type_contiguous mpi_bcast mpi_bcast proc~tem_load_weights->mpi_bcast mpi_type_size mpi_type_size proc~tem_load_weights->mpi_type_size mpi_type_commit mpi_type_commit proc~tem_load_weights->mpi_type_commit mpi_file_get_size mpi_file_get_size proc~tem_load_weights->mpi_file_get_size mpi_file_read_all mpi_file_read_all proc~tem_load_weights->mpi_file_read_all mpi_file_close mpi_file_close proc~tem_load_weights->mpi_file_close proc~check_mpi_error check_mpi_error proc~tem_load_weights->proc~check_mpi_error mpi_file_open mpi_file_open proc~tem_load_weights->mpi_file_open mpi_file_set_view mpi_file_set_view proc~tem_load_weights->mpi_file_set_view mpi_type_free mpi_type_free proc~tem_load_weights->mpi_type_free proc~tem_levelof tem_LevelOf proc~tem_load_weights->proc~tem_levelof 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~~tem_load_weights~~CalledByGraph proc~tem_load_weights tem_load_weights proc~getelemweight getElemWeight proc~getelemweight->proc~tem_load_weights proc~load_tem load_tem proc~load_tem->proc~tem_load_weights 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_weights( me, weights, success )
    ! -------------------------------------------------------------------- !
    type(treelmesh_type), intent(in) :: me
    real(kind=rk), intent(out) :: weights(me%nElems)
    logical, intent(out) :: success
    ! -------------------------------------------------------------------- !
    integer(kind=MPI_OFFSET_KIND) :: displacement
    integer(kind=MPI_OFFSET_KIND) :: filebytes
    integer :: fh, ftype, iError, typesize
    integer :: iostatus( MPI_STATUS_SIZE )
    integer :: level
    integer :: iElem
    logical :: ex
    character(len=PathLen) :: filename
    character(len=4) :: EndianSuffix
    ! -------------------------------------------------------------------- !

    EndianSuffix = tem_create_endianSuffix()
    if (trim(me%weights_file) /= '') then
      filename = trim(me%weights_file)//EndianSuffix
      ! check if a corresponding weight file exists
      if (me%global%myPart == 0) then
        inquire(file=trim(filename), exist=ex)
      end if

      call MPI_Bcast(ex, 1, MPI_LOGICAL, 0, me%global%comm, iError)
    else
      ex = .false.
    end if

    if (ex) then
      ! Found a weights file, which is used to read a weight for each
      ! element.
      write(logUnit(3),*) 'Loading Weights from file: '//trim(filename)
      ! Open the binary file for MPI I/O (Write)
      call MPI_File_open( me%global%comm, trim(filename),            &
        &                 MPI_MODE_RDONLY, MPI_INFO_NULL, fh, iError )
      call check_mpi_error(iError,'file_open in load_weights')

      ! Create a MPI Subarray  as ftype for file view
      call MPI_Type_contiguous( me%nElems, rk_mpi , ftype, iError )
      call check_mpi_error(iError,'type ftype in load_weights')
      call MPI_Type_commit( ftype, iError )
      call check_mpi_error(iError,'commit ftype in load_weights')

      !get size of etype
      call MPI_Type_size(rk_mpi, typesize, iError )
      call check_mpi_error(iError,'typesize in load_weights')

      call MPI_File_get_size(fh, filebytes, iError)

      ! Check whether the weight file has as many values as there are
      ! elements in the mesh.
      if (filebytes /= typesize*me%global%nElems) then

        write(logunit(3),*) 'Mesh file has wrong number of elements!'
        write(logunit(3),*) 'NO BALANCING!'
        success = .false.

      else

        ! File has correct number of elements matching the mesh.
        ! calculate displacement
        displacement = me%elemoffset * typesize * 1_MPI_OFFSET_KIND

        ! Set the view for each process on the file above
        call MPI_File_set_view( fh, displacement, rk_mpi, ftype, &
          &                     "native", MPI_INFO_NULL, iError  )
        call check_mpi_error(iError,'set_view in load_weights')

        ! Read data from the file
        call MPI_File_read_all( fh, weights, me%nElems, rk_mpi, iostatus, &
          &                     iError                                    )
        call check_mpi_error(iError,'read_all in load_weights')

        success = .true.

      end if

      !Free the MPI_Datatypes which were created and close the file
      call MPI_Type_free(ftype, iError)
      call check_mpi_error(iError,'free ftype in load_weights')
      call MPI_File_close(fh,    iError)
      call check_mpi_error(iError,'close file in load_weights')

    else if (maxval(me%levelweight) > 0.0_rk) then
      write(logunit(3),*) 'Using levelwise weights.'
      do iElem=1,me%nElems
        level = tem_levelOf( me%treeID(iElem) )
        weights(iElem) = me%levelWeight(level)
      end do
      success = .true.
    else
      write(logunit(3),*) 'NO BALANCING!'
      success = .false.
    end if
    write(logUnit(3),*) 'Done loading weights.'

  end subroutine tem_load_weights