tem_restart_closeWrite Subroutine

public subroutine tem_restart_closeWrite(me, tree, timing, varSys, subTree)

Complete a number of empty writes (due to higher amount of mpi_file_writes from other processes to finalize the write process), close the restart dump file and write the last header.

Arguments

Type IntentOptional Attributes Name
type(tem_restart_type) :: me

the restart information

type(treelmesh_type) :: tree

mesh, provided in treelm format

type(tem_time_type), intent(in) :: timing

The timing object holding the simulation time

type(tem_varSys_type), intent(in) :: varSys

global variable system defined in solver

type(tem_subTree_type), optional :: subTree

optional subTree of given tree


Source Code

  subroutine tem_restart_closeWrite( me, tree, timing, varSys, subTree )
    ! -------------------------------------------------------------------- !
    !> the restart information
    type(tem_restart_type) :: me
    !> mesh, provided in treelm format
    type(treelmesh_type) :: tree
    !> The timing object holding the simulation time
    type(tem_time_type), intent(in) :: timing
    !> global variable system defined in solver
    type(tem_varSys_type), intent(in) :: varSys
    !> optional subTree of given tree
    type(tem_subTree_type), optional :: subTree
    ! -------------------------------------------------------------------- !
    ! variables to catch possible MPI I/O errors
    integer :: iError
    integer :: iChunk ! chunk counter
    ! temporary buffer array for empty writes
    real(kind=rk) :: chunk(1)
    ! -------------------------------------------------------------------- !
    ! Check, if the number of calls to mpi_file_write_all corresponds to
    ! the maximum number throughout all processes
    ! if local process has performed less writes, we need to complete
    ! as many empty writes as the difference to maxnChunks
    me%nChunkElems = 0
    do iChunk = me%write_file%nChunks+1, me%write_file%maxnChunks
      call tem_restart_writeData( restart = me, chunk = chunk )
    end do

    ! now close the binary file
    call MPI_File_close(me%binaryUnit, iError)
    call check_mpi_error( iError,'File close in tem_restart_closeWrite')
    call MPI_Barrier( me%comm%comm, iError )

    ! Now write out the 'last' header which points to the last successful
    ! restart file.
    call tem_restart_writeHeader( me         = me,      &
      &                           tree       = tree,    &
      &                           subTree    = subTree, &
      &                           timing     = timing,  &
      &                           varSys     = varSys,  &
      &                           lastHeader = .true.   )

  end subroutine tem_restart_closeWrite