Initialize the global restart data type and prepare for the restart output
This routine is called as well for read as for write restart tasks It must be called after tem_load_restart and before the first call to any of tem_restart_* This routine is called only when restart is read from the restart table. If restart is performed from the initial condition table, then different routines (which are basically wrapped around by this one) are invoked
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_restart_type), | intent(inout) | :: | me |
The restart object to initialize. |
||
type(tem_solveHead_type), | intent(in), | optional | :: | solver |
Details about the solver. |
|
type(tem_varMap_type), | intent(in) | :: | varMap |
Description of each variable system. This is ignored, if the data is provided by reading a restart. Contains position of variables to dump in restart file in global variable system for a scheme |
||
type(treelmesh_type), | intent(in) | :: | tree |
Mesh, provided in treelm format |
||
type(tem_subTree_type), | intent(in), | optional | :: | subTree |
optional subTree of the given tree |
|
integer, | intent(in), | optional | :: | nDofs_write |
number of degrees of freedom for each variable of the equation system |
|
integer, | intent(in), | optional | :: | chunkSize |
use predefined chunkSize |
|
integer, | intent(in), | optional | :: | solSpec_unit |
Solver specific unit for restart header |
subroutine tem_init_restart( me, solver, varMap, tree, subTree, &
& nDofs_write, chunkSize, solSpec_unit )
! -------------------------------------------------------------------- !
!> The restart object to initialize.
type(tem_restart_type), intent(inout) :: me
!> Details about the solver.
type(tem_solveHead_type), optional, intent(in) :: solver
!> Description of each variable system.
!! This is ignored, if the data is provided by reading
!! a restart.
! type(tem_varSys_type), intent(in) :: varSys
!> Contains position of variables to dump in restart file in global
!! variable system for a scheme
type(tem_varMap_type), intent(in) :: varMap
!> Mesh, provided in treelm format
type(treelmesh_type), intent(in) :: tree
!> optional subTree of the given tree
type(tem_subTree_type), optional, intent(in) :: subTree
!> number of degrees of freedom for each variable of the equation system
integer, optional, intent(in) :: nDofs_write
!> use predefined chunkSize
integer, optional, intent(in) :: chunkSize
!> Solver specific unit for restart header
integer, optional, intent(in) :: solSpec_unit
! -------------------------------------------------------------------- !
integer :: comm, rank, comm_size, locElems
integer(kind=long_k) :: globElems, elemOff
integer :: read_stat
logical :: nUnitOpened
character(len=320) :: solve_line
! -------------------------------------------------------------------- !
! set information about variables to be dumped in restart format.
! Do this irrespective of read or write restart
me%varMap = varMap
! Get the total number of chunks necessary to write data-set through IO
! buffer to disk.
if ( present( subTree ) ) then
comm = subTree%global%comm
rank = subTree%global%myPart
comm_size = subTree%global%nParts
globElems = subTree%global%nElems
elemOff = subTree%elemOffset
locElems = int(subTree%nElems)
else
comm = tree%global%comm
rank = tree%global%myPart
comm_size = tree%global%nParts
globElems = tree%global%nElems
elemOff = tree%elemOffset
locElems = int(tree%nElems)
end if
! Invoke the routine to allocate various variables
call tem_init_restart_alloc( me = me, &
& comm = comm, &
& rank = rank, &
& comm_size = comm_size, &
& solver = solver, &
& nDofs_write = nDofs_write )
! Loop over all the systems to create MPI types for reading restart
! Invoke the routine which creates MPI Types
call tem_init_restart_create_types( me = me, &
& elemOff = elemOff, &
& locElems = locElems )
call tem_restart_getTotalChunks( restart = me, &
& nElems = locElems, &
& comm = comm, &
& chunkSize = chunkSize )
! In root of this restart type,
! if write restart is active and solSpec_unit is present
! and opened then copy the content in solSpec_unit to internal
! restart scratch unit
me%solSpec_unit = -1
if (me%comm%rank == 0 .and. me%controller%writeRestart) then
if (present(solSpec_unit)) then
if (solSpec_unit>0) then
write(dbgUnit(10),*) 'Writing solver specific info in restart:'
inquire(unit=solSpec_unit, opened=nUnitOpened)
if (nUnitOpened) then
me%solSpec_unit = newunit()
open(unit=me%solSpec_unit, status='scratch')
rewind(solSpec_unit)
do
read(solSpec_unit,'(a)', iostat=read_stat) solve_line
if (read_stat /= 0) EXIT
write(dbgUnit(10),*) trim(solve_line)
write(me%solSpec_unit,'(a)') trim(solve_line)
end do
end if !unitOpened
end if !solSpecUnit>0
end if !present solSpecUnit
end if !root process and active writeRestart
end subroutine tem_init_restart