This commit implements a better way (less I/O operations) restoring of the parameters during the job restart. Some parameters, such as the size of the block, or the equation system, cannot be changed during the restart. Such parameters need to be restored from the restart snapshots. The new way restores these parameters on the MPI master process only, updates their values if they were changed, and distributes them to other MPI processes. This way the number of I/O operations is kept to the minimum (only one process access one file). Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
676 lines
21 KiB
Fortran
676 lines
21 KiB
Fortran
!!******************************************************************************
|
|
!!
|
|
!! This file is part of the AMUN source code, a program to perform
|
|
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
|
|
!! adaptive mesh.
|
|
!!
|
|
!! Copyright (C) 2017-2023 Grzegorz Kowal <grzegorz@amuncode.org>
|
|
!!
|
|
!! This program is free software: you can redistribute it and/or modify
|
|
!! it under the terms of the GNU General Public License as published by
|
|
!! the Free Software Foundation, either version 3 of the License, or
|
|
!! (at your option) any later version.
|
|
!!
|
|
!! This program is distributed in the hope that it will be useful,
|
|
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
!! GNU General Public License for more details.
|
|
!!
|
|
!! You should have received a copy of the GNU General Public License
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
|
!!
|
|
!!*****************************************************************************
|
|
!!
|
|
!! module: SYSTEM
|
|
!!
|
|
!! This module handles all operations related to the simulated system, i.e.
|
|
!! the initialization of system parameters and related modules (equations,
|
|
!! mesh, evolution, schemes, time statistics), the system execution and
|
|
!! termination.
|
|
!!
|
|
!!*****************************************************************************
|
|
!
|
|
module system
|
|
|
|
implicit none
|
|
|
|
! this flags determines if the system is being resumed
|
|
!
|
|
logical, save :: resumed = .false.
|
|
|
|
#ifdef SIGNALS
|
|
! the flag indicating that the system evolution should be interrupted
|
|
!
|
|
integer, save :: quit = 0
|
|
#endif /* SIGNALS */
|
|
|
|
! the number of resumed run
|
|
!
|
|
integer, save :: nrun = 0
|
|
|
|
! the number of iterations executed during the evolution
|
|
!
|
|
integer, save :: nsteps = 1
|
|
|
|
! the system's permament parameters
|
|
!
|
|
character(len=64) , save :: name = "none"
|
|
character(len=32) , save :: eqsys = "hydrodynamic"
|
|
character(len=32) , save :: eos = "adiabatic"
|
|
character(len=32) , save :: rngtype = "same"
|
|
integer , save :: ncells = 8
|
|
integer , save :: maxlev = 1
|
|
integer , dimension(3) , save :: bdims = 1
|
|
real(kind=8), dimension(3,2), save :: dbnds
|
|
|
|
! the remaining system's parameters
|
|
!
|
|
integer , save :: nghosts = 2
|
|
integer , save :: nmax = huge(1)
|
|
real(kind=8) , save :: tmax = 0.0d+00
|
|
real(kind=8) , save :: trun = 9.9d+99
|
|
real(kind=8) , save :: tsav = 3.0d+01
|
|
|
|
! the minimum workspace size
|
|
!
|
|
integer , save :: nwork = 4194304 ! 32MiB
|
|
|
|
private
|
|
public :: initialize_system, finalize_system, prepare_system, evolve_system
|
|
public :: print_system_info
|
|
public :: resumed, name, nrun, eqsys, eos, ncells, nghosts, maxlev, nsteps
|
|
public :: bdims, dbnds, rngtype, tsav, nwork
|
|
#ifdef SIGNALS
|
|
public :: quit
|
|
#endif /* SIGNALS */
|
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
!
|
|
contains
|
|
!
|
|
!===============================================================================
|
|
!!
|
|
!!*** PUBLIC SUBROUTINES *****************************************************
|
|
!!
|
|
!===============================================================================
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine INITIALIZE_SYSTEM:
|
|
! ----------------------------
|
|
!
|
|
! Subroutine initializes this module and all related modules.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! verbose - the verbose flag;
|
|
! status - the subroutine call status;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine initialize_system(verbose, status)
|
|
|
|
use blocks , only : initialize_blocks
|
|
use boundaries , only : initialize_boundaries
|
|
use coordinates , only : initialize_coordinates, nn => bcells
|
|
use equations , only : initialize_equations, nf, nv
|
|
use evolution , only : initialize_evolution, registers
|
|
use helpers , only : print_message
|
|
use interpolations, only : initialize_interpolations
|
|
use mesh , only : initialize_mesh
|
|
use mpitools , only : check_status
|
|
use problems , only : initialize_problems
|
|
use random , only : reset_seeds
|
|
use statistics , only : initialize_statistics
|
|
use user_problem , only : initialize_user_problem
|
|
|
|
implicit none
|
|
|
|
logical, intent(in) :: verbose
|
|
integer, intent(out) :: status
|
|
|
|
character(len=*), parameter :: loc = 'SYSTEM::initialize_system()'
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
status = 0
|
|
|
|
call restore_parameters(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not restore the system parameters!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call reset_seeds(rngtype)
|
|
|
|
call initialize_equations(eqsys, eos, verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module EQUATIONS!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_interpolations(nghosts, verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module INTERPOLATIONS!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_coordinates(ncells, nghosts, maxlev, bdims, &
|
|
dbnds(1,1), dbnds(1,2), &
|
|
dbnds(2,1), dbnds(2,2), &
|
|
dbnds(3,1), dbnds(3,2), verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module COORDINATES!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_evolution(verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module EVOLUTION!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_blocks(nv, nf, nn, registers, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module BLOCKS!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_problems(name, nrun, verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module PROBLEMS!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_user_problem(name, nrun, verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module USER_PROBLEM!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_boundaries(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module BOUNDARIES!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_mesh(verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module MESH!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
call initialize_statistics(verbose, nrun, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not initialize module STATISTICS!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
nwork = 3 * NDIMS * nf * (ncells + 2 * nghosts)**NDIMS
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine initialize_system
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine FINALIZE_SYSTEM:
|
|
! --------------------------
|
|
!
|
|
! Subroutine frees all the memory used by this module and finalizes related
|
|
! modules.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! verbose - the verbose flag;
|
|
! status - the subroutine call status;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine finalize_system(verbose, status)
|
|
|
|
use blocks , only : finalize_blocks
|
|
use boundaries , only : finalize_boundaries
|
|
use coordinates , only : finalize_coordinates
|
|
use equations , only : finalize_equations
|
|
use evolution , only : finalize_evolution
|
|
use helpers , only : print_message
|
|
use interpolations, only : finalize_interpolations
|
|
use mesh , only : finalize_mesh
|
|
use problems , only : finalize_problems
|
|
use statistics , only : finalize_statistics
|
|
use user_problem , only : finalize_user_problem
|
|
|
|
|
|
implicit none
|
|
|
|
logical, intent(in) :: verbose
|
|
integer, intent(out) :: status
|
|
|
|
character(len=*), parameter :: loc = 'SYSTEM::finalize_system()'
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
status = 0
|
|
|
|
call finalize_statistics(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module STATISTICS!")
|
|
|
|
call finalize_mesh(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module MESH!")
|
|
|
|
call finalize_boundaries(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module BOUNDARIES!")
|
|
|
|
call finalize_user_problem(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module USER_PROBLEM!")
|
|
|
|
call finalize_problems(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module PROBLEMS!")
|
|
|
|
call finalize_blocks(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module BLOCKS!")
|
|
|
|
call finalize_evolution(verbose, status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module EVOLUTION!")
|
|
|
|
call finalize_coordinates(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module COORDINATES!")
|
|
|
|
call finalize_interpolations(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module INTERPOLATIONS!")
|
|
|
|
call finalize_equations(status)
|
|
if (status /=0) &
|
|
call print_message(loc, "Could not finalize module EQUATIONS!")
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine finalize_system
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine PRINT_SYSTEM_INFO:
|
|
! ----------------------------
|
|
!
|
|
! Subroutine prints all system related information.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! verbose - the verbose flag;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine print_system_info(verbose)
|
|
|
|
use boundaries , only : print_boundaries
|
|
use coordinates , only : print_coordinates
|
|
use equations , only : print_equations
|
|
use evolution , only : print_evolution
|
|
use helpers , only : print_section, print_parameter
|
|
use interpolations , only : print_interpolations
|
|
use io , only : print_io
|
|
use mesh , only : print_mesh
|
|
|
|
implicit none
|
|
|
|
logical, intent(in) :: verbose
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
call print_section(verbose, "Problem")
|
|
call print_parameter(verbose, "simulated problem", name)
|
|
call print_parameter(verbose, "resumed", resumed, "yes")
|
|
|
|
call print_equations(verbose)
|
|
call print_evolution(verbose)
|
|
call print_interpolations(verbose)
|
|
call print_coordinates(verbose)
|
|
call print_boundaries(verbose)
|
|
call print_mesh(verbose)
|
|
call print_io(verbose)
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine print_system_info
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine PREPARE_SYSTEM:
|
|
! -------------------------
|
|
!
|
|
! Subroutine prepares the system by generating the mesh and the initial
|
|
! setup or by resuming it from the restart snapshot.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! status - the subroutine call status;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine prepare_system(status)
|
|
|
|
use boundaries, only : boundary_variables
|
|
use evolution , only : initialize_time_step, time
|
|
use helpers , only : print_message
|
|
use io , only : read_restart_snapshot
|
|
use mesh , only : generate_mesh
|
|
|
|
implicit none
|
|
|
|
integer, intent(out) :: status
|
|
|
|
character(len=*), parameter :: loc = 'SYSTEM::prepare_system()'
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
status = 0
|
|
|
|
if (resumed) then
|
|
|
|
call read_restart_snapshot(status)
|
|
if (status /= 0) then
|
|
call print_message(loc, "Could not read the restart snapshot!")
|
|
return
|
|
end if
|
|
|
|
call boundary_variables(time, 0.0d+00, status)
|
|
if (status /= 0) then
|
|
call print_message(loc, "Could not update variable boundaries!")
|
|
return
|
|
end if
|
|
|
|
else
|
|
|
|
call generate_mesh(status)
|
|
if (status /= 0) then
|
|
call print_message(loc, "Could not generate the initial mesh!")
|
|
return
|
|
end if
|
|
|
|
call boundary_variables(0.0d+00, 0.0d+00, status)
|
|
if (status /= 0) then
|
|
call print_message(loc, "Could not update variable boundaries!")
|
|
return
|
|
end if
|
|
|
|
call initialize_time_step()
|
|
|
|
call store_system(status)
|
|
if (status /= 0) &
|
|
call print_message(loc, "Could not store the initial system state!")
|
|
|
|
end if
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine prepare_system
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine EVOLVE_SYSTEM:
|
|
! ------------------------
|
|
!
|
|
! Subroutine performed the time evolution of the system.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! verbose - the verbose flag;
|
|
! status - the subroutine call status;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine evolve_system(verbose, status)
|
|
|
|
use blocks , only : get_nleafs
|
|
use evolution , only : step, time, dt
|
|
use evolution , only : error_control, errtol, nrejections
|
|
use evolution , only : advance, new_time_step
|
|
use helpers , only : print_message
|
|
use io , only : update_dtp, write_restart_snapshot
|
|
use iso_fortran_env, only : output_unit
|
|
use mpitools , only : check_status
|
|
use timers , only : get_timer_total
|
|
|
|
implicit none
|
|
|
|
logical, intent(in) :: verbose
|
|
integer, intent(out) :: status
|
|
|
|
logical :: proceed = .true.
|
|
integer :: ed = 9999, eh = 23, em = 59, es = 59
|
|
real(kind=8) :: thrs = 0.0d+00, tprn = 0.0d+00, ss = 0.0d+00
|
|
|
|
real(kind=8), dimension(8) :: tsim = 0.0d+00, tclc = 0.0d+00
|
|
|
|
character(len=132) :: hfmt, sfmt
|
|
|
|
character(len=*), parameter :: loc = 'SYSTEM::evolve_system()'
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
status = 0
|
|
|
|
hfmt = "('" // achar(10) // "',1x,'Evolving the system:','" // achar(10) //&
|
|
"',5x,'step',4x,'time',10x,'timestep',2x,"
|
|
sfmt = "(i9,2(1x,1es13.6),"
|
|
if (error_control) then
|
|
hfmt = adjustl(trim(hfmt)) // "3x,'err/tol',2x,'nrejs',"
|
|
sfmt = adjustl(trim(sfmt)) // "1x,1es9.2,1x,i6,"
|
|
end if
|
|
hfmt = adjustl(trim(hfmt)) // "4x,'blocks',4x,'completed in')"
|
|
sfmt = adjustl(trim(sfmt)) // &
|
|
"1x,i9,2x,1i4.1,'d',1i2.2,'h',1i2.2,'m',1i2.2,'s'," // &
|
|
"15x,'" // achar(13) // "')"
|
|
|
|
if (verbose) then
|
|
|
|
tsim(:) = time
|
|
tclc(:) = get_timer_total()
|
|
|
|
write(output_unit, hfmt)
|
|
if (error_control) then
|
|
write(output_unit, sfmt, advance="no") step, time, dt, &
|
|
errtol, nrejections, &
|
|
get_nleafs(), ed, eh, em, es
|
|
else
|
|
write(output_unit, sfmt, advance="no") step, time, dt, &
|
|
get_nleafs(), ed, eh, em, es
|
|
end if
|
|
flush(output_unit)
|
|
|
|
end if
|
|
|
|
proceed = (nsteps <= nmax) .and. (time < tmax) .and. proceed
|
|
|
|
do while(proceed)
|
|
|
|
call advance(status)
|
|
if (status /= 0) &
|
|
call print_message(loc, "Advancing the system failed!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
time = time + dt
|
|
step = step + 1
|
|
nsteps = nsteps + 1
|
|
|
|
call update_dtp()
|
|
|
|
call new_time_step()
|
|
|
|
call store_system(status)
|
|
if (status /= 0) &
|
|
call print_message(loc, "Could not store the current system state!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
tclc(1) = get_timer_total()
|
|
thrs = tclc(1) / 3.6d+03
|
|
|
|
call write_restart_snapshot(thrs, name, nrun, status)
|
|
if (status /= 0) &
|
|
call print_message(loc, "Could not store the restart snapshot!")
|
|
if (check_status(status /= 0)) return
|
|
|
|
proceed = (nsteps <= nmax) .and. (time < tmax) .and. (thrs <= trun)
|
|
#ifdef SIGNALS
|
|
proceed = proceed .and. .not. check_status(quit /= 0)
|
|
#endif /* SIGNALS */
|
|
|
|
if (verbose) then
|
|
|
|
tsim(1) = time
|
|
tsim = cshift(tsim, 1)
|
|
tclc = cshift(tclc, 1)
|
|
|
|
if (time >= tmax .or. (tclc(8) - tprn) >= 1.0d+00) then
|
|
|
|
ss = (tclc(8) - tclc(1)) * ((tmax - time) / &
|
|
max(1.0d-08 * tmax, time - tsim(1)))
|
|
es = int(max(0.0d+00, min(8.63999999d+08, ss)), kind=4)
|
|
ed = es / 86400
|
|
es = mod(es, 86400)
|
|
eh = es / 3600
|
|
es = mod(es, 3600)
|
|
em = es / 60
|
|
es = mod(es, 60)
|
|
|
|
if (error_control) then
|
|
write(output_unit, sfmt, advance="no") step, time, dt, &
|
|
errtol, nrejections, &
|
|
get_nleafs(), ed, eh, em, es
|
|
else
|
|
write(output_unit, sfmt, advance="no") step, time, dt, &
|
|
get_nleafs(), ed, eh, em, es
|
|
end if
|
|
flush(output_unit)
|
|
|
|
tprn = tclc(8)
|
|
end if
|
|
end if
|
|
|
|
end do ! main loop
|
|
|
|
if (verbose) write(*,*)
|
|
|
|
call write_restart_snapshot(1.0d+16, name, nrun, status)
|
|
if (status /= 0) &
|
|
call print_message(loc, "Could not store the restart snapshot!")
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine evolve_system
|
|
!
|
|
!===============================================================================
|
|
!!
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
!!
|
|
!===============================================================================
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine STORE_SYSTEM:
|
|
! -----------------------
|
|
!
|
|
! Subroutine stores system's progress: statistics and snapshots.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! status - the subroutine call status;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine store_system(status)
|
|
|
|
use evolution , only : time
|
|
use io , only : write_snapshot
|
|
use statistics , only : store_statistics
|
|
use user_problem, only : user_statistics
|
|
|
|
implicit none
|
|
|
|
integer, intent(out) :: status
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
status = 0
|
|
|
|
call store_statistics()
|
|
call user_statistics(time)
|
|
call write_snapshot(name)
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine store_system
|
|
!
|
|
!===============================================================================
|
|
!
|
|
! subroutine RESTORE_PARAMETERS:
|
|
! -----------------------------
|
|
!
|
|
! Subroutine reads the parameters from the parameters' file, then
|
|
! restores them if the job is resumed.
|
|
!
|
|
! Arguments:
|
|
!
|
|
! status - the subroutine call status;
|
|
!
|
|
!===============================================================================
|
|
!
|
|
subroutine restore_parameters(status)
|
|
|
|
use io , only : restart_from_snapshot, restart_snapshot_number
|
|
use io , only : restore_snapshot_parameters
|
|
use parameters, only : get_parameter
|
|
|
|
implicit none
|
|
|
|
integer, intent(out) :: status
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
status = 0
|
|
|
|
resumed = restart_from_snapshot()
|
|
|
|
nrun = max(1, restart_snapshot_number() + 1)
|
|
|
|
dbnds(:,1) = 0.0d+00
|
|
dbnds(:,2) = 1.0d+00
|
|
|
|
if (resumed) then
|
|
call restore_snapshot_parameters(status)
|
|
if (status /=0) return
|
|
end if
|
|
|
|
call get_parameter("problem", name)
|
|
call get_parameter("equation_system" , eqsys)
|
|
call get_parameter("equation_of_state", eos)
|
|
call get_parameter("ncells", ncells)
|
|
call get_parameter("maxlev", maxlev)
|
|
call get_parameter("xblocks", bdims(1))
|
|
call get_parameter("yblocks", bdims(2))
|
|
#if NDIMS == 3
|
|
call get_parameter("zblocks", bdims(3))
|
|
#endif /* NDIMS == 3 */
|
|
call get_parameter("xmin", dbnds(1,1))
|
|
call get_parameter("xmax", dbnds(1,2))
|
|
call get_parameter("ymin", dbnds(2,1))
|
|
call get_parameter("ymax", dbnds(2,2))
|
|
#if NDIMS == 3
|
|
call get_parameter("zmin", dbnds(3,1))
|
|
call get_parameter("zmax", dbnds(3,2))
|
|
#endif /* NDIMS == 3 */
|
|
|
|
call get_parameter("rngtype", rngtype)
|
|
call get_parameter("tmax" , tmax)
|
|
call get_parameter("trun" , trun)
|
|
call get_parameter("tsav" , tsav)
|
|
|
|
trun = trun - tsav / 6.0d+01
|
|
|
|
!-------------------------------------------------------------------------------
|
|
!
|
|
end subroutine restore_parameters
|
|
|
|
!===============================================================================
|
|
!
|
|
end module system
|