MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_input_output.t
Go to the documentation of this file.
1 !> Module for reading input and writing output
3 
4  implicit none
5  public
6 
7  !> Version number of the .dat file output
8  integer, parameter :: version_number = 5
9 
10  !> List of compatible versions
11  integer, parameter :: compatible_versions(3) = [3, 4, 5]
12 
13  !> number of w found in dat files
14  integer :: nw_found
15 
16  !> tag for MPI message
17  integer, private :: itag
18 
19  !> coefficient for rk2_alfa
20  double precision, private :: rk2_alfa
21 
22  !> whether a manually inserted snapshot is saved
23  logical :: save_now
24 
25  ! Formats used in output
26  character(len=*), parameter :: fmt_r = 'es16.8' ! Default precision
27  character(len=*), parameter :: fmt_r2 = 'es10.2' ! Two digits
28  character(len=*), parameter :: fmt_i = 'i8' ! Integer format
29 
30  ! public methods
31  public :: count_ix
32  public :: create_output_file
34  public :: block_shape_io
35 
36 
37 contains
38 
39  !> Read the command line arguments passed to amrvac
40  subroutine read_arguments()
42  use mod_slice, only: slicenext
43 
44  integer :: len, stat, n, i, ipars
45  integer, parameter :: max_files = 20 ! Maximum number of par files
46  integer :: n_par_files
47  character(len=max_files*std_len) :: all_par_files
48  character(len=std_len) :: tmp_files(max_files), arg
49  logical :: unknown_arg, help, morepars
50 
51  if (mype == 0) then
52  print *, '-----------------------------------------------------------------------------'
53  print *, '-----------------------------------------------------------------------------'
54  print *, '| __ __ ____ ___ _ __ __ ______ ___ ____ |'
55  print *, '| | \/ | _ \_ _| / \ | \/ | _ \ \ / / \ / ___| |'
56  print *, '| | |\/| | |_) | |_____ / _ \ | |\/| | |_) \ \ / / _ \| | |'
57  print *, '| | | | | __/| |_____/ ___ \| | | | _ < \ V / ___ \ |___ |'
58  print *, '| |_| |_|_| |___| /_/ \_\_| |_|_| \_\ \_/_/ \_\____| |'
59  print *, '-----------------------------------------------------------------------------'
60  print *, '-----------------------------------------------------------------------------'
61  end if
62 
63  ! =============== Fortran 2003 command line reading ================
64 
65  ! Default command line arguments
66  all_par_files="amrvac.par"
68  snapshotnext=-1
69  slicenext=-1
70  collapsenext=-1
71  help=.false.
72  convert=.false.
73  resume_previous_run=.false.
74 
75  ! Argument 0 is program name, so we start from one
76  i = 1
77  unknown_arg=.false.
78  DO
79  CALL get_command_argument(i, arg)
80  IF (len_trim(arg) == 0) EXIT
81  select case(arg)
82  case("-i")
83  i = i+1
84  CALL get_command_argument(i, arg)
85  !if(mype==0)print *,'found argument=',arg
86  all_par_files=trim(arg)
87  morepars=.true.
88  ipars=1
89  do while (ipars<max_files.and.morepars)
90  CALL get_command_argument(i+ipars,arg)
91  !if(mype==0)print *,'found argument=',arg
92  if (index(trim(arg),"-")==1.or.len_trim(arg)==0) then
93  morepars=.false.
94  else
95  ipars=ipars+1
96  all_par_files=trim(all_par_files)//" "//trim(arg)
97  !if(mype==0)print *,'now all_par_files=',TRIM(all_par_files)
98  endif
99  end do
100  !if(mype==0)print *,'-i identified ',ipars, ' par-file arguments'
101  i=i+ipars-1
102  !if(mype==0)print *,'-i arguments passed to all_par_files=',TRIM(all_par_files)
103  case("-if")
104  i = i+1
105  CALL get_command_argument(i, arg)
106  restart_from_file=trim(arg)
107  !if(mype==0)print *,'-if has argument=',TRIM(arg),' passed to restart_from_file=',TRIM(restart_from_file)
108  case("-slicenext")
109  i = i+1
110  CALL get_command_argument(i, arg)
111  read(arg,*,iostat=stat) slicenext
112  !if(mype==0)print *,'-slicenext has argument=',arg,' passed to slicenext=',slicenext
113  case("-collapsenext")
114  i = i+1
115  CALL get_command_argument(i, arg)
116  read(arg,*,iostat=stat) collapsenext
117  !if(mype==0)print *,'-collapsenext has argument=',arg,' passed to collapsenext=',collapsenext
118  case("-snapshotnext")
119  i = i+1
120  CALL get_command_argument(i, arg)
121  read(arg,*,iostat=stat) snapshotnext
122  !if(mype==0)print *,'-snapshotnext has argument=',arg,' passed to snapshotnext=',snapshotnext
123  case("-resume")
124  resume_previous_run=.true.
125  !if(mype==0)print *,'resume specified: resume_previous_run=T'
126  case("-convert")
127  convert=.true.
128  if(mype==0)print *,'convert specified: convert=T'
129  case("--help","-help")
130  help=.true.
131  EXIT
132  case default
133  unknown_arg=.true.
134  help=.true.
135  EXIT
136  end select
137  i = i+1
138  END DO
139 
140  if (unknown_arg) then
141  print*,"======================================="
142  print*,"Error: Command line argument ' ",trim(arg)," ' not recognized"
143  print*,"======================================="
144  help=.true.
145  end if
146 
147  ! Show the usage if the help flag was given, or no par file was specified
148  if (help) then
149  if (mype == 0) then
150  print *, 'Usage example:'
151  print *, 'mpirun -np 4 ./amrvac -i file.par [file2.par ...]'
152  print *, ' (later .par files override earlier ones)'
153  print *, ''
154  print *, 'Optional arguments:'
155  print *, '-convert Convert snapshot files'
156  print *, '-if file0001.dat Use this snapshot to restart from'
157  print *, ' (you can modify e.g. output names)'
158  print *, '-resume Automatically resume previous run'
159  print *, ' (useful for long runs on HPC systems)'
160  print *, '-snapshotnext N Manual index for next snapshot'
161  print *, '-slicenext N Manual index for next slice output'
162  print *, '-collapsenext N Manual index for next collapsed output'
163  print *, ''
164  end if
165  call mpi_finalize(ierrmpi)
166  stop
167  end if
168 
169  ! Split the input files, in case multiple were given
170  call get_fields_string(all_par_files, " ,'"""//char(9), max_files, &
171  tmp_files, n_par_files)
172 
173  allocate(par_files(n_par_files))
174  par_files = tmp_files(1:n_par_files)
175 
176  end subroutine read_arguments
177 
178  !> Read in the user-supplied parameter-file
179  subroutine read_par_files()
182  use mod_small_values
183  use mod_limiter
184  use mod_slice
185  use mod_geometry
186  use mod_source
187 
188  logical :: fileopen, file_exists
189  integer :: i, j, k, ifile, io_state
190  integer :: iB, isave, iw, level, idim, islice
191  integer :: nx_vec(^ND), block_nx_vec(^ND)
192  integer :: my_unit, iostate
193  integer :: ilev
194  double precision :: dx_vec(^ND)
195 
196  character :: c_ndim
197  character(len=80) :: fmt_string
198  character(len=std_len) :: err_msg, restart_from_file_arg
199  character(len=std_len) :: basename_full, basename_prev, dummy_file
200  character(len=std_len), dimension(:), allocatable :: &
201  typeboundary_min^D, typeboundary_max^D
202  character(len=std_len), allocatable :: limiter(:)
203  character(len=std_len), allocatable :: gradient_limiter(:)
204  character(len=name_len) :: stretch_dim(ndim)
205  !> How to apply dimensional splitting to the source terms, see
206  !> @ref discretization.md
207  character(len=std_len) :: typesourcesplit
208  !> Which flux scheme of spatial discretization to use (per grid level)
209  character(len=std_len), allocatable :: flux_scheme(:)
210  !> Which type of the maximal bound speed of Riemann fan to use
211  character(len=std_len) :: typeboundspeed
212  !> Which time stepper to use
213  character(len=std_len) :: time_stepper
214  !> Which time integrator to use
215  character(len=std_len) :: time_integrator
216  !> type of curl operator
217  character(len=std_len) :: typecurl
218  !> Limiter used for prolongation to refined grids and ghost cells
219  character(len=std_len) :: typeprolonglimit
220  !> How to compute the CFL-limited time step.
221  !> Options are 'maxsum': max(sum(c/dx)); 'summax': sum(max(c/dx)) and
222  !> 'minimum: max(c/dx), where the summations loop over the grid dimensions and
223  !> c is the velocity. The default 'maxsum' is the conventiontal way of
224  !> computing CFL-limited time steps.
225  character(len=std_len) :: typecourant
226 
227  double precision, dimension(nsavehi) :: tsave_log, tsave_dat, tsave_slice, &
228  tsave_collapsed, tsave_custom
229  double precision :: dtsave_log, dtsave_dat, dtsave_slice, &
230  dtsave_collapsed, dtsave_custom
231  integer :: ditsave_log, ditsave_dat, ditsave_slice, &
232  ditsave_collapsed, ditsave_custom
233  double precision :: tsavestart_log, tsavestart_dat, tsavestart_slice, &
234  tsavestart_collapsed, tsavestart_custom
235  integer :: windex, ipower
236  double precision :: sizeuniformpart^D
237  double precision :: im_delta,im_nu,rka54,rka51,rkb54,rka55
238 
239  namelist /filelist/ base_filename,restart_from_file, &
247 
248  namelist /savelist/ tsave,itsave,dtsave,ditsave,nslices,slicedir, &
250  tsave_log, tsave_dat, tsave_slice, tsave_collapsed, tsave_custom, &
251  dtsave_log, dtsave_dat, dtsave_slice, dtsave_collapsed, dtsave_custom, &
252  ditsave_log, ditsave_dat, ditsave_slice, ditsave_collapsed, ditsave_custom,&
253  tsavestart_log, tsavestart_dat, tsavestart_slice, tsavestart_collapsed,&
254  tsavestart_custom, tsavestart
255 
258 
259  namelist /methodlist/ time_stepper,time_integrator, &
260  source_split_usr,typesourcesplit,&
261  dimsplit,typedimsplit,flux_scheme,&
262  limiter,gradient_limiter,cada3_radius,&
263  loglimit,typeboundspeed, h_correction,&
265  typegrad,typediv,typecurl,&
267  flatcd,flatsh,&
272  schmid_rad^d
273 
274  namelist /boundlist/ nghostcells,ghost_copy,&
276 
277  namelist /meshlist/ refine_max_level,nbufferx^d,refine_threshold,&
279  stretch_dim, stretch_uncentered, &
283  typeprolonglimit, &
285  namelist /paramlist/ courantpar, dtpar, dtdiffpar, &
286  typecourant, slowsteps
287 
288  namelist /emissionlist/ filename_euv,wavelength,resolution_euv,&
293 
294  ! default maximum number of grid blocks in a processor
295  max_blocks=4000
296 
297  ! allocate cell size of all levels
298  allocate(dx(ndim,nlevelshi))
299  {allocate(dg^d(nlevelshi))\}
300  {allocate(ng^d(nlevelshi))\}
301 
302  ! default block size excluding ghost cells
303  {block_nx^d = 16\}
304 
305  ! default resolution of level-1 mesh (full domain)
306  {domain_nx^d = 32\}
307 
308  !! default number of ghost-cell layers at each boundary of a block
309  ! this is now done when the variable is defined in mod_global_parameters
310  ! the physics modules might set this variable in their init subroutine called earlier
311  !nghostcells = 2
312 
313  ! Allocate boundary conditions arrays in new and old style
314  {
315  allocate(typeboundary_min^d(nwfluxbc))
316  allocate(typeboundary_max^d(nwfluxbc))
317  typeboundary_min^d = undefined
318  typeboundary_max^d = undefined
319  }
320 
321  allocate(typeboundary(nwflux+nwaux,2*ndim))
322  typeboundary=0
323 
324  ! not save physical boundary in dat files by default
325  save_physical_boundary = .false.
326 
327  internalboundary = .false.
328 
329  ! defaults for specific options
330  typegrad = 'central'
331  typediv = 'central'
332  typecurl = 'central'
333 
334  ! defaults for smallest physical values allowed
335  small_temperature = 0.d0
336  small_pressure = 0.d0
337  small_density = 0.d0
338 
339  allocate(small_values_fix_iw(nw))
340  small_values_fix_iw(:) = .true.
341 
342  ! defaults for convert behavior
343 
344  ! store the -if value from argument in command line
345  restart_from_file_arg = restart_from_file
346  nwauxio = 0
347  nocartesian = .false.
348  saveprim = .false.
349  autoconvert = .false.
350  convert_type = 'vtuBCCmpi'
351  slice_type = 'vtuCC'
352  collapse_type = 'vti'
353  allocate(w_write(nw))
354  w_write(1:nw) = .true.
355  allocate(writelevel(nlevelshi))
356  writelevel(1:nlevelshi) = .true.
357  writespshift(1:ndim,1:2) = zero
358  level_io = -1
359  level_io_min = 1
361  ! endianness: littleendian (default) is 1, bigendian otherwise
362  type_endian = 1
363 
364  ! normalization of primitive variables: only for output
365  ! note that length_convert_factor is for length
366  ! this scaling is optional, and must be set consistently if used
367  allocate(w_convert_factor(nw))
368  w_convert_factor(:) = 1.0d0
369  time_convert_factor = 1.0d0
370  length_convert_factor = 1.0d0
371 
372  ! AMR related defaults
373  refine_max_level = 1
374  {nbufferx^d = 0\}
375  allocate(refine_threshold(nlevelshi))
376  refine_threshold(1:nlevelshi) = 0.1d0
377  allocate(derefine_ratio(nlevelshi))
378  derefine_ratio(1:nlevelshi) = 1.0d0/8.0d0
379  typeprolonglimit = 'default'
380  refine_criterion = 3
381  allocate(w_refine_weight(nw+1))
382  w_refine_weight = 0.d0
383  allocate(logflag(nw+1))
384  logflag = .false.
385  allocate(amr_wavefilter(nlevelshi))
386  amr_wavefilter(1:nlevelshi) = 1.0d-2
387  tfixgrid = bigdouble
388  itfixgrid = biginteger
389  ditregrid = 1
390 
391  ! Grid stretching defaults
392  stretch_uncentered = .true.
393  stretch_dim(1:ndim) = undefined
394  qstretch_baselevel(1:ndim) = bigdouble
396 
397  ! IO defaults
398  it_init = 0
399  it_max = biginteger
400  time_init = 0.d0
401  time_max = bigdouble
402  wall_time_max = bigdouble
403  final_dt_reduction=.true.
404  final_dt_exit=.false.
405  dtmin = 1.0d-10
406  nslices = 0
407  collapse = .false.
408  collapselevel = 1
409  time_between_print = 30.0d0 ! Print status every 30 seconds
410 
411  do ifile=1,nfile
412  do isave=1,nsavehi
413  tsave(isave,ifile) = bigdouble ! global_time of saves into the output files
414  itsave(isave,ifile) = biginteger ! it of saves into the output files
415  end do
416  dtsave(ifile) = bigdouble ! time between saves
417  ditsave(ifile) = biginteger ! timesteps between saves
418  isavet(ifile) = 1 ! index for saves by global_time
419  isaveit(ifile) = 1 ! index for saves by it
420  tsavestart(ifile) = 0.0d0
421  end do
422 
423  tsave_log = bigdouble
424  tsave_dat = bigdouble
425  tsave_slice = bigdouble
426  tsave_collapsed = bigdouble
427  tsave_custom = bigdouble
428 
429  dtsave_log = bigdouble
430  dtsave_dat = bigdouble
431  dtsave_slice = bigdouble
432  dtsave_collapsed = bigdouble
433  dtsave_custom = bigdouble
434 
435  ditsave_log = biginteger
436  ditsave_dat = biginteger
437  ditsave_slice = biginteger
438  ditsave_collapsed = biginteger
439  ditsave_custom = biginteger
440 
441  tsavestart_log = bigdouble
442  tsavestart_dat = bigdouble
443  tsavestart_slice = bigdouble
444  tsavestart_collapsed = bigdouble
445  tsavestart_custom = bigdouble
446 
447  typefilelog = 'default'
448 
449  ! defaults for input
450  reset_time = .false.
451  reset_it = .false.
452  firstprocess = .false.
453  reset_grid = .false.
454  base_filename = 'data'
455  usr_filename = ''
456 
457 
458  ! Defaults for discretization methods
459  typeaverage = 'default'
460  tvdlfeps = one
461  nxdiffusehllc = 0
462  flathllc = .false.
463  slowsteps = -1
464  courantpar = 0.8d0
465  typecourant = 'maxsum'
466  dimsplit = .false.
467  typedimsplit = 'default'
468  if(physics_type=='mhd') then
469  cada3_radius = 0.5d0
470  else
471  cada3_radius = 0.1d0
472  end if
473  {schmid_rad^d = 1.d0\}
474  typetvd = 'roe'
475  typeboundspeed = 'Einfeldt'
476  source_split_usr= .false.
477  time_stepper = 'twostep'
478  time_integrator = 'default'
479  angmomfix = .false.
480  ! default PC or explicit midpoint, hence alfa=0.5
481  rk2_alfa = half
482  ! default IMEX-RK22Ln hence lambda = 1 - 1/sqrt(2)
483  imex222_lambda = 1.0d0 - 1.0d0 / dsqrt(2.0d0)
484  ! default SSPRK(3,3) or Gottlieb-Shu 1998 for threestep
485  ! default SSPRK(4,3) or Spireti-Ruuth for fourstep
486  ! default SSPRK(5,4) using Gottlieb coeffs
487  ssprk_order = 3
488  ! default RK3 butcher table: Heun 3rd order
489  rk3_switch = 3
490  ! default IMEX threestep is IMEX_ARK(232)
491  imex_switch = 1
492 
493  ! Defaults for synthesing emission
494  los_theta = 0.d0
495  los_phi = 0.d0
496  image_rotate = 0.d0
497  x_origin = 0.d0
498  big_image = .false.
499  location_slit = 0.d0
500  resolution_euv = 'instrument'
501  resolution_sxr = 'instrument'
502  resolution_spectrum = 'instrument'
503 
504  allocate(flux_scheme(nlevelshi),typepred1(nlevelshi),flux_method(nlevelshi))
505  allocate(limiter(nlevelshi),gradient_limiter(nlevelshi))
506  do level=1,nlevelshi
507  flux_scheme(level) = 'tvdlf'
508  typepred1(level) = 0
509  limiter(level) = 'minmod'
510  gradient_limiter(level) = 'minmod'
511  end do
512 
513  flatcd = .false.
514  flatsh = .false.
515  typesourcesplit = 'sfs'
516  allocate(loglimit(nw))
517  loglimit(1:nw) = .false.
518 
519  allocate(typeentropy(nw))
520 
521  do iw=1,nw
522  typeentropy(iw)='nul' ! Entropy fix type
523  end do
524 
525  dtdiffpar = 0.5d0
526  dtpar = -1.d0
527 
528  ! problem setup defaults
529  iprob = 1
530 
531  ! end defaults
532 
533  ! Initialize Kronecker delta, and Levi-Civita tensor
534  do i=1,3
535  do j=1,3
536  if(i==j)then
537  kr(i,j)=1
538  else
539  kr(i,j)=0
540  endif
541  do k=1,3
542  if(i==j.or.j==k.or.k==i)then
543  lvc(i,j,k)=0
544  else if(i+1==j.or.i-2==j)then
545  lvc(i,j,k)=1
546  else
547  lvc(i,j,k)=-1
548  endif
549  enddo
550  enddo
551  enddo
552 
553  ! These are used to construct file and log names from multiple par files
554  basename_full = ''
555  basename_prev = ''
556 
557  do i = 1, size(par_files)
558  if (mype == 0) print *, "Reading " // trim(par_files(i))
559 
560  ! Check whether the file exists
561  inquire(file=trim(par_files(i)), exist=file_exists)
562 
563  if (.not. file_exists) then
564  write(err_msg, *) "The parameter file " // trim(par_files(i)) // &
565  " does not exist"
566  call mpistop(trim(err_msg))
567  end if
568 
569  open(unitpar, file=trim(par_files(i)), status='old')
570 
571  ! Try to read in the namelists. They can be absent or in a different
572  ! order, since we rewind before each read.
573  rewind(unitpar)
574  read(unitpar, filelist, end=101)
575 
576 101 rewind(unitpar)
577  read(unitpar, savelist, end=102)
578 
579 102 rewind(unitpar)
580  read(unitpar, stoplist, end=103)
581 
582 103 rewind(unitpar)
583  read(unitpar, methodlist, end=104)
584 
585 104 rewind(unitpar)
586  read(unitpar, boundlist, end=105)
587 
588 105 rewind(unitpar)
589  read(unitpar, meshlist, end=106)
590 
591 106 rewind(unitpar)
592  read(unitpar, paramlist, end=107)
593 
594 107 rewind(unitpar)
595  read(unitpar, emissionlist, end=108)
596 
597 108 close(unitpar)
598 
599  ! Append the log and file names given in the par files
600  if (base_filename /= basename_prev) &
601  basename_full = trim(basename_full) // trim(base_filename)
602  basename_prev = base_filename
603  end do
604 
605  base_filename = basename_full
606 
607  ! Check whether output directory is writable
608  if(mype==0) then
609  dummy_file = trim(base_filename)//"DUMMY"
610  open(newunit=my_unit, file=trim(dummy_file), iostat=iostate)
611  if (iostate /= 0) then
612  call mpistop("Can't write to output directory (" // &
613  trim(base_filename) // ")")
614  else
615  close(my_unit, status='delete')
616  end if
617  end if
618 
619  ! restart filename from command line overwrites the one in par file
620  if(restart_from_file_arg /= undefined) &
621  restart_from_file=restart_from_file_arg
622 
623  ! Root process will search snapshot
624  if (mype == 0) then
625  if(restart_from_file == undefined) then
626  do index_latest_data = -1, 9998
627  ! Check if the next snapshot is missing
628  if (.not. snapshot_exists(index_latest_data+1)) exit
629  end do
630 
631  if (index_latest_data == -1) then
632  ! If initial data is missing (e.g. moved due to lack of space),
633  ! search file with highest index
634  do index_latest_data = 9999, 0, -1
635  if (snapshot_exists(index_latest_data)) exit
636  end do
637  end if
638  else
639  ! get index of the given data restarted from
640  index_latest_data=get_snapshot_index(trim(restart_from_file))
641  end if
642  end if
643  call mpi_bcast(index_latest_data, 1, mpi_integer, 0, icomm, ierrmpi)
644 
645  if (resume_previous_run) then
646  if (index_latest_data == -1) then
647  if(mype==0) write(*,*) "No snapshots found to resume from, start a new run..."
648  else
649  ! Set file name to restart from
650  write(restart_from_file, "(a,i4.4,a)") trim(base_filename),index_latest_data, ".dat"
651  end if
652  end if
653 
654  if (restart_from_file == undefined) then
655  snapshotnext = 0
656  slicenext = 0
657  collapsenext = 0
658  if (firstprocess) &
659  call mpistop("Please restart from a snapshot when firstprocess=T")
660  if (convert) &
661  call mpistop('Change convert to .false. for a new run!')
662  end if
663 
664  if (small_pressure < 0.d0) call mpistop("small_pressure should be positive.")
665  if (small_density < 0.d0) call mpistop("small_density should be positive.")
666  ! Give priority to non-zero small temperature
667  if (small_temperature>0.d0) small_pressure=small_density*small_temperature
668 
669  if(convert) autoconvert=.false.
670 
671  where (tsave_log < bigdouble) tsave(:, 1) = tsave_log
672  where (tsave_dat < bigdouble) tsave(:, 2) = tsave_dat
673  where (tsave_slice < bigdouble) tsave(:, 3) = tsave_slice
674  where (tsave_collapsed < bigdouble) tsave(:, 4) = tsave_collapsed
675  where (tsave_custom < bigdouble) tsave(:, 5) = tsave_custom
676 
677  if (dtsave_log < bigdouble) dtsave(1) = dtsave_log
678  if (dtsave_dat < bigdouble) dtsave(2) = dtsave_dat
679  if (dtsave_slice < bigdouble) dtsave(3) = dtsave_slice
680  if (dtsave_collapsed < bigdouble) dtsave(4) = dtsave_collapsed
681  if (dtsave_custom < bigdouble) dtsave(5) = dtsave_custom
682 
683  if (tsavestart_log < bigdouble) tsavestart(1) = tsavestart_log
684  if (tsavestart_dat < bigdouble) tsavestart(2) = tsavestart_dat
685  if (tsavestart_slice < bigdouble) tsavestart(3) = tsavestart_slice
686  if (tsavestart_collapsed < bigdouble) tsavestart(4) = tsavestart_collapsed
687  if (tsavestart_custom < bigdouble) tsavestart(5) = tsavestart_custom
688 
689  if (ditsave_log < bigdouble) ditsave(1) = ditsave_log
690  if (ditsave_dat < bigdouble) ditsave(2) = ditsave_dat
691  if (ditsave_slice < bigdouble) ditsave(3) = ditsave_slice
692  if (ditsave_collapsed < bigdouble) ditsave(4) = ditsave_collapsed
693  if (ditsave_custom < bigdouble) ditsave(5) = ditsave_custom
694  ! convert hours to seconds for ending wall time
695  if (wall_time_max < bigdouble) wall_time_max=wall_time_max*3600.d0
696 
697  if (mype == 0) then
698  write(unitterm, *) ''
699  write(unitterm, *) 'Output type | tsavestart | dtsave | ditsave | itsave(1) | tsave(1)'
700  write(fmt_string, *) '(A12," | ",E9.3E2," | ",E9.3E2," | ",I6," | "'//&
701  ',I6, " | ",E9.3E2)'
702  end if
703 
704  do ifile = 1, nfile
705  if (mype == 0) write(unitterm, fmt_string) trim(output_names(ifile)), &
706  tsavestart(ifile), dtsave(ifile), ditsave(ifile), itsave(1, ifile), tsave(1, ifile)
707  end do
708 
709  if (mype == 0) write(unitterm, *) ''
710 
711  do islice=1,nslices
712  if(slicedir(islice) > ndim) &
713  write(uniterr,*)'Warning in read_par_files: ', &
714  'Slice ', islice,' direction',slicedir(islice),'larger than ndim=',ndim
715  if(slicedir(islice) < 1) &
716  write(uniterr,*)'Warning in read_par_files: ', &
717  'Slice ', islice,' direction',slicedir(islice),'too small, should be [',1,ndim,']'
718  end do
719 
720  if(it_max==biginteger .and. time_max==bigdouble.and.mype==0) write(uniterr,*) &
721  'Warning in read_par_files: it_max or time_max not given!'
722 
723  select case (typecourant)
724  case ('maxsum')
725  type_courant=type_maxsum
726  case ('summax')
727  type_courant=type_summax
728  case ('minimum')
729  type_courant=type_minimum
730  case default
731  write(unitterm,*)'Unknown typecourant=',typecourant
732  call mpistop("Error from read_par_files: no such typecourant!")
733  end select
734 
735  do level=1,nlevelshi
736  select case (flux_scheme(level))
737  case ('hll')
738  flux_method(level)=fs_hll
739  case ('hllc')
740  flux_method(level)=fs_hllc
741  case ('hlld')
742  flux_method(level)=fs_hlld
743  case ('hllcd')
744  flux_method(level)=fs_hllcd
745  case ('tvdlf')
746  flux_method(level)=fs_tvdlf
747  case ('tvdmu')
748  flux_method(level)=fs_tvdmu
749  case ('tvd')
750  flux_method(level)=fs_tvd
751  case ('cd')
752  flux_method(level)=fs_cd
753  case ('cd4')
754  flux_method(level)=fs_cd4
755  case ('fd')
756  flux_method(level)=fs_fd
757  case ('source')
758  flux_method(level)=fs_source
759  case ('nul','null')
760  flux_method(level)=fs_nul
761  case default
762  call mpistop("unkown or bad flux scheme")
763  end select
764  if(flux_scheme(level)=='tvd'.and.time_stepper/='onestep') &
765  call mpistop(" tvd is onestep method, reset time_stepper='onestep'")
766  if(flux_scheme(level)=='tvd')then
767  if(mype==0.and.(.not.dimsplit)) write(unitterm,*) &
768  'Warning: setting dimsplit=T for tvd, as used for level=',level
769  dimsplit=.true.
770  endif
771  if(flux_scheme(level)=='hlld'.and.physics_type/='mhd' .and. physics_type/='twofl') &
772  call mpistop("Cannot use hlld flux if not using MHD or 2FL only charges physics!")
773 
774  if(flux_scheme(level)=='hllc'.and.physics_type=='mf') &
775  call mpistop("Cannot use hllc flux if using magnetofriction physics!")
776 
777  if(flux_scheme(level)=='tvd'.and.physics_type=='mf') &
778  call mpistop("Cannot use tvd flux if using magnetofriction physics!")
779 
780  if(flux_scheme(level)=='tvdmu'.and.physics_type=='mf') &
781  call mpistop("Cannot use tvdmu flux if using magnetofriction physics!")
782 
783  if (typepred1(level)==0) then
784  select case (flux_scheme(level))
785  case ('cd')
786  typepred1(level)=fs_cd
787  case ('cd4')
788  typepred1(level)=fs_cd4
789  case ('fd')
790  typepred1(level)=fs_fd
791  case ('tvdlf','tvdmu')
792  typepred1(level)=fs_hancock
793  case ('hll')
794  typepred1(level)=fs_hll
795  case ('hllc')
796  typepred1(level)=fs_hllc
797  case ('hllcd')
798  typepred1(level)=fs_hllcd
799  case ('hlld')
800  typepred1(level)=fs_hlld
801  case ('nul','source','tvd')
802  typepred1(level)=fs_nul
803  case default
804  call mpistop("No default predictor for this full step")
805  end select
806  end if
807  end do
808 
809  ! finite difference scheme fd need global maximal speed
810  if(any(flux_scheme=='fd')) need_global_cmax=.true.
811  if(any(limiter=='schmid1')) need_global_a2max=.true.
812 
813  ! initialize type_curl
814  select case (typecurl)
815  case ("central")
816  type_curl=central
817  case ("Gaussbased")
818  type_curl=gaussbased
819  case ("Stokesbased")
820  type_curl=stokesbased
821  case default
822  write(unitterm,*) "typecurl=",typecurl
823  call mpistop("unkown type of curl operator in read_par_files")
824  end select
825 
826  ! initialize types of time stepper and time integrator
827  select case (time_stepper)
828  case ("onestep")
829  t_stepper=onestep
830  nstep=1
831  if (time_integrator=='default') then
832  time_integrator="Forward_Euler"
833  end if
834  select case (time_integrator)
835  case ("Forward_Euler")
836  t_integrator=forward_euler
837  case ("IMEX_Euler")
838  t_integrator=imex_euler
839  case ("IMEX_SP")
840  t_integrator=imex_sp
841  case default
842  write(unitterm,*) "time_integrator=",time_integrator,"time_stepper=",time_stepper
843  call mpistop("unkown onestep time_integrator in read_par_files")
844  end select
845  use_imex_scheme=(t_integrator==imex_euler.or.t_integrator==imex_sp)
846  case ("twostep")
847  t_stepper=twostep
848  nstep=2
849  if (time_integrator=='default') then
850  time_integrator="Predictor_Corrector"
851  endif
852  select case (time_integrator)
853  case ("Predictor_Corrector")
854  t_integrator=predictor_corrector
855  case ("RK2_alfa")
856  t_integrator=rk2_alf
857  case ("ssprk2")
858  t_integrator=ssprk2
859  case ("IMEX_Midpoint")
860  t_integrator=imex_midpoint
861  case ("IMEX_Trapezoidal")
862  t_integrator=imex_trapezoidal
863  case ("IMEX_222")
864  t_integrator=imex_222
865  case default
866  write(unitterm,*) "time_integrator=",time_integrator,"time_stepper=",time_stepper
867  call mpistop("unkown twostep time_integrator in read_par_files")
868  end select
869  use_imex_scheme=(t_integrator==imex_midpoint.or.t_integrator==imex_trapezoidal&
870  .or.t_integrator==imex_222)
871  if (t_integrator==rk2_alf) then
872  if(rk2_alfa<smalldouble.or.rk2_alfa>one)call mpistop("set rk2_alfa within [0,1]")
873  rk_a21=rk2_alfa
874  rk_b2=1.0d0/(2.0d0*rk2_alfa)
875  rk_b1=1.0d0-rk_b2
876  endif
877  case ("threestep")
878  t_stepper=threestep
879  nstep=3
880  if (time_integrator=='default') then
881  time_integrator='ssprk3'
882  endif
883  select case (time_integrator)
884  case ("ssprk3")
885  t_integrator=ssprk3
886  case ("RK3_BT")
887  t_integrator=rk3_bt
888  case ("IMEX_ARS3")
889  t_integrator=imex_ars3
890  case ("IMEX_232")
891  t_integrator=imex_232
892  case ("IMEX_CB3a")
893  t_integrator=imex_cb3a
894  case default
895  write(unitterm,*) "time_integrator=",time_integrator,"time_stepper=",time_stepper
896  call mpistop("unkown threestep time_integrator in read_par_files")
897  end select
898  if(t_integrator==rk3_bt) then
899  select case(rk3_switch)
900  case(1)
901  ! we code up Ralston 3rd order here
902  rk3_a21=1.0d0/2.0d0
903  rk3_a31=0.0d0
904  rk3_a32=3.0d0/4.0d0
905  rk3_b1=2.0d0/9.0d0
906  rk3_b2=1.0d0/3.0d0
907  case(2)
908  ! we code up RK-Wray 3rd order here
909  rk3_a21=8.0d0/15.0d0
910  rk3_a31=1.0d0/4.0d0
911  rk3_a32=5.0d0/12.0d0
912  rk3_b1=1.0d0/4.0d0
913  rk3_b2=0.0d0
914  case(3)
915  ! we code up Heun 3rd order here
916  rk3_a21=1.0d0/3.0d0
917  rk3_a31=0.0d0
918  rk3_a32=2.0d0/3.0d0
919  rk3_b1=1.0d0/4.0d0
920  rk3_b2=0.0d0
921  case(4)
922  ! we code up Nystrom 3rd order here
923  rk3_a21=2.0d0/3.0d0
924  rk3_a31=0.0d0
925  rk3_a32=2.0d0/3.0d0
926  rk3_b1=1.0d0/4.0d0
927  rk3_b2=3.0d0/8.0d0
928  case default
929  call mpistop("Unknown rk3_switch")
930  end select
931  ! the rest is fixed from above
932  rk3_b3=1.0d0-rk3_b1-rk3_b2
933  rk3_c2=rk3_a21
934  rk3_c3=rk3_a31+rk3_a32
935  endif
936  if(t_integrator==ssprk3) then
937  select case(ssprk_order)
938  case(3) ! this is SSPRK(3,3) Gottlieb-Shu
939  rk_beta11=1.0d0
940  rk_beta22=1.0d0/4.0d0
941  rk_beta33=2.0d0/3.0d0
942  rk_alfa21=3.0d0/4.0d0
943  rk_alfa31=1.0d0/3.0d0
944  rk_c2=1.0d0
945  rk_c3=1.0d0/2.0d0
946  case(2) ! this is SSP(3,2)
947  rk_beta11=1.0d0/2.0d0
948  rk_beta22=1.0d0/2.0d0
949  rk_beta33=1.0d0/3.0d0
950  rk_alfa21=0.0d0
951  rk_alfa31=1.0d0/3.0d0
952  rk_c2=1.0d0/2.0d0
953  rk_c3=1.0d0
954  case default
955  call mpistop("Unknown ssprk3_order")
956  end select
957  rk_alfa22=1.0d0-rk_alfa21
958  rk_alfa33=1.0d0-rk_alfa31
959  endif
960  if(t_integrator==imex_ars3) then
961  ars_gamma=(3.0d0+dsqrt(3.0d0))/6.0d0
962  endif
963  if(t_integrator==imex_232) then
964  select case(imex_switch)
965  case(1) ! this is IMEX_ARK(232)
966  im_delta=1.0d0-1.0d0/dsqrt(2.0d0)
967  im_nu=(3.0d0+2.0d0*dsqrt(2.0d0))/6.0d0
968  imex_a21=2.0d0*im_delta
969  imex_a31=1.0d0-im_nu
970  imex_a32=im_nu
971  imex_b1=1.0d0/(2.0d0*dsqrt(2.0d0))
972  imex_b2=1.0d0/(2.0d0*dsqrt(2.0d0))
973  imex_ha21=im_delta
974  imex_ha22=im_delta
975  case(2) ! this is IMEX_SSP(232)
976  ! doi 10.1002/2017MS001065 Rokhzadi et al
977  imex_a21=0.711664700366941d0
978  imex_a31=0.077338168947683d0
979  imex_a32=0.917273367886007d0
980  imex_b1=0.398930808264688d0
981  imex_b2=0.345755244189623d0
982  imex_ha21=0.353842865099275d0
983  imex_ha22=0.353842865099275d0
984  case default
985  call mpistop("Unknown imex_siwtch")
986  end select
987  imex_c2=imex_a21
988  imex_c3=imex_a31+imex_a32
989  imex_b3=1.0d0-imex_b1-imex_b2
990  endif
991  if(t_integrator==imex_cb3a) then
992  imex_c2 = 0.8925502329346865
993  imex_a22 = imex_c2
994  imex_ha21 = imex_c2
995  imex_c3 = imex_c2 / (6.0d0*imex_c2**2 - 3.0d0*imex_c2 + 1.0d0)
996  imex_ha32 = imex_c3
997  imex_b2 = (3.0d0*imex_c2 - 1.0d0) / (6.0d0*imex_c2**2)
998  imex_b3 = (6.0d0*imex_c2**2 - 3.0d0*imex_c2 + 1.0d0) / (6.0d0*imex_c2**2)
999  imex_a33 = (1.0d0/6.0d0 - imex_b2*imex_c2**2 - imex_b3*imex_c2*imex_c3) / (imex_b3*(imex_c3-imex_c2))
1000  imex_a32 = imex_c3 - imex_a33
1001  ! if (mype == 0) then
1002  ! write(*,*) "================================="
1003  ! write(*,*) "Asserting the order conditions..."
1004  ! ! First order convergence: OK
1005  ! ! Second order convergence
1006  ! write(*,*) -1.0d0/2.0d0 + imex_b2*imex_c2 + imex_b3*imex_c3
1007  ! ! Third order convergence
1008  ! write(*,*) -1.0d0/3.0d0 + imex_b2*imex_c2**2 + imex_b3*imex_c3**2
1009  ! write(*,*) -1.0d0/6.0d0 + imex_b3*imex_ha32*imex_c2
1010  ! write(*,*) -1.0d0/6.0d0 + imex_b2*imex_a22*imex_c2 + imex_b3*imex_a32*imex_c2 + imex_b3*imex_a33*imex_c3
1011  ! write(*,*) "================================="
1012  ! end if
1013  end if
1014  use_imex_scheme=(t_integrator==imex_ars3.or.t_integrator==imex_232.or.t_integrator==imex_cb3a)
1015  case ("fourstep")
1016  t_stepper=fourstep
1017  nstep=4
1018  if (time_integrator=='default') then
1019  time_integrator="ssprk4"
1020  endif
1021  select case (time_integrator)
1022  case ("ssprk4")
1023  t_integrator=ssprk4
1024  case ("rk4")
1025  t_integrator=rk4
1026  case default
1027  write(unitterm,*) "time_integrator=",time_integrator,"time_stepper=",time_stepper
1028  call mpistop("unkown fourstep time_integrator in read_par_files")
1029  end select
1030  if(t_integrator==ssprk4) then
1031  select case(ssprk_order)
1032  case(3) ! this is SSPRK(4,3) Spireti-Ruuth
1033  rk_beta11=1.0d0/2.0d0
1034  rk_beta22=1.0d0/2.0d0
1035  rk_beta33=1.0d0/6.0d0
1036  rk_beta44=1.0d0/2.0d0
1037  rk_alfa21=0.0d0
1038  rk_alfa31=2.0d0/3.0d0
1039  rk_alfa41=0.0d0
1040  rk_c2=1.0d0/2.0d0
1041  rk_c3=1.0d0
1042  rk_c4=1.0d0/2.0d0
1043  case(2) ! this is SSP(4,2)
1044  rk_beta11=1.0d0/3.0d0
1045  rk_beta22=1.0d0/3.0d0
1046  rk_beta33=1.0d0/3.0d0
1047  rk_beta44=1.0d0/4.0d0
1048  rk_alfa21=0.0d0
1049  rk_alfa31=0.0d0
1050  rk_alfa41=1.0d0/4.0d0
1051  rk_c2=1.0d0/3.0d0
1052  rk_c3=2.0d0/3.0d0
1053  rk_c4=1.0d0
1054  case default
1055  call mpistop("Unknown ssprk_order")
1056  end select
1057  rk_alfa22=1.0d0-rk_alfa21
1058  rk_alfa33=1.0d0-rk_alfa31
1059  rk_alfa44=1.0d0-rk_alfa41
1060  endif
1061  case ("fivestep")
1062  t_stepper=fivestep
1063  nstep=5
1064  if (time_integrator=='default') then
1065  time_integrator="ssprk5"
1066  end if
1067  select case (time_integrator)
1068  case ("ssprk5")
1069  t_integrator=ssprk5
1070  case default
1071  write(unitterm,*) "time_integrator=",time_integrator,"time_stepper=",time_stepper
1072  call mpistop("unkown fivestep time_integrator in read_par_files")
1073  end select
1074  if(t_integrator==ssprk5) then
1075  select case(ssprk_order)
1076  ! we use ssprk_order to intercompare the different coefficient choices
1077  case(3) ! From Gottlieb 2005
1078  rk_beta11=0.391752226571890d0
1079  rk_beta22=0.368410593050371d0
1080  rk_beta33=0.251891774271694d0
1081  rk_beta44=0.544974750228521d0
1082  rk_beta54=0.063692468666290d0
1083  rk_beta55=0.226007483236906d0
1084  rk_alfa21=0.444370493651235d0
1085  rk_alfa31=0.620101851488403d0
1086  rk_alfa41=0.178079954393132d0
1087  rk_alfa53=0.517231671970585d0
1088  rk_alfa54=0.096059710526147d0
1089 
1090  rk_alfa22=0.555629506348765d0
1091  rk_alfa33=0.379898148511597d0
1092  rk_alfa44=0.821920045606868d0
1093  rk_alfa55=0.386708617503269d0
1094  rk_alfa22=1.0d0-rk_alfa21
1095  rk_alfa33=1.0d0-rk_alfa31
1096  rk_alfa44=1.0d0-rk_alfa41
1097  rk_alfa55=1.0d0-rk_alfa53-rk_alfa54
1098  rk_c2=rk_beta11
1099  rk_c3=rk_alfa22*rk_c2+rk_beta22
1100  rk_c4=rk_alfa33*rk_c3+rk_beta33
1101  rk_c5=rk_alfa44*rk_c4+rk_beta44
1102  case(2) ! From Spireti-Ruuth
1103  rk_beta11=0.39175222700392d0
1104  rk_beta22=0.36841059262959d0
1105  rk_beta33=0.25189177424738d0
1106  rk_beta44=0.54497475021237d0
1107  !rk_beta54=0.06369246925946d0
1108  !rk_beta55=0.22600748319395d0
1109  rk_alfa21=0.44437049406734d0
1110  rk_alfa31=0.62010185138540d0
1111  rk_alfa41=0.17807995410773d0
1112  rk_alfa53=0.51723167208978d0
1113  !rk_alfa54=0.09605971145044d0
1114 
1115  rk_alfa22=1.0d0-rk_alfa21
1116  rk_alfa33=1.0d0-rk_alfa31
1117  rk_alfa44=1.0d0-rk_alfa41
1118 
1119  rka51=0.00683325884039d0
1120  rka54=0.12759831133288d0
1121  rkb54=0.08460416338212d0
1122  rk_beta54=rkb54-rk_beta44*rka51/rk_alfa41
1123  rk_alfa54=rka54-rk_alfa44*rka51/rk_alfa41
1124 
1125  rk_alfa55=1.0d0-rk_alfa53-rk_alfa54
1126  rk_c2=rk_beta11
1127  rk_c3=rk_alfa22*rk_c2+rk_beta22
1128  rk_c4=rk_alfa33*rk_c3+rk_beta33
1129  rk_c5=rk_alfa44*rk_c4+rk_beta44
1130  rk_beta55=1.0d0-rk_beta54-rk_alfa53*rk_c3-rk_alfa54*rk_c4-rk_alfa55*rk_c5
1131  case default
1132  call mpistop("Unknown ssprk_order")
1133  end select
1134  ! the following combinations must be unity
1135  !print *,rk_beta55+rk_beta54+rk_alfa53*rk_c3+rk_alfa54*rk_c4+rk_alfa55*rk_c5
1136  !print *,rk_alfa22+rk_alfa21
1137  !print *,rk_alfa33+rk_alfa31
1138  !print *,rk_alfa44+rk_alfa41
1139  !print *,rk_alfa55+rk_alfa53+rk_alfa54
1140  endif
1141  use_imex_scheme=.false.
1142  case default
1143  call mpistop("Unknown time_stepper in read_par_files")
1144  end select
1145 
1146  do i = 1, ndim
1147  select case (stretch_dim(i))
1148  case (undefined, 'none')
1149  stretch_type(i) = stretch_none
1150  stretched_dim(i) = .false.
1151  case ('uni','uniform')
1152  stretch_type(i) = stretch_uni
1153  stretched_dim(i) = .true.
1154  case ('symm', 'symmetric')
1155  stretch_type(i) = stretch_symm
1156  stretched_dim(i) = .true.
1157  case default
1158  stretch_type(i) = stretch_none
1159  stretched_dim(i) = .false.
1160  if (mype == 0) print *, 'Got stretch_type = ', stretch_type(i)
1161  call mpistop('Unknown stretch type')
1162  end select
1163  end do
1164 
1165  ! Harmonize the parameters for dimensional splitting and source splitting
1166  if(typedimsplit =='default'.and. dimsplit) typedimsplit='xyyx'
1167  if(typedimsplit =='default'.and..not.dimsplit) typedimsplit='unsplit'
1168  dimsplit = typedimsplit /='unsplit'
1169 
1170  ! initialize types of split-source addition
1171  select case (typesourcesplit)
1172  case ('sfs')
1173  sourcesplit=sourcesplit_sfs
1174  case ('sf')
1175  sourcesplit=sourcesplit_sf
1176  case ('ssf')
1177  sourcesplit=sourcesplit_ssf
1178  case ('ssfss')
1179  sourcesplit=sourcesplit_ssfss
1180  case default
1181  write(unitterm,*)'No such typesourcesplit=',typesourcesplit
1182  call mpistop("Error: Unknown typesourcesplit!")
1183  end select
1184 
1185  if(coordinate==-1) then
1186  coordinate=cartesian
1187  if(mype==0) then
1188  write(*,*) 'Warning: coordinate system is not specified!'
1189  write(*,*) 'call set_coordinate_system in usr_init in mod_usr.t'
1190  write(*,*) 'Now use Cartesian coordinate'
1191  end if
1192  end if
1193 
1194  if(coordinate==cartesian) then
1195  slab=.true.
1196  slab_uniform=.true.
1197  if(any(stretched_dim)) then
1198  coordinate=cartesian_stretched
1199  slab_uniform=.false.
1200  end if
1201  else
1202  slab=.false.
1203  slab_uniform=.false.
1204  end if
1205 
1206  if(coordinate==spherical) then
1207  if(dimsplit) then
1208  if(mype==0)print *,'Warning: spherical symmetry needs dimsplit=F, resetting'
1209  dimsplit=.false.
1210  end if
1211  end if
1212 
1213  if (ndim==1) dimsplit=.false.
1214 
1215  ! type limiter of prolongation
1216  select case(typeprolonglimit)
1217  case('unlimit')
1218  ! unlimited
1219  prolong_limiter=1
1220  case('minmod')
1221  prolong_limiter=2
1222  case('woodward')
1223  prolong_limiter=3
1224  case('koren')
1225  prolong_limiter=4
1226  case default
1227  prolong_limiter=0
1228  end select
1229 
1230  ! Type limiter is of integer type for performance
1231  allocate(type_limiter(nlevelshi))
1232  allocate(type_gradient_limiter(nlevelshi))
1233 
1234  do level=1,nlevelshi
1235  type_limiter(level) = limiter_type(limiter(level))
1236  type_gradient_limiter(level) = limiter_type(gradient_limiter(level))
1237  end do
1238 
1239  if (any(limiter(1:nlevelshi)== 'ppm')&
1240  .and.(flatsh.and.physics_type=='rho')) then
1241  call mpistop(" PPM with flatsh=.true. can not be used with physics_type='rho'!")
1242  end if
1243 
1244  ! Copy boundary conditions to typeboundary, which is used internally
1245  {
1246  do iw=1,nwfluxbc
1247  select case(typeboundary_min^d(iw))
1248  case("special")
1249  typeboundary(iw,2*^d-1)=bc_special
1250  case("cont")
1251  typeboundary(iw,2*^d-1)=bc_cont
1252  case("symm")
1253  typeboundary(iw,2*^d-1)=bc_symm
1254  case("asymm")
1255  typeboundary(iw,2*^d-1)=bc_asymm
1256  case("periodic")
1257  typeboundary(iw,2*^d-1)=bc_periodic
1258  case("aperiodic")
1259  typeboundary(iw,2*^d-1)=bc_aperiodic
1260  case("noinflow")
1261  typeboundary(iw,2*^d-1)=bc_noinflow
1262  case("pole")
1263  typeboundary(iw,2*^d-1)=12
1264  case("bc_data")
1265  typeboundary(iw,2*^d-1)=bc_data
1266  case("character")
1267  typeboundary(iw,2*^d-1)=bc_character
1268  case default
1269  write (unitterm,*) "Undefined boundarytype found in read_par_files", &
1270  typeboundary_min^d(iw),"for variable iw=",iw," and side iB=",2*^d-1
1271  end select
1272  end do
1273  do iw=1,nwfluxbc
1274  select case(typeboundary_max^d(iw))
1275  case("special")
1276  typeboundary(iw,2*^d)=bc_special
1277  case("cont")
1278  typeboundary(iw,2*^d)=bc_cont
1279  case("symm")
1280  typeboundary(iw,2*^d)=bc_symm
1281  case("asymm")
1282  typeboundary(iw,2*^d)=bc_asymm
1283  case("periodic")
1284  typeboundary(iw,2*^d)=bc_periodic
1285  case("aperiodic")
1286  typeboundary(iw,2*^d)=bc_aperiodic
1287  case("noinflow")
1288  typeboundary(iw,2*^d)=bc_noinflow
1289  case("pole")
1290  typeboundary(iw,2*^d)=12
1291  case("bc_data")
1292  typeboundary(iw,2*^d)=bc_data
1293  case("bc_character")
1294  typeboundary(iw,2*^d)=bc_character
1295  case default
1296  write (unitterm,*) "Undefined boundarytype found in read_par_files", &
1297  typeboundary_max^d(iw),"for variable iw=",iw," and side iB=",2*^d
1298  end select
1299  end do
1300  }
1301 
1302  ! psi, tracers take the same boundary type as the first variable
1303  if (nwfluxbc<nwflux) then
1304  do iw = nwfluxbc+1, nwflux
1305  typeboundary(iw,:) = typeboundary(1, :)
1306  end do
1307  end if
1308  ! auxiliary variables take the same boundary type as the first variable
1309  if (nwaux>0) then
1310  do iw = nwflux+1, nwflux+nwaux
1311  typeboundary(iw,:) = typeboundary(1, :)
1312  end do
1313  end if
1314 
1315  if (any(typeboundary == 0)) then
1316  call mpistop("Not all boundary conditions have been defined")
1317  end if
1318 
1319  do idim=1,ndim
1320  periodb(idim)=(any(typeboundary(:,2*idim-1:2*idim)==bc_periodic))
1321  aperiodb(idim)=(any(typeboundary(:,2*idim-1:2*idim)==bc_aperiodic))
1322  if (periodb(idim).or.aperiodb(idim)) then
1323  do iw=1,nwflux
1324  if (typeboundary(iw,2*idim-1) .ne. typeboundary(iw,2*idim)) &
1325  call mpistop("Wrong counterpart in periodic boundary")
1326 
1327  if (typeboundary(iw,2*idim-1) /= bc_periodic .and. &
1328  typeboundary(iw,2*idim-1) /= bc_aperiodic) then
1329  call mpistop("Each dimension should either have all "//&
1330  "or no variables periodic, some can be aperiodic")
1331  end if
1332  end do
1333  end if
1334  end do
1335  {^nooned
1336  do idim=1,ndim
1337  if(any(typeboundary(:,2*idim-1)==12)) then
1338  if(any(typeboundary(:,2*idim-1)/=12)) typeboundary(:,2*idim-1)=12
1339  if(phys_energy) then
1340  windex=2
1341  else
1342  windex=1
1343  end if
1344  typeboundary(:,2*idim-1)=bc_symm
1345  if(physics_type/='rho') then
1346  select case(coordinate)
1347  case(cylindrical)
1348  typeboundary(phi_+1,2*idim-1)=bc_asymm
1349  if(physics_type=='mhd') typeboundary(ndir+windex+phi_,2*idim-1)=bc_asymm
1350  case(spherical)
1351  typeboundary(3:ndir+1,2*idim-1)=bc_asymm
1352  if(physics_type=='mhd') typeboundary(ndir+windex+2:ndir+windex+ndir,2*idim-1)=bc_asymm
1353  case default
1354  call mpistop('Pole is in cylindrical, polar, spherical coordinates!')
1355  end select
1356  end if
1357  end if
1358  if(any(typeboundary(:,2*idim)==12)) then
1359  if(any(typeboundary(:,2*idim)/=12)) typeboundary(:,2*idim)=12
1360  if(phys_energy) then
1361  windex=2
1362  else
1363  windex=1
1364  end if
1365  typeboundary(:,2*idim)=bc_symm
1366  if(physics_type/='rho') then
1367  select case(coordinate)
1368  case(cylindrical)
1369  typeboundary(phi_+1,2*idim)=bc_asymm
1370  if(physics_type=='mhd') typeboundary(ndir+windex+phi_,2*idim)=bc_asymm
1371  case(spherical)
1372  typeboundary(3:ndir+1,2*idim)=bc_asymm
1373  if(physics_type=='mhd') typeboundary(ndir+windex+2:ndir+windex+ndir,2*idim)=bc_asymm
1374  case default
1375  call mpistop('Pole is in cylindrical, polar, spherical coordinates!')
1376  end select
1377  end if
1378  end if
1379  end do
1380  }
1381 
1382  if(any(limiter(1:nlevelshi)=='mp5')) then
1383  nghostcells=max(nghostcells,3)
1384  end if
1385 
1386  if(any(limiter(1:nlevelshi)=='weno5')) then
1387  nghostcells=max(nghostcells,3)
1388  end if
1389 
1390  if(any(limiter(1:nlevelshi)=='weno5nm')) then
1391  nghostcells=max(nghostcells,3)
1392  end if
1393 
1394  if(any(limiter(1:nlevelshi)=='wenoz5')) then
1395  nghostcells=max(nghostcells,3)
1396  end if
1397 
1398  if(any(limiter(1:nlevelshi)=='wenoz5nm')) then
1399  nghostcells=max(nghostcells,3)
1400  end if
1401 
1402  if(any(limiter(1:nlevelshi)=='wenozp5')) then
1403  nghostcells=max(nghostcells,3)
1404  end if
1405 
1406  if(any(limiter(1:nlevelshi)=='wenozp5nm')) then
1407  nghostcells=max(nghostcells,3)
1408  end if
1409 
1410  if(any(limiter(1:nlevelshi)=='teno5ad')) then
1411  nghostcells=max(nghostcells,3)
1412  end if
1413 
1414  if(any(limiter(1:nlevelshi)=='weno5cu6')) then
1415  nghostcells=max(nghostcells,3)
1416  end if
1417 
1418  if(any(limiter(1:nlevelshi)=='ppm')) then
1419  nghostcells=max(nghostcells,4)
1420  end if
1421 
1422  if(any(limiter(1:nlevelshi)=='weno7')) then
1423  nghostcells=max(nghostcells,4)
1424  end if
1425 
1426  if(any(limiter(1:nlevelshi)=='mpweno7')) then
1427  nghostcells=max(nghostcells,4)
1428  end if
1429 
1430  ! If a wider stencil is used, extend the number of ghost cells
1431  nghostcells = nghostcells + phys_wider_stencil
1432 
1433  ! prolongation in AMR for constrained transport MHD needs even number ghosts
1434  if(stagger_grid .and. refine_max_level>1 .and. mod(nghostcells,2)/=0) then
1435  nghostcells=nghostcells+1
1436  end if
1437 
1438  select case (coordinate)
1439  {^nooned
1440  case (spherical)
1441  xprob^lim^de=xprob^lim^de*two*dpi;
1442  \}
1443  case (cylindrical)
1444  {
1445  if (^d==phi_) then
1446  xprob^lim^d=xprob^lim^d*two*dpi;
1447  end if
1448  \}
1449  end select
1450 
1451  ! full block size including ghostcells
1452  {ixghi^d = block_nx^d + 2*nghostcells\}
1453  {ixgshi^d = ixghi^d\}
1454 
1455  nx_vec = [{domain_nx^d|, }]
1456  block_nx_vec = [{block_nx^d|, }]
1457 
1458  if (any(nx_vec < 4) .or. any(mod(nx_vec, 2) == 1)) &
1459  call mpistop('Grid size (domain_nx^D) has to be even and >= 4')
1460 
1461  if (any(block_nx_vec < 4) .or. any(mod(block_nx_vec, 2) == 1)) &
1462  call mpistop('Block size (block_nx^D) has to be even and >= 4')
1463 
1464  { if(mod(domain_nx^d,block_nx^d)/=0) &
1465  call mpistop('Grid (domain_nx^D) and block (block_nx^D) must be consistent') \}
1466 
1467  if(refine_max_level>nlevelshi.or.refine_max_level<1)then
1468  write(unitterm,*)'Error: refine_max_level',refine_max_level,'>nlevelshi ',nlevelshi
1469  call mpistop("Reset nlevelshi and recompile!")
1470  endif
1471 
1472  if (any(stretched_dim)) then
1473  allocate(qstretch(0:nlevelshi,1:ndim),dxfirst(0:nlevelshi,1:ndim),&
1474  dxfirst_1mq(0:nlevelshi,1:ndim),dxmid(0:nlevelshi,1:ndim))
1475  allocate(nstretchedblocks(1:nlevelshi,1:ndim))
1476  qstretch(0:nlevelshi,1:ndim)=0.0d0
1477  dxfirst(0:nlevelshi,1:ndim)=0.0d0
1478  nstretchedblocks(1:nlevelshi,1:ndim)=0
1479  {if (stretch_type(^d) == stretch_uni) then
1480  ! first some sanity checks
1481  if(qstretch_baselevel(^d)<1.0d0.or.qstretch_baselevel(^d)==bigdouble) then
1482  if(mype==0) then
1483  write(*,*) 'stretched grid needs finite qstretch_baselevel>1'
1484  write(*,*) 'will try default value for qstretch_baselevel in dimension', ^d
1485  endif
1486  if(xprobmin^d>smalldouble)then
1487  qstretch_baselevel(^d)=(xprobmax^d/xprobmin^d)**(1.d0/dble(domain_nx^d))
1488  else
1489  call mpistop("can not set qstretch_baselevel automatically")
1490  endif
1491  endif
1492  if(mod(block_nx^d,2)==1) &
1493  call mpistop("stretched grid needs even block size block_nxD")
1494  if(mod(domain_nx^d/block_nx^d,2)/=0) &
1495  call mpistop("number level 1 blocks in D must be even")
1496  qstretch(1,^d)=qstretch_baselevel(^d)
1497  dxfirst(1,^d)=(xprobmax^d-xprobmin^d) &
1498  *(1.0d0-qstretch(1,^d))/(1.0d0-qstretch(1,^d)**domain_nx^d)
1499  qstretch(0,^d)=qstretch(1,^d)**2
1500  dxfirst(0,^d)=dxfirst(1,^d)*(1.0d0+qstretch(1,^d))
1501  if(refine_max_level>1)then
1502  do ilev=2,refine_max_level
1503  qstretch(ilev,^d)=dsqrt(qstretch(ilev-1,^d))
1504  dxfirst(ilev,^d)=dxfirst(ilev-1,^d) &
1505  /(1.0d0+dsqrt(qstretch(ilev-1,^d)))
1506  enddo
1507  endif
1508  endif \}
1509  if(mype==0) then
1510  {if(stretch_type(^d) == stretch_uni) then
1511  write(*,*) 'Stretched dimension ', ^d
1512  write(*,*) 'Using stretched grid with qs=',qstretch(0:refine_max_level,^d)
1513  write(*,*) ' and first cell sizes=',dxfirst(0:refine_max_level,^d)
1514  endif\}
1515  end if
1516  {if(stretch_type(^d) == stretch_symm) then
1517  if(mype==0) then
1518  write(*,*) 'will apply symmetric stretch in dimension', ^d
1519  endif
1520  if(mod(block_nx^d,2)==1) &
1521  call mpistop("stretched grid needs even block size block_nxD")
1522  ! checks on the input variable nstretchedblocks_baselevel
1523  if(nstretchedblocks_baselevel(^d)==0) &
1524  call mpistop("need finite even number of stretched blocks at baselevel")
1525  if(mod(nstretchedblocks_baselevel(^d),2)==1) &
1526  call mpistop("need even number of stretched blocks at baselevel")
1527  if(qstretch_baselevel(^d)<1.0d0.or.qstretch_baselevel(^d)==bigdouble) &
1528  call mpistop('stretched grid needs finite qstretch_baselevel>1')
1529  ! compute stretched part to ensure uniform center
1530  ipower=(nstretchedblocks_baselevel(^d)/2)*block_nx^d
1531  if(nstretchedblocks_baselevel(^d)==domain_nx^d/block_nx^d)then
1532  xstretch^d=0.5d0*(xprobmax^d-xprobmin^d)
1533  else
1534  xstretch^d=(xprobmax^d-xprobmin^d) &
1535  /(2.0d0+dble(domain_nx^d-nstretchedblocks_baselevel(^d)*block_nx^d) &
1536  *(1.0d0-qstretch_baselevel(^d))/(1.0d0-qstretch_baselevel(^d)**ipower))
1537  endif
1538  if(xstretch^d>(xprobmax^d-xprobmin^d)*0.5d0) &
1539  call mpistop(" stretched grid part should not exceed full domain")
1540  dxfirst(1,^d)=xstretch^d*(1.0d0-qstretch_baselevel(^d)) &
1541  /(1.0d0-qstretch_baselevel(^d)**ipower)
1542  nstretchedblocks(1,^d)=nstretchedblocks_baselevel(^d)
1543  qstretch(1,^d)=qstretch_baselevel(^d)
1544  qstretch(0,^d)=qstretch(1,^d)**2
1545  dxfirst(0,^d)=dxfirst(1,^d)*(1.0d0+qstretch(1,^d))
1546  dxmid(1,^d)=dxfirst(1,^d)
1547  dxmid(0,^d)=dxfirst(1,^d)*2.0d0
1548  if(refine_max_level>1)then
1549  do ilev=2,refine_max_level
1550  nstretchedblocks(ilev,^d)=2*nstretchedblocks(ilev-1,^d)
1551  qstretch(ilev,^d)=dsqrt(qstretch(ilev-1,^d))
1552  dxfirst(ilev,^d)=dxfirst(ilev-1,^d) &
1553  /(1.0d0+dsqrt(qstretch(ilev-1,^d)))
1554  dxmid(ilev,^d)=dxmid(ilev-1,^d)/2.0d0
1555  enddo
1556  endif
1557  ! sanity check on total domain size:
1558  sizeuniformpart^d=dxfirst(1,^d) &
1559  *(domain_nx^d-nstretchedblocks_baselevel(^d)*block_nx^d)
1560  if(mype==0) then
1561  print *,'uniform part of size=',sizeuniformpart^d
1562  print *,'setting of domain is then=',2*xstretch^d+sizeuniformpart^d
1563  print *,'versus=',xprobmax^d-xprobmin^d
1564  endif
1565  if(dabs(xprobmax^d-xprobmin^d-2*xstretch^d-sizeuniformpart^d)>smalldouble) then
1566  call mpistop('mismatch in domain size!')
1567  endif
1568  endif \}
1569  dxfirst_1mq(0:refine_max_level,1:ndim)=dxfirst(0:refine_max_level,1:ndim) &
1570  /(1.0d0-qstretch(0:refine_max_level,1:ndim))
1571  end if
1572 
1573  dx_vec = [{xprobmax^d-xprobmin^d|, }] / nx_vec
1574 
1575  if (mype==0) then
1576  write(c_ndim, '(I1)') ^nd
1577  write(unitterm, '(A30,' // c_ndim // '(I0," "))') &
1578  ' Domain size (cells): ', nx_vec
1579  write(unitterm, '(A30,' // c_ndim // '(E9.3," "))') &
1580  ' Level one dx: ', dx_vec
1581  end if
1582 
1583  if (any(dx_vec < smalldouble)) &
1584  call mpistop("Incorrect domain size (too small grid spacing)")
1585 
1586  dx(:, 1) = dx_vec
1587 
1588  if(sum(w_refine_weight(:))==0) w_refine_weight(1) = 1.d0
1589  if(dabs(sum(w_refine_weight(:))-1.d0)>smalldouble) then
1590  write(unitterm,*) "Sum of all elements in w_refine_weight be 1.d0"
1591  call mpistop("Reset w_refine_weight so the sum is 1.d0")
1592  end if
1593 
1594  select case (typeboundspeed)
1595  case('Einfeldt')
1596  boundspeed=1
1597  case('cmaxmean')
1598  boundspeed=2
1599  case('cmaxleftright')
1600  boundspeed=3
1601  case default
1602  call mpistop("set typeboundspeed='Einfieldt' or 'cmaxmean' or 'cmaxleftright'")
1603  end select
1604 
1605  if (mype==0) write(unitterm, '(A30)', advance='no') 'Refine estimation: '
1606 
1607  select case (refine_criterion)
1608  case (0)
1609  if (mype==0) write(unitterm, '(A)') "user defined"
1610  case (1)
1611  if (mype==0) write(unitterm, '(A)') "relative error"
1612  case (2)
1613  if (mype==0) write(unitterm, '(A)') "Lohner's original scheme"
1614  case (3)
1615  if (mype==0) write(unitterm, '(A)') "Lohner's scheme"
1616  case default
1617  call mpistop("Unknown error estimator, change refine_criterion")
1618  end select
1619 
1620  if (tfixgrid<bigdouble/2.0d0) then
1621  if(mype==0)print*,'Warning, at time=',tfixgrid,'the grid will be fixed'
1622  end if
1623  if (itfixgrid<biginteger/2) then
1624  if(mype==0)print*,'Warning, at iteration=',itfixgrid,'the grid will be fixed'
1625  end if
1626  if (ditregrid>1) then
1627  if(mype==0)print*,'Note, Grid is reconstructed once every',ditregrid,'iterations'
1628  end if
1629 
1630 
1631  do islice=1,nslices
1632  select case(slicedir(islice))
1633  {case(^d)
1634  if(slicecoord(islice)<xprobmin^d.or.slicecoord(islice)>xprobmax^d) &
1635  write(uniterr,*)'Warning in read_par_files: ', &
1636  'Slice ', islice, ' coordinate',slicecoord(islice),'out of bounds for dimension ',slicedir(islice)
1637  \}
1638  end select
1639  end do
1640 
1641  if (mype==0) then
1642  write(unitterm, '(A30,A,A)') 'restart_from_file: ', ' ', trim(restart_from_file)
1643  write(unitterm, '(A30,L1)') 'converting: ', convert
1644  write(unitterm, '(A)') ''
1645  endif
1646 
1647  deallocate(flux_scheme)
1648 
1649  end subroutine read_par_files
1650 
1651  !> Routine to find entries in a string
1652  subroutine get_fields_string(line, delims, n_max, fields, n_found, fully_read)
1653  !> The line from which we want to read
1654  character(len=*), intent(in) :: line
1655  !> A string with delimiters. For example delims = " ,'"""//char(9)
1656  character(len=*), intent(in) :: delims
1657  !> Maximum number of entries to read in
1658  integer, intent(in) :: n_max
1659  !> Number of entries found
1660  integer, intent(inout) :: n_found
1661  !> Fields in the strings
1662  character(len=*), intent(inout) :: fields(n_max)
1663  logical, intent(out), optional :: fully_read
1664 
1665  integer :: ixs_start(n_max)
1666  integer :: ixs_end(n_max)
1667  integer :: ix, ix_prev
1668 
1669  ix_prev = 0
1670  n_found = 0
1671 
1672  do while (n_found < n_max)
1673  ! Find the starting point of the next entry (a non-delimiter value)
1674  ix = verify(line(ix_prev+1:), delims)
1675  if (ix == 0) exit
1676 
1677  n_found = n_found + 1
1678  ixs_start(n_found) = ix_prev + ix ! This is the absolute position in 'line'
1679 
1680  ! Get the end point of the current entry (next delimiter index minus one)
1681  ix = scan(line(ixs_start(n_found)+1:), delims) - 1
1682 
1683  if (ix == -1) then ! If there is no last delimiter,
1684  ixs_end(n_found) = len(line) ! the end of the line is the endpoint
1685  else
1686  ixs_end(n_found) = ixs_start(n_found) + ix
1687  end if
1688 
1689  fields(n_found) = line(ixs_start(n_found):ixs_end(n_found))
1690  ix_prev = ixs_end(n_found) ! We continue to search from here
1691  end do
1692 
1693  if (present(fully_read)) then
1694  ix = verify(line(ix_prev+1:), delims)
1695  fully_read = (ix == 0) ! Are there only delimiters?
1696  end if
1697 
1698  end subroutine get_fields_string
1699 
1700  subroutine saveamrfile(ifile)
1701 
1704  use mod_particles, only: write_particles_snapshot
1705  use mod_slice, only: write_slice
1706  use mod_collapse, only: write_collapsed
1707  integer:: ifile
1708 
1709  select case (ifile)
1710  case (fileout_)
1711  ! Write .dat snapshot
1712  call write_snapshot()
1713 
1714  ! Generate formatted output (e.g., VTK)
1716 
1717  if(use_particles) call write_particles_snapshot()
1718 
1720  case (fileslice_)
1721  call write_slice
1722  case (filecollapse_)
1723  call write_collapsed
1724  case (filelog_)
1725  select case (typefilelog)
1726  case ('default')
1727  call printlog_default
1728  case ('regression_test')
1730  case ('special')
1731  if (.not. associated(usr_print_log)) then
1732  call mpistop("usr_print_log not defined")
1733  else
1734  call usr_print_log()
1735  end if
1736  case default
1737  call mpistop("Error in SaveFile: Unknown typefilelog")
1738  end select
1739  case (fileanalysis_)
1740  if (associated(usr_write_analysis)) then
1741  call usr_write_analysis()
1742  end if
1743  case default
1744  write(*,*) 'No save method is defined for ifile=',ifile
1745  call mpistop("")
1746  end select
1747 
1748  ! opedit: Flush stdout and stderr from time to time.
1749  flush(unit=unitterm)
1750 
1751  end subroutine saveamrfile
1752 
1753  !> Standard method for creating a new output file
1754  subroutine create_output_file(fh, ix, extension, suffix)
1756  integer, intent(out) :: fh !< File handle
1757  integer, intent(in) :: ix !< Index of file
1758  character(len=*), intent(in) :: extension !< Extension of file
1759  character(len=*), intent(in), optional :: suffix !< Optional suffix
1760  character(len=std_len) :: filename
1761  integer :: amode
1762 
1763  if (ix >= 10000) then
1764  call mpistop("Number of output files is limited to 10000 (0...9999)")
1765  end if
1766 
1767  if (present(suffix)) then
1768  write(filename,"(a,a,i4.4,a)") trim(base_filename), &
1769  trim(suffix), ix, extension
1770  else
1771  write(filename,"(a,i4.4,a)") trim(base_filename), ix, extension
1772  end if
1773 
1774  ! MPI cannot easily replace existing files
1775  open(unit=unitsnapshot,file=filename,status='replace')
1776  close(unitsnapshot, status='delete')
1777 
1778  amode = ior(mpi_mode_create, mpi_mode_wronly)
1779  call mpi_file_open(mpi_comm_self,filename,amode, &
1780  mpi_info_null, fh, ierrmpi)
1781 
1782  if (ierrmpi /= 0) then
1783  print *, "Error, cannot create file ", trim(filename)
1784  call mpistop("Fatal error")
1785  end if
1786 
1787  end subroutine create_output_file
1788 
1789  ! Check if a snapshot exists
1790  logical function snapshot_exists(ix)
1792  integer, intent(in) :: ix !< Index of snapshot
1793  character(len=std_len) :: filename
1794 
1795  write(filename, "(a,i4.4,a)") trim(base_filename), ix, ".dat"
1796  inquire(file=trim(filename), exist=snapshot_exists)
1797  end function snapshot_exists
1798 
1799  integer function get_snapshot_index(filename)
1800  character(len=*), intent(in) :: filename
1801  integer :: i
1802 
1803  ! Try to parse index in restart_from_file string (e.g. basename0000.dat)
1804  i = len_trim(filename) - 7
1805  read(filename(i:i+3), '(I4)') get_snapshot_index
1806  end function get_snapshot_index
1807 
1808 
1809  !> Write header for a snapshot, generalize cons_wnames and nw
1810  !>
1811  !> If you edit the header, don't forget to update: snapshot_write_header(),
1812  !> snapshot_read_header(), doc/fileformat.md, tools/python/dat_reader.py
1813  subroutine snapshot_write_header1(fh, offset_tree, offset_block, dataset_names, nw_vars)
1814  use mod_forest
1815  use mod_physics
1817  use mod_slice, only: slicenext
1818  integer, intent(in) :: fh !< File handle
1819  integer(kind=MPI_OFFSET_KIND), intent(in) :: offset_tree !< Offset of tree info
1820  integer(kind=MPI_OFFSET_KIND), intent(in) :: offset_block !< Offset of block data
1821  character(len=*), intent(in) :: dataset_names(:)
1822  integer, intent(in) :: nw_vars
1823  integer, dimension(MPI_STATUS_SIZE) :: st
1824  integer :: iw, er
1825 
1826  character(len=name_len) :: dname
1827 
1828  call mpi_file_write(fh, version_number, 1, mpi_integer, st, er)
1829  call mpi_file_write(fh, int(offset_tree), 1, mpi_integer, st, er)
1830  call mpi_file_write(fh, int(offset_block), 1, mpi_integer, st, er)
1831  call mpi_file_write(fh, nw_vars, 1, mpi_integer, st, er)
1832  call mpi_file_write(fh, ndir, 1, mpi_integer, st, er)
1833  call mpi_file_write(fh, ndim, 1, mpi_integer, st, er)
1834  call mpi_file_write(fh, levmax, 1, mpi_integer, st, er)
1835  call mpi_file_write(fh, nleafs, 1, mpi_integer, st, er)
1836  call mpi_file_write(fh, nparents, 1, mpi_integer, st, er)
1837  call mpi_file_write(fh, it, 1, mpi_integer, st, er)
1838  ! Note: It is nice when this double has an even number of 4 byte
1839  ! integers before it (for alignment)
1840  call mpi_file_write(fh, global_time, 1, mpi_double_precision, st, er)
1841 
1842  call mpi_file_write(fh, [ xprobmin^d ], ndim, &
1843  mpi_double_precision, st, er)
1844  call mpi_file_write(fh, [ xprobmax^d ], ndim, &
1845  mpi_double_precision, st, er)
1846  call mpi_file_write(fh, [ domain_nx^d ], ndim, &
1847  mpi_integer, st, er)
1848  call mpi_file_write(fh, [ block_nx^d ], ndim, &
1849  mpi_integer, st, er)
1850 
1851  ! Periodicity (assume all variables are periodic if one is)
1852  call mpi_file_write(fh, periodb, ndim, mpi_logical, st, er)
1853 
1854  ! Geometry
1855  call mpi_file_write(fh, geometry_name(1:name_len), &
1856  name_len, mpi_character, st, er)
1857 
1858  ! Write stagger grid mark
1859  call mpi_file_write(fh, stagger_grid, 1, mpi_logical, st, er)
1860 
1861  do iw = 1, nw_vars
1862  ! using directly trim(adjustl((dataset_names(iw)))) in MPI_FILE_WRITE call
1863  ! does not work, there will be trailing characters
1864  dname = trim(adjustl((dataset_names(iw))))
1865  call mpi_file_write(fh, dname, name_len, mpi_character, st, er)
1866  end do
1867 
1868  ! Physics related information
1869  call mpi_file_write(fh, physics_type, name_len, mpi_character, st, er)
1870 
1871  ! Format:
1872  ! integer :: n_par
1873  ! double precision :: values(n_par)
1874  ! character(n_par * name_len) :: names
1875  call phys_write_info(fh)
1876 
1877  ! Write snapshotnext etc., which is useful for restarting.
1878  ! Note we add one, since snapshotnext is updated *after* this procedure
1879  if(pass_wall_time.or.save_now) then
1880  call mpi_file_write(fh, snapshotnext, 1, mpi_integer, st, er)
1881  else
1882  call mpi_file_write(fh, snapshotnext+1, 1, mpi_integer, st, er)
1883  end if
1884  call mpi_file_write(fh, slicenext, 1, mpi_integer, st, er)
1885  call mpi_file_write(fh, collapsenext, 1, mpi_integer, st, er)
1886 
1887  end subroutine snapshot_write_header1
1888 
1889  !> Write header for a snapshot
1890  !>
1891  !> If you edit the header, don't forget to update: snapshot_write_header(),
1892  !> snapshot_read_header(), doc/fileformat.md, tools/python/dat_reader.py
1893  subroutine snapshot_write_header(fh, offset_tree, offset_block)
1894  use mod_forest
1895  use mod_physics
1897  use mod_slice, only: slicenext
1898  integer, intent(in) :: fh !< File handle
1899  integer(kind=MPI_OFFSET_KIND), intent(in) :: offset_tree !< Offset of tree info
1900  integer(kind=MPI_OFFSET_KIND), intent(in) :: offset_block !< Offset of block data
1901  call snapshot_write_header1(fh, offset_tree, offset_block, cons_wnames, nw)
1902  end subroutine snapshot_write_header
1903 
1904  !> Read header for a snapshot
1905  !>
1906  !> If you edit the header, don't forget to update: snapshot_write_header(),
1907  !> snapshot_read_header(), doc/fileformat.md, tools/python/dat_reader.py
1908  subroutine snapshot_read_header(fh, offset_tree, offset_block)
1909  use mod_forest
1911  use mod_physics, only: physics_type
1912  use mod_slice, only: slicenext
1913  integer, intent(in) :: fh !< File handle
1914  integer(MPI_OFFSET_KIND), intent(out) :: offset_tree !< Offset of tree info
1915  integer(MPI_OFFSET_KIND), intent(out) :: offset_block !< Offset of block data
1916  integer :: i, version
1917  integer :: ibuf(ndim), iw
1918  double precision :: rbuf(ndim)
1919  integer, dimension(MPI_STATUS_SIZE) :: st
1920  character(len=name_len), allocatable :: var_names(:), param_names(:)
1921  double precision, allocatable :: params(:)
1922  character(len=name_len) :: phys_name, geom_name
1923  integer :: er, n_par, tmp_int
1924  logical :: stagger_mark_dat, periodic(ndim)
1925 
1926  ! Version number
1927  call mpi_file_read(fh, version, 1, mpi_integer, st, er)
1928  if (all(compatible_versions /= version)) then
1929  call mpistop("Incompatible file version (maybe old format?)")
1930  end if
1931 
1932  ! offset_tree
1933  call mpi_file_read(fh, ibuf(1), 1, mpi_integer, st, er)
1934  offset_tree = ibuf(1)
1935 
1936  ! offset_block
1937  call mpi_file_read(fh, ibuf(1), 1, mpi_integer, st, er)
1938  offset_block = ibuf(1)
1939 
1940  ! nw
1941  call mpi_file_read(fh, ibuf(1), 1, mpi_integer, st, er)
1942  nw_found=ibuf(1)
1943  if (nw /= ibuf(1)) then
1944  write(*,*) "nw=",nw," and nw found in restart file=",ibuf(1)
1945  write(*,*) "Please be aware of changes in w at restart."
1946  !call mpistop("currently, changing nw at restart is not allowed")
1947  end if
1948 
1949  ! ndir
1950  call mpi_file_read(fh, ibuf(1), 1, mpi_integer, st, er)
1951  if (ibuf(1) /= ndir) then
1952  write(*,*) "ndir in restart file = ",ibuf(1)
1953  write(*,*) "ndir = ",ndir
1954  call mpistop("reset ndir to ndir in restart file")
1955  end if
1956 
1957  ! ndim
1958  call mpi_file_read(fh, ibuf(1), 1, mpi_integer, st, er)
1959  if (ibuf(1) /= ndim) then
1960  write(*,*) "ndim in restart file = ",ibuf(1)
1961  write(*,*) "ndim = ",ndim
1962  call mpistop("reset ndim to ndim in restart file")
1963  end if
1964 
1965  ! levmax
1966  call mpi_file_read(fh, ibuf(1), 1, mpi_integer, st, er)
1967  if (ibuf(1) > refine_max_level) then
1968  write(*,*) "number of levels in restart file = ",ibuf(1)
1969  write(*,*) "refine_max_level = ",refine_max_level
1970  call mpistop("refine_max_level < num. levels in restart file")
1971  end if
1972 
1973  ! nleafs
1974  call mpi_file_read(fh, nleafs, 1, mpi_integer, st, er)
1975 
1976  ! nparents
1977  call mpi_file_read(fh, nparents, 1, mpi_integer, st, er)
1978 
1979  ! it
1980  call mpi_file_read(fh, it, 1, mpi_integer, st, er)
1981 
1982  ! global time
1983  call mpi_file_read(fh, global_time, 1, mpi_double_precision, st, er)
1984 
1985  ! xprobmin^D
1986  call mpi_file_read(fh,rbuf(1:ndim),ndim,mpi_double_precision,st,er)
1987  if (maxval(abs(rbuf(1:ndim) - [ xprobmin^d ])) > 0) then
1988  write(*,*) "Error: xprobmin differs from restart data: ", rbuf(1:ndim)
1989  call mpistop("change xprobmin^D in par file")
1990  end if
1991 
1992  ! xprobmax^D
1993  call mpi_file_read(fh,rbuf(1:ndim),ndim,mpi_double_precision,st,er)
1994  if (maxval(abs(rbuf(1:ndim) - [ xprobmax^d ])) > 0) then
1995  write(*,*) "Error: xprobmax differs from restart data: ", rbuf(1:ndim)
1996  call mpistop("change xprobmax^D in par file")
1997  end if
1998 
1999  ! domain_nx^D
2000  call mpi_file_read(fh,ibuf(1:ndim), ndim, mpi_integer,st,er)
2001  if (any(ibuf(1:ndim) /= [ domain_nx^d ])) then
2002  write(*,*) "Error: mesh size differs from restart data: ", ibuf(1:ndim)
2003  call mpistop("change domain_nx^D in par file")
2004  end if
2005 
2006  ! block_nx^D
2007  call mpi_file_read(fh,ibuf(1:ndim), ndim, mpi_integer,st,er)
2008  if (any(ibuf(1:ndim) /= [ block_nx^d ])) then
2009  write(*,*) "Error: block size differs from restart data:", ibuf(1:ndim)
2010  call mpistop("change block_nx^D in par file")
2011  end if
2012 
2013  ! From version 5, read more info about the grid
2014  if (version > 4) then
2015  call mpi_file_read(fh, periodic, ndim, mpi_logical, st, er)
2016  if ({periodic(^d) .and. .not.periodb(^d) .or. .not.periodic(^d) .and. periodb(^d)| .or. }) &
2017  call mpistop("change in periodicity in par file")
2018 
2019  call mpi_file_read(fh, geom_name, name_len, mpi_character, st, er)
2020 
2021  if (geom_name /= geometry_name(1:name_len)) then
2022  write(*,*) "type of coordinates in data is: ", geom_name
2023  call mpistop("select the correct coordinates in mod_usr.t file")
2024  end if
2025 
2026  call mpi_file_read(fh, stagger_mark_dat, 1, mpi_logical, st, er)
2027  if (stagger_grid .and. .not. stagger_mark_dat .or. .not.stagger_grid.and.stagger_mark_dat) then
2028  write(*,*) "Error: stagger grid flag differs from restart data:", stagger_mark_dat
2029  call mpistop("change parameter to use stagger grid")
2030  end if
2031  end if
2032 
2033  ! From version 4 onwards, the later parts of the header must be present
2034  if (version > 3) then
2035  ! w_names (not used here)
2036  allocate(var_names(nw_found))
2037  do iw = 1, nw_found
2038  call mpi_file_read(fh, var_names(iw), name_len, mpi_character, st, er)
2039  end do
2040 
2041  ! Physics related information
2042  call mpi_file_read(fh, phys_name, name_len, mpi_character, st, er)
2043 
2044  if (phys_name /= physics_type) then
2045 ! call mpistop("Cannot restart with a different physics type")
2046  end if
2047 
2048  call mpi_file_read(fh, n_par, 1, mpi_integer, st, er)
2049  allocate(params(n_par))
2050  allocate(param_names(n_par))
2051  call mpi_file_read(fh, params, n_par, mpi_double_precision, st, er)
2052  call mpi_file_read(fh, param_names, name_len * n_par, mpi_character, st, er)
2053 
2054  ! Read snapshotnext etc. for restarting
2055  call mpi_file_read(fh, tmp_int, 1, mpi_integer, st, er)
2056 
2057  ! Only set snapshotnext if the user hasn't specified it
2058  if (snapshotnext == -1) snapshotnext = tmp_int
2059 
2060  call mpi_file_read(fh, tmp_int, 1, mpi_integer, st, er)
2061  if (slicenext == -1) slicenext = tmp_int
2062 
2063  call mpi_file_read(fh, tmp_int, 1, mpi_integer, st, er)
2064  if (collapsenext == -1) collapsenext = tmp_int
2065  else
2066  ! Guess snapshotnext from file name if not set
2067  if (snapshotnext == -1) &
2069  ! Set slicenext and collapsenext if not set
2070  if (slicenext == -1) slicenext = 0
2071  if (collapsenext == -1) collapsenext = 0
2072  end if
2073 
2074  ! Still used in convert
2076 
2077  end subroutine snapshot_read_header
2078 
2079  !> Compute number of elements in index range
2080  pure integer function count_ix(ixO^L)
2081  integer, intent(in) :: ixo^l
2082 
2083  count_ix = product([ ixomax^d ] - [ ixomin^d ] + 1)
2084  end function count_ix
2085 
2086  !> Determine the shape of a block for output (whether to include ghost cells,
2087  !> and on which sides)
2088  subroutine block_shape_io(igrid, n_ghost, ixO^L, n_values)
2090 
2091  integer, intent(in) :: igrid
2092  !> nghost(1:ndim) contains the number of ghost cells on the block's minimum
2093  !> boundaries, and nghost(ndim+1:2*ndim) on the block's maximum boundaries
2094  integer, intent(out) :: n_ghost(2*ndim)
2095  !> Index range on output block
2096  integer, intent(out) :: ixo^l
2097  !> Number of cells/values in output
2098  integer, intent(out) :: n_values
2099 
2100  integer :: idim
2101 
2102  n_ghost(:) = 0
2103 
2104  if(save_physical_boundary) then
2105  do idim=1,ndim
2106  ! Include ghost cells on lower boundary
2107  if(ps(igrid)%is_physical_boundary(2*idim-1)) n_ghost(idim)=nghostcells
2108  ! Include ghost cells on upper boundary
2109  if(ps(igrid)%is_physical_boundary(2*idim)) n_ghost(ndim+idim)=nghostcells
2110  end do
2111  end if
2112 
2113  {ixomin^d = ixmlo^d - n_ghost(^d)\}
2114  {ixomax^d = ixmhi^d + n_ghost(ndim+^d)\}
2115 
2116  n_values = count_ix(ixo^l) * nw
2117 
2118  end subroutine block_shape_io
2119 
2120  subroutine write_snapshot
2121  use mod_forest
2123  use mod_physics
2124 
2125  integer :: file_handle, igrid, Morton_no, iwrite
2126  integer :: ipe, ix_buffer(2*ndim+1), n_values
2127  integer :: ixO^L, n_ghost(2*ndim)
2128  integer :: ixOs^L,n_values_stagger
2129  integer :: iorecvstatus(MPI_STATUS_SIZE)
2130  integer :: ioastatus(MPI_STATUS_SIZE)
2131  integer :: igrecvstatus(MPI_STATUS_SIZE)
2132  integer :: istatus(MPI_STATUS_SIZE)
2133  type(tree_node), pointer :: pnode
2134  integer(kind=MPI_OFFSET_KIND) :: offset_tree_info
2135  integer(kind=MPI_OFFSET_KIND) :: offset_block_data
2136  integer(kind=MPI_OFFSET_KIND) :: offset_offsets
2137  double precision, allocatable :: w_buffer(:)
2138 
2139  integer, allocatable :: block_ig(:, :)
2140  integer, allocatable :: block_lvl(:)
2141  integer(kind=MPI_OFFSET_KIND), allocatable :: block_offset(:)
2142 
2143  call mpi_barrier(icomm, ierrmpi)
2144 
2145  ! Allocate send/receive buffer
2146  n_values = count_ix(ixg^ll) * nw
2147  if(stagger_grid) then
2148  n_values = n_values + count_ix(ixgs^ll) * nws
2149  end if
2150  allocate(w_buffer(n_values))
2151 
2152  ! Allocate arrays with information about grid blocks
2153  allocate(block_ig(ndim, nleafs))
2154  allocate(block_lvl(nleafs))
2155  allocate(block_offset(nleafs+1))
2156 
2157  ! master processor
2158  if (mype==0) then
2159  call create_output_file(file_handle, snapshotnext, ".dat")
2160 
2161  ! Don't know offsets yet, we will write header again later
2162  offset_tree_info = -1
2163  offset_block_data = -1
2164  call snapshot_write_header(file_handle, offset_tree_info, &
2165  offset_block_data)
2166 
2167  call mpi_file_get_position(file_handle, offset_tree_info, ierrmpi)
2168 
2169  call write_forest(file_handle)
2170 
2171  ! Collect information about the spatial index (ig^D) and refinement level
2172  ! of leaves
2173  do morton_no = morton_start(0), morton_stop(npe-1)
2174  igrid = sfc(1, morton_no)
2175  ipe = sfc(2, morton_no)
2176  pnode => igrid_to_node(igrid, ipe)%node
2177 
2178  block_ig(:, morton_no) = [ pnode%ig^d ]
2179  block_lvl(morton_no) = pnode%level
2180  block_offset(morton_no) = 0 ! Will be determined later
2181  end do
2182 
2183  call mpi_file_write(file_handle, block_lvl, size(block_lvl), &
2184  mpi_integer, istatus, ierrmpi)
2185 
2186  call mpi_file_write(file_handle, block_ig, size(block_ig), &
2187  mpi_integer, istatus, ierrmpi)
2188 
2189  ! Block offsets are currently unknown, but will be overwritten later
2190  call mpi_file_get_position(file_handle, offset_offsets, ierrmpi)
2191  call mpi_file_write(file_handle, block_offset(1:nleafs), nleafs, &
2192  mpi_offset, istatus, ierrmpi)
2193 
2194  call mpi_file_get_position(file_handle, offset_block_data, ierrmpi)
2195 
2196  ! Check whether data was written as expected
2197  if (offset_block_data - offset_tree_info /= &
2198  (nleafs + nparents) * size_logical + &
2199  nleafs * ((1+ndim) * size_int + 2 * size_int)) then
2200  if (mype == 0) then
2201  print *, "Warning: MPI_OFFSET type /= 8 bytes"
2202  print *, "This *could* cause problems when reading .dat files"
2203  end if
2204  end if
2205 
2206  block_offset(1) = offset_block_data
2207  iwrite = 0
2208  end if
2209 
2210  do morton_no=morton_start(mype), morton_stop(mype)
2211  igrid = sfc_to_igrid(morton_no)
2212  itag = morton_no
2213 
2214  call block_shape_io(igrid, n_ghost, ixo^l, n_values)
2215  if(stagger_grid) then
2216  w_buffer(1:n_values) = pack(ps(igrid)%w(ixo^s, 1:nw), .true.)
2217  {ixosmin^d = ixomin^d -1\}
2218  {ixosmax^d = ixomax^d \}
2219  n_values_stagger= count_ix(ixos^l)*nws
2220  w_buffer(n_values+1:n_values+n_values_stagger) = pack(ps(igrid)%ws(ixos^s, 1:nws), .true.)
2221  n_values=n_values+n_values_stagger
2222  else
2223  w_buffer(1:n_values) = pack(ps(igrid)%w(ixo^s, 1:nw), .true.)
2224  end if
2225  ix_buffer(1) = n_values
2226  ix_buffer(2:) = n_ghost
2227 
2228  if (mype /= 0) then
2229  call mpi_send(ix_buffer, 2*ndim+1, &
2230  mpi_integer, 0, itag, icomm, ierrmpi)
2231  call mpi_send(w_buffer, n_values, &
2232  mpi_double_precision, 0, itag, icomm, ierrmpi)
2233  else
2234  iwrite = iwrite+1
2235  call mpi_file_write(file_handle, ix_buffer(2:), &
2236  2*ndim, mpi_integer, istatus, ierrmpi)
2237  call mpi_file_write(file_handle, w_buffer, &
2238  n_values, mpi_double_precision, istatus, ierrmpi)
2239 
2240  ! Set offset of next block
2241  block_offset(iwrite+1) = block_offset(iwrite) + &
2242  int(n_values, mpi_offset_kind) * size_double + &
2243  2 * ndim * size_int
2244  end if
2245  end do
2246 
2247  ! Write data communicated from other processors
2248  if (mype == 0) then
2249  do ipe = 1, npe-1
2250  do morton_no=morton_start(ipe), morton_stop(ipe)
2251  iwrite=iwrite+1
2252  itag=morton_no
2253 
2254  call mpi_recv(ix_buffer, 2*ndim+1, mpi_integer, ipe, itag, icomm,&
2255  igrecvstatus, ierrmpi)
2256  n_values = ix_buffer(1)
2257 
2258  call mpi_recv(w_buffer, n_values, mpi_double_precision,&
2259  ipe, itag, icomm, iorecvstatus, ierrmpi)
2260 
2261  call mpi_file_write(file_handle, ix_buffer(2:), &
2262  2*ndim, mpi_integer, istatus, ierrmpi)
2263  call mpi_file_write(file_handle, w_buffer, &
2264  n_values, mpi_double_precision, istatus, ierrmpi)
2265 
2266  ! Set offset of next block
2267  block_offset(iwrite+1) = block_offset(iwrite) + &
2268  int(n_values, mpi_offset_kind) * size_double + &
2269  2 * ndim * size_int
2270  end do
2271  end do
2272 
2273  ! Write block offsets (now we know them)
2274  call mpi_file_seek(file_handle, offset_offsets, mpi_seek_set, ierrmpi)
2275  call mpi_file_write(file_handle, block_offset(1:nleafs), nleafs, &
2276  mpi_offset, istatus, ierrmpi)
2277 
2278  ! Write header again, now with correct offsets
2279  call mpi_file_seek(file_handle, 0_mpi_offset_kind, mpi_seek_set, ierrmpi)
2280  call snapshot_write_header(file_handle, offset_tree_info, &
2281  offset_block_data)
2282 
2283  call mpi_file_close(file_handle, ierrmpi)
2284  end if
2285 
2286  call mpi_barrier(icomm, ierrmpi)
2287 
2288  end subroutine write_snapshot
2289 
2290  !> Routine to read in snapshots (.dat files). When it cannot recognize the
2291  !> file version, it will automatically try the 'old' reader.
2292  subroutine read_snapshot
2293  use mod_usr_methods, only: usr_transform_w
2294  use mod_forest
2296  use mod_slice, only: slicenext
2297 
2298  integer :: ix_buffer(2*ndim+1), n_values, n_values_stagger
2299  integer :: ixO^L, ixOs^L
2300  integer :: file_handle, amode, igrid, Morton_no, iread
2301  integer :: istatus(MPI_STATUS_SIZE)
2302  integer :: iorecvstatus(MPI_STATUS_SIZE)
2303  integer :: ipe,inrecv,nrecv, file_version
2304  logical :: fexist
2305  integer(MPI_OFFSET_KIND) :: offset_tree_info
2306  integer(MPI_OFFSET_KIND) :: offset_block_data
2307  double precision, allocatable :: w_buffer(:)
2308  double precision, dimension(:^D&,:), allocatable :: w
2309  double precision :: ws(ixGs^T,1:ndim)
2310 
2311  if (mype==0) then
2312  inquire(file=trim(restart_from_file), exist=fexist)
2313  if (.not.fexist) call mpistop(trim(restart_from_file)//" not found!")
2314 
2315  call mpi_file_open(mpi_comm_self,restart_from_file,mpi_mode_rdonly, &
2316  mpi_info_null,file_handle,ierrmpi)
2317  call mpi_file_read(file_handle, file_version, 1, mpi_integer, &
2318  istatus, ierrmpi)
2319  end if
2320 
2321  call mpi_bcast(file_version,1,mpi_integer,0,icomm,ierrmpi)
2322 
2323  if (all(compatible_versions /= file_version)) then
2324  if (mype == 0) print *, "Unknown version, trying old snapshot reader..."
2325  call mpi_file_close(file_handle,ierrmpi)
2326  call read_snapshot_old()
2327 
2328  ! Guess snapshotnext from file name if not set
2329  if (snapshotnext == -1) &
2331  ! Set slicenext and collapsenext if not set
2332  if (slicenext == -1) slicenext = 0
2333  if (collapsenext == -1) collapsenext = 0
2334 
2335  ! Still used in convert
2337 
2338  return ! Leave this routine
2339  else if (mype == 0) then
2340  call mpi_file_seek(file_handle, 0_mpi_offset_kind, mpi_seek_set, ierrmpi)
2341  call snapshot_read_header(file_handle, offset_tree_info, &
2342  offset_block_data)
2343  end if
2344 
2345  ! Share information about restart file
2346  call mpi_bcast(nw_found,1,mpi_integer,0,icomm,ierrmpi)
2347  call mpi_bcast(nleafs,1,mpi_integer,0,icomm,ierrmpi)
2348  call mpi_bcast(nparents,1,mpi_integer,0,icomm,ierrmpi)
2349  call mpi_bcast(it,1,mpi_integer,0,icomm,ierrmpi)
2350  call mpi_bcast(global_time,1,mpi_double_precision,0,icomm,ierrmpi)
2351 
2352  call mpi_bcast(snapshotnext,1,mpi_integer,0,icomm,ierrmpi)
2353  call mpi_bcast(slicenext,1,mpi_integer,0,icomm,ierrmpi)
2354  call mpi_bcast(collapsenext,1,mpi_integer,0,icomm,ierrmpi)
2355 
2356  ! Allocate send/receive buffer
2357  n_values = count_ix(ixg^ll) * nw_found
2358  if(stagger_grid) then
2359  n_values = n_values + count_ix(ixgs^ll) * nws
2360  end if
2361  allocate(w_buffer(n_values))
2362  allocate(w(ixg^t,1:nw_found))
2363 
2365 
2366  if (mype == 0) then
2367  call mpi_file_seek(file_handle, offset_tree_info, &
2368  mpi_seek_set, ierrmpi)
2369  end if
2370 
2371  call read_forest(file_handle)
2372 
2373  do morton_no=morton_start(mype),morton_stop(mype)
2374  igrid=sfc_to_igrid(morton_no)
2375  call alloc_node(igrid)
2376  end do
2377 
2378  if (mype==0) then
2379  call mpi_file_seek(file_handle, offset_block_data, mpi_seek_set, ierrmpi)
2380 
2381  iread = 0
2382  do ipe = 0, npe-1
2383  do morton_no=morton_start(ipe),morton_stop(ipe)
2384  iread=iread+1
2385  itag=morton_no
2386 
2387  call mpi_file_read(file_handle,ix_buffer(1:2*ndim), 2*ndim, &
2388  mpi_integer, istatus,ierrmpi)
2389 
2390  ! Construct ixO^L array from number of ghost cells
2391  {ixomin^d = ixmlo^d - ix_buffer(^d)\}
2392  {ixomax^d = ixmhi^d + ix_buffer(ndim+^d)\}
2393  n_values = count_ix(ixo^l) * nw_found
2394  if(stagger_grid) then
2395  {ixosmin^d = ixomin^d - 1\}
2396  {ixosmax^d = ixomax^d\}
2397  n_values_stagger = n_values
2398  n_values = n_values + count_ix(ixos^l) * nws
2399  end if
2400 
2401  call mpi_file_read(file_handle, w_buffer, n_values, &
2402  mpi_double_precision, istatus, ierrmpi)
2403 
2404  if (mype == ipe) then ! Root task
2405  igrid=sfc_to_igrid(morton_no)
2406  block=>ps(igrid)
2407  if(stagger_grid) then
2408  w(ixo^s, 1:nw_found) = reshape(w_buffer(1:n_values_stagger), &
2409  shape(w(ixo^s, 1:nw_found)))
2410  ps(igrid)%ws(ixos^s,1:nws)=reshape(w_buffer(n_values_stagger+1:n_values), &
2411  shape(ws(ixos^s, 1:nws)))
2412  else
2413  w(ixo^s, 1:nw_found) = reshape(w_buffer(1:n_values), &
2414  shape(w(ixo^s, 1:nw_found)))
2415  end if
2416  if (nw_found<nw) then
2417  if (associated(usr_transform_w)) then
2418  call usr_transform_w(ixg^ll,ixm^ll,nw_found,w,ps(igrid)%x,ps(igrid)%w)
2419  else
2420  ps(igrid)%w(ixo^s,1:nw_found)=w(ixo^s,1:nw_found)
2421  end if
2422  else if (nw_found>nw) then
2423  if (associated(usr_transform_w)) then
2424  call usr_transform_w(ixg^ll,ixm^ll,nw_found,w,ps(igrid)%x,ps(igrid)%w)
2425  else
2426  ps(igrid)%w(ixo^s,1:nw)=w(ixo^s,1:nw)
2427  end if
2428  else
2429  ps(igrid)%w(ixo^s,1:nw)=w(ixo^s,1:nw)
2430  end if
2431  else
2432  call mpi_send([ ixo^l, n_values ], 2*ndim+1, &
2433  mpi_integer, ipe, itag, icomm, ierrmpi)
2434  call mpi_send(w_buffer, n_values, &
2435  mpi_double_precision, ipe, itag, icomm, ierrmpi)
2436  end if
2437  end do
2438  end do
2439 
2440  call mpi_file_close(file_handle,ierrmpi)
2441 
2442  else ! mype > 0
2443 
2444  do morton_no=morton_start(mype),morton_stop(mype)
2445  igrid=sfc_to_igrid(morton_no)
2446  block=>ps(igrid)
2447  itag=morton_no
2448 
2449  call mpi_recv(ix_buffer, 2*ndim+1, mpi_integer, 0, itag, icomm,&
2450  iorecvstatus, ierrmpi)
2451  {ixomin^d = ix_buffer(^d)\}
2452  {ixomax^d = ix_buffer(ndim+^d)\}
2453  n_values = ix_buffer(2*ndim+1)
2454 
2455  call mpi_recv(w_buffer, n_values, mpi_double_precision,&
2456  0, itag, icomm, iorecvstatus, ierrmpi)
2457 
2458  if(stagger_grid) then
2459  n_values_stagger = count_ix(ixo^l) * nw_found
2460  {ixosmin^d = ixomin^d - 1\}
2461  {ixosmax^d = ixomax^d\}
2462  w(ixo^s, 1:nw_found) = reshape(w_buffer(1:n_values_stagger), &
2463  shape(w(ixo^s, 1:nw_found)))
2464  ps(igrid)%ws(ixos^s,1:nws)=reshape(w_buffer(n_values_stagger+1:n_values), &
2465  shape(ws(ixos^s, 1:nws)))
2466  else
2467  w(ixo^s, 1:nw_found) = reshape(w_buffer(1:n_values), &
2468  shape(w(ixo^s, 1:nw_found)))
2469  end if
2470  if (nw_found<nw) then
2471  if (associated(usr_transform_w)) then
2472  call usr_transform_w(ixg^ll,ixm^ll,nw_found,w,ps(igrid)%x,ps(igrid)%w)
2473  else
2474  ps(igrid)%w(ixo^s,1:nw_found)=w(ixo^s,1:nw_found)
2475  end if
2476  else if (nw_found>nw) then
2477  if (associated(usr_transform_w)) then
2478  call usr_transform_w(ixg^ll,ixm^ll,nw_found,w,ps(igrid)%x,ps(igrid)%w)
2479  else
2480  ps(igrid)%w(ixo^s,1:nw)=w(ixo^s,1:nw)
2481  end if
2482  else
2483  ps(igrid)%w(ixo^s,1:nw)=w(ixo^s,1:nw)
2484  end if
2485  end do
2486  end if
2487 
2488  call mpi_barrier(icomm,ierrmpi)
2489 
2490  end subroutine read_snapshot
2491 
2492  subroutine read_snapshot_old()
2493  use mod_forest
2495 
2496  double precision :: wio(ixG^T,1:nw)
2497  integer :: fh, igrid, Morton_no, iread
2498  integer :: levmaxini, ndimini, ndirini
2499  integer :: nwini, neqparini, nxini^D
2500  integer(kind=MPI_OFFSET_KIND) :: offset
2501  integer :: istatus(MPI_STATUS_SIZE)
2502  integer, allocatable :: iorecvstatus(:,:)
2503  integer :: ipe,inrecv,nrecv
2504  integer :: sendini(7+^ND)
2505  character(len=80) :: filename
2506  logical :: fexist
2507  double precision :: eqpar_dummy(100)
2508 
2509  if (mype==0) then
2510  call mpi_file_open(mpi_comm_self,trim(restart_from_file), &
2511  mpi_mode_rdonly,mpi_info_null,fh,ierrmpi)
2512 
2513  offset=-int(7*size_int+size_double,kind=mpi_offset_kind)
2514  call mpi_file_seek(fh,offset,mpi_seek_end,ierrmpi)
2515 
2516  call mpi_file_read(fh,nleafs,1,mpi_integer,istatus,ierrmpi)
2518  call mpi_file_read(fh,levmaxini,1,mpi_integer,istatus,ierrmpi)
2519  call mpi_file_read(fh,ndimini,1,mpi_integer,istatus,ierrmpi)
2520  call mpi_file_read(fh,ndirini,1,mpi_integer,istatus,ierrmpi)
2521  call mpi_file_read(fh,nwini,1,mpi_integer,istatus,ierrmpi)
2522  call mpi_file_read(fh,neqparini,1,mpi_integer,istatus,ierrmpi)
2523  call mpi_file_read(fh,it,1,mpi_integer,istatus,ierrmpi)
2524  call mpi_file_read(fh,global_time,1,mpi_double_precision,istatus,ierrmpi)
2525 
2526  ! check if settings are suitable for restart
2527  if (levmaxini>refine_max_level) then
2528  write(*,*) "number of levels in restart file = ",levmaxini
2529  write(*,*) "refine_max_level = ",refine_max_level
2530  call mpistop("refine_max_level < number of levels in restart file")
2531  end if
2532  if (ndimini/=ndim) then
2533  write(*,*) "ndim in restart file = ",ndimini
2534  write(*,*) "ndim = ",ndim
2535  call mpistop("reset ndim to ndim in restart file")
2536  end if
2537  if (ndirini/=ndir) then
2538  write(*,*) "ndir in restart file = ",ndirini
2539  write(*,*) "ndir = ",ndir
2540  call mpistop("reset ndir to ndir in restart file")
2541  end if
2542  if (nw/=nwini) then
2543  write(*,*) "nw=",nw," and nw in restart file=",nwini
2544  call mpistop("currently, changing nw at restart is not allowed")
2545  end if
2546 
2547  offset=offset-int(ndimini*size_int+neqparini*size_double,kind=mpi_offset_kind)
2548  call mpi_file_seek(fh,offset,mpi_seek_end,ierrmpi)
2549 
2550  {call mpi_file_read(fh,nxini^d,1,mpi_integer,istatus,ierrmpi)\}
2551  if (ixghi^d/=nxini^d+2*nghostcells|.or.) then
2552  write(*,*) "Error: reset resolution to ",nxini^d+2*nghostcells
2553  call mpistop("change with setamrvac")
2554  end if
2555 
2556  call mpi_file_read(fh,eqpar_dummy,neqparini, &
2557  mpi_double_precision,istatus,ierrmpi)
2558  end if
2559 
2560  ! broadcast the global parameters first
2561  if (npe>1) then
2562  if (mype==0) then
2563  sendini=(/nleafs,levmaxini,ndimini,ndirini,nwini,neqparini,it ,^d&nxini^d /)
2564  end if
2565  call mpi_bcast(sendini,7+^nd,mpi_integer,0,icomm,ierrmpi)
2566  nleafs=sendini(1);levmaxini=sendini(2);ndimini=sendini(3);
2567  ndirini=sendini(4);nwini=sendini(5);
2568  neqparini=sendini(6);it=sendini(7);
2569  nxini^d=sendini(7+^d);
2571  call mpi_bcast(global_time,1,mpi_double_precision,0,icomm,ierrmpi)
2572  end if
2573 
2574  if (mype == 0) then
2575  offset = int(size_block_io,kind=mpi_offset_kind) * &
2576  int(nleafs,kind=mpi_offset_kind)
2577  call mpi_file_seek(fh,offset,mpi_seek_set,ierrmpi)
2578  end if
2579 
2580  call read_forest(fh)
2581 
2582  do morton_no=morton_start(mype),morton_stop(mype)
2583  igrid=sfc_to_igrid(morton_no)
2584  call alloc_node(igrid)
2585  end do
2586 
2587  if (mype==0)then
2588  iread=0
2589 
2590  do morton_no=morton_start(0),morton_stop(0)
2591  igrid=sfc_to_igrid(morton_no)
2592  iread=iread+1
2593  offset=int(size_block_io,kind=mpi_offset_kind) &
2594  *int(morton_no-1,kind=mpi_offset_kind)
2595  call mpi_file_read_at(fh,offset,ps(igrid)%w,1,type_block_io, &
2596  istatus,ierrmpi)
2597  end do
2598  if (npe>1) then
2599  do ipe=1,npe-1
2600  do morton_no=morton_start(ipe),morton_stop(ipe)
2601  iread=iread+1
2602  itag=morton_no
2603  offset=int(size_block_io,kind=mpi_offset_kind)&
2604  *int(morton_no-1,kind=mpi_offset_kind)
2605  call mpi_file_read_at(fh,offset,wio,1,type_block_io,&
2606  istatus,ierrmpi)
2607  call mpi_send(wio,1,type_block_io,ipe,itag,icomm,ierrmpi)
2608  end do
2609  end do
2610  end if
2611  call mpi_file_close(fh,ierrmpi)
2612  else
2613  nrecv=(morton_stop(mype)-morton_start(mype)+1)
2614  allocate(iorecvstatus(mpi_status_size,nrecv))
2615  inrecv=0
2616  do morton_no=morton_start(mype),morton_stop(mype)
2617  igrid=sfc_to_igrid(morton_no)
2618  itag=morton_no
2619  inrecv=inrecv+1
2620  call mpi_recv(ps(igrid)%w,1,type_block_io,0,itag,icomm,&
2621  iorecvstatus(:,inrecv),ierrmpi)
2622  end do
2623  deallocate(iorecvstatus)
2624  end if
2625 
2626  call mpi_barrier(icomm,ierrmpi)
2627 
2628  end subroutine read_snapshot_old
2629 
2630  !> Write volume-averaged values and other information to the log file
2631  subroutine printlog_default
2632 
2633  use mod_timing
2636 
2637  logical :: fileopen
2638  integer :: i, iw, level
2639  double precision :: wmean(1:nw), total_volume
2640  double precision :: volume_coverage(refine_max_level)
2641  integer :: nx^D, nc, ncells, dit
2642  double precision :: dtTimeLast, now, cellupdatesPerSecond
2643  double precision :: activeBlocksPerCore, wctPerCodeTime, timeToFinish
2644  character(len=40) :: fmt_string
2645  character(len=80) :: filename
2646  character(len=2048) :: line
2647  logical, save :: opened = .false.
2648  integer :: amode, istatus(MPI_STATUS_SIZE)
2649  integer, parameter :: my_unit = 20
2650 
2651  ! Compute the volume-average of w**1 = w
2652  call get_volume_average(1, wmean, total_volume)
2653 
2654  ! Compute the volume coverage
2655  call get_volume_coverage(volume_coverage)
2656 
2657  if (mype == 0) then
2658 
2659  ! To compute cell updates per second, we do the following:
2660  nx^d=ixmhi^d-ixmlo^d+1;
2661  nc={nx^d*}
2662  ncells = nc * nleafs_active
2663 
2664  ! assumes the number of active leafs haven't changed since last compute.
2665  now = mpi_wtime()
2666  dit = it - ittimelast
2667  dttimelast = now - timelast
2668  ittimelast = it
2669  timelast = now
2670  cellupdatespersecond = dble(ncells) * dble(nstep) * &
2671  dble(dit) / (dttimelast * dble(npe))
2672 
2673  ! blocks per core:
2674  activeblockspercore = dble(nleafs_active) / dble(npe)
2675 
2676  ! Wall clock time per code time unit in seconds:
2677  wctpercodetime = dttimelast / max(dit * dt, epsilon(1.0d0))
2678 
2679  ! Wall clock time to finish in hours:
2680  timetofinish = (time_max - global_time) * wctpercodetime / 3600.0d0
2681 
2682  ! On first entry, open the file and generate the header
2683  if (.not. opened) then
2684 
2685  filename = trim(base_filename) // ".log"
2686 
2687  ! Delete the log when not doing a restart run
2688  if (restart_from_file == undefined) then
2689  open(unit=my_unit,file=trim(filename),status='replace')
2690  close(my_unit, status='delete')
2691  end if
2692 
2693  amode = ior(mpi_mode_create,mpi_mode_wronly)
2694  amode = ior(amode,mpi_mode_append)
2695 
2696  call mpi_file_open(mpi_comm_self, filename, amode, &
2697  mpi_info_null, log_fh, ierrmpi)
2698 
2699  opened = .true.
2700 
2701  ! Start of file headern
2702  line = "it global_time dt"
2703  do level=1,nw
2704  i = len_trim(line) + 2
2705  write(line(i:),"(a,a)") trim(cons_wnames(level)), " "
2706  end do
2707 
2708  ! Volume coverage per level
2709  do level = 1, refine_max_level
2710  i = len_trim(line) + 2
2711  write(line(i:), "(a,i0)") "c", level
2712  end do
2713 
2714  ! Cell counts per level
2715  do level=1,refine_max_level
2716  i = len_trim(line) + 2
2717  write(line(i:), "(a,i0)") "n", level
2718  end do
2719 
2720  ! Rest of file header
2721  line = trim(line) // " | Xload Xmemory 'Cell_Updates /second/core'"
2722  line = trim(line) // " 'Active_Blocks/Core' 'Wct Per Code Time [s]'"
2723  line = trim(line) // " 'TimeToFinish [hrs]'"
2724 
2725  ! Only write header if not restarting
2726  if (restart_from_file == undefined .or. reset_time) then
2727  call mpi_file_write(log_fh, trim(line) // new_line('a'), &
2728  len_trim(line)+1, mpi_character, istatus, ierrmpi)
2729  end if
2730  end if
2731 
2732  ! Construct the line to be added to the log
2733 
2734  fmt_string = '(' // fmt_i // ',2' // fmt_r // ')'
2735  write(line, fmt_string) it, global_time, dt
2736  i = len_trim(line) + 2
2737 
2738  write(fmt_string, '(a,i0,a)') '(', nw, fmt_r // ')'
2739  write(line(i:), fmt_string) wmean(1:nw)
2740  i = len_trim(line) + 2
2741 
2742  write(fmt_string, '(a,i0,a)') '(', refine_max_level, fmt_r // ')'
2743  write(line(i:), fmt_string) volume_coverage(1:refine_max_level)
2744  i = len_trim(line) + 2
2745 
2746  write(fmt_string, '(a,i0,a)') '(', refine_max_level, fmt_i // ')'
2747  write(line(i:), fmt_string) nleafs_level(1:refine_max_level)
2748  i = len_trim(line) + 2
2749 
2750  fmt_string = '(a,6' // fmt_r2 // ')'
2751  write(line(i:), fmt_string) '| ', xload, xmemory, cellupdatespersecond, &
2752  activeblockspercore, wctpercodetime, timetofinish
2753 
2754  call mpi_file_write(log_fh, trim(line) // new_line('a') , &
2755  len_trim(line)+1, mpi_character, istatus, ierrmpi)
2756  end if
2757 
2758  end subroutine printlog_default
2759 
2760  !> Print a log that can be used to check whether the code still produces the
2761  !> same output (regression test)
2764 
2765  integer, parameter :: n_modes = 2
2766  integer, parameter :: my_unit = 123
2767  character(len=40) :: fmt_string
2768  logical, save :: file_open = .false.
2769  integer :: power
2770  double precision :: modes(nw, n_modes), volume
2771 
2772  do power = 1, n_modes
2773  call get_volume_average(power, modes(:, power), volume)
2774  end do
2775 
2776  if (mype == 0) then
2777  if (.not. file_open) then
2778  open(my_unit, file = trim(base_filename) // ".log")
2779  file_open = .true.
2780 
2781  write(my_unit, *) "# time mean(w) mean(w**2)"
2782  end if
2783 
2784  write(fmt_string, "(a,i0,a)") "(", nw * n_modes + 1, fmt_r // ")"
2785  write(my_unit, fmt_string) global_time, modes
2786  end if
2787  end subroutine printlog_regression_test
2788 
2789  !> Compute mean(w**power) over the leaves of the grid. The first mode
2790  !> (power=1) corresponds to the mean, the second to the mean squared values
2791  !> and so on.
2792  subroutine get_volume_average(power, mode, volume)
2794 
2795  integer, intent(in) :: power !< Which mode to compute
2796  double precision, intent(out) :: mode(nw) !< The computed mode
2797  double precision, intent(out) :: volume !< The total grid volume
2798  integer :: iigrid, igrid, iw
2799  double precision :: wsum(nw+1)
2800  double precision :: dsum_recv(1:nw+1)
2801 
2802  wsum(:) = 0
2803 
2804  ! Loop over all the grids
2805  do iigrid = 1, igridstail
2806  igrid = igrids(iigrid)
2807 
2808  ! Store total volume in last element
2809  wsum(nw+1) = wsum(nw+1) + sum(ps(igrid)%dvolume(ixm^t))
2810 
2811  ! Compute the modes of the cell-centered variables, weighted by volume
2812  do iw = 1, nw
2813  wsum(iw) = wsum(iw) + &
2814  sum(ps(igrid)%dvolume(ixm^t)*ps(igrid)%w(ixm^t,iw)**power)
2815  end do
2816  end do
2817 
2818  ! Make the information available on all tasks
2819  call mpi_allreduce(wsum, dsum_recv, nw+1, mpi_double_precision, &
2820  mpi_sum, icomm, ierrmpi)
2821 
2822  ! Set the volume and the average
2823  volume = dsum_recv(nw+1)
2824  mode = dsum_recv(1:nw) / volume
2825 
2826  end subroutine get_volume_average
2827 
2828  !> Compute how much of the domain is covered by each grid level. This routine
2829  !> does not take a non-Cartesian geometry into account.
2830  subroutine get_volume_coverage(vol_cov)
2832 
2833  double precision, intent(out) :: vol_cov(1:refine_max_level)
2834  double precision :: dsum_recv(1:refine_max_level)
2835  integer :: iigrid, igrid, iw, level
2836 
2837  ! First determine the total 'flat' volume in each level
2838  vol_cov(1:refine_max_level)=zero
2839 
2840  do iigrid = 1, igridstail
2841  igrid = igrids(iigrid);
2842  level = node(plevel_,igrid)
2843  vol_cov(level) = vol_cov(level)+ &
2844  {(rnode(rpxmax^d_,igrid)-rnode(rpxmin^d_,igrid))|*}
2845  end do
2846 
2847  ! Make the information available on all tasks
2848  call mpi_allreduce(vol_cov, dsum_recv, refine_max_level, mpi_double_precision, &
2849  mpi_sum, icomm, ierrmpi)
2850 
2851  ! Normalize
2852  vol_cov = dsum_recv / sum(dsum_recv)
2853  end subroutine get_volume_coverage
2854 
2855  !> Compute the volume average of func(w) over the leaves of the grid.
2856  subroutine get_volume_average_func(func, f_avg, volume)
2858 
2859  interface
2860  pure function func(w_vec, w_size) result(val)
2861  integer, intent(in) :: w_size
2862  double precision, intent(in) :: w_vec(w_size)
2863  double precision :: val
2864  end function func
2865  end interface
2866  double precision, intent(out) :: f_avg !< The volume average of func
2867  double precision, intent(out) :: volume !< The total grid volume
2868  integer :: iigrid, igrid, i^D
2869  double precision :: wsum(2)
2870  double precision :: dsum_recv(2)
2871 
2872  wsum(:) = 0
2873 
2874  ! Loop over all the grids
2875  do iigrid = 1, igridstail
2876  igrid = igrids(iigrid)
2877 
2878  ! Store total volume in last element
2879  wsum(2) = wsum(2) + sum(ps(igrid)%dvolume(ixm^t))
2880 
2881  ! Compute the modes of the cell-centered variables, weighted by volume
2882  {do i^d = ixmlo^d, ixmhi^d\}
2883  wsum(1) = wsum(1) + ps(igrid)%dvolume(i^d) * &
2884  func(ps(igrid)%w(i^d, :), nw)
2885  {end do\}
2886  end do
2887 
2888  ! Make the information available on all tasks
2889  call mpi_allreduce(wsum, dsum_recv, 2, mpi_double_precision, &
2890  mpi_sum, icomm, ierrmpi)
2891 
2892  ! Set the volume and the average
2893  volume = dsum_recv(2)
2894  f_avg = dsum_recv(1) / volume
2895 
2896  end subroutine get_volume_average_func
2897 
2898  !> Compute global maxima of iw variables over the leaves of the grid.
2899  subroutine get_global_maxima(wmax)
2901 
2902  double precision, intent(out) :: wmax(nw) !< The global maxima
2903 
2904  integer :: iigrid, igrid, iw
2905  double precision :: wmax_mype(nw),wmax_recv(nw)
2906 
2907  wmax_mype(1:nw) = -bigdouble
2908 
2909  ! Loop over all the grids
2910  do iigrid = 1, igridstail
2911  igrid = igrids(iigrid)
2912  do iw = 1, nw
2913  wmax_mype(iw)=max(wmax_mype(iw),maxval(ps(igrid)%w(ixm^t,iw)))
2914  end do
2915  end do
2916 
2917  ! Make the information available on all tasks
2918  call mpi_allreduce(wmax_mype, wmax_recv, nw, mpi_double_precision, &
2919  mpi_max, icomm, ierrmpi)
2920 
2921  wmax(1:nw)=wmax_recv(1:nw)
2922 
2923  end subroutine get_global_maxima
2924 
2925  !> Compute global minima of iw variables over the leaves of the grid.
2926  subroutine get_global_minima(wmin)
2928 
2929  double precision, intent(out) :: wmin(nw) !< The global maxima
2930 
2931  integer :: iigrid, igrid, iw
2932  double precision :: wmin_mype(nw),wmin_recv(nw)
2933 
2934  wmin_mype(1:nw) = bigdouble
2935 
2936  ! Loop over all the grids
2937  do iigrid = 1, igridstail
2938  igrid = igrids(iigrid)
2939  do iw = 1, nw
2940  wmin_mype(iw)=min(wmin_mype(iw),minval(ps(igrid)%w(ixm^t,iw)))
2941  end do
2942  end do
2943 
2944  ! Make the information available on all tasks
2945  call mpi_allreduce(wmin_mype, wmin_recv, nw, mpi_double_precision, &
2946  mpi_min, icomm, ierrmpi)
2947 
2948  wmin(1:nw)=wmin_recv(1:nw)
2949 
2950  end subroutine get_global_minima
2951 
2952 
2953 
2954 end module mod_input_output
subroutine alloc_node(igrid)
allocate arrays on igrid node
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:194
subroutine generate_plotfile
Definition: convert.t:3
subroutine read_forest(file_handle)
Definition: forest.t:318
subroutine write_forest(file_handle)
Definition: forest.t:283
Collapses D-dimensional output to D-1 view by line-of-sight integration.
Definition: mod_collapse.t:4
subroutine write_collapsed
Definition: mod_collapse.t:11
Module with basic grid data structures.
Definition: mod_forest.t:2
integer, dimension(:), allocatable, save nleafs_level
How many leaves are present per refinement level.
Definition: mod_forest.t:81
integer, dimension(:), allocatable, save sfc_to_igrid
Go from a Morton number to an igrid index (for a single processor)
Definition: mod_forest.t:53
integer nleafs_active
Definition: mod_forest.t:78
integer, dimension(:), allocatable, save morton_start
First Morton number per processor.
Definition: mod_forest.t:62
integer, save nleafs
Number of leaf block.
Definition: mod_forest.t:76
integer, dimension(:), allocatable, save morton_stop
Last Morton number per processor.
Definition: mod_forest.t:65
integer, dimension(:,:), allocatable, save sfc
Array to go from a Morton number to an igrid and processor index. Sfc(1:3, MN) contains [igrid,...
Definition: mod_forest.t:43
integer, save nparents
Number of parent blocks.
Definition: mod_forest.t:73
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Definition: mod_forest.t:32
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
This module contains definitions of global parameters and variables and some generic functions/subrou...
character(len=std_len), dimension(:), allocatable typeentropy
Which type of entropy fix to use with Riemann-type solvers.
double precision, dimension(:), allocatable w_convert_factor
Conversion factors the primitive variables.
type(state), pointer block
Block pointer for using one block and its previous state.
double precision xload
Stores the memory and load imbalance, used in printlog.
integer nstep
How many sub-steps the time integrator takes.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
integer it_max
Stop the simulation after this many time steps have been taken.
logical internalboundary
if there is an internal boundary
character(len=std_len) filename_sxr
Base file name for synthetic SXR emission output.
integer spectrum_wl
wave length for spectrum
logical nocartesian
IO switches for conversion.
integer, dimension(:), allocatable typepred1
The spatial discretization for the predictor step when using a two step PC method.
double precision dtdiffpar
For resistive MHD, the time step is also limited by the diffusion time: .
character(len=std_len) typegrad
logical reset_it
If true, reset iteration count to 0.
double precision small_pressure
integer ixghi
Upper index of grid block arrays.
character(len=std_len) geometry_name
integer, dimension(3, 3, 3) lvc
Levi-Civita tensor.
logical source_split_usr
Use split or unsplit way to add user's source terms, default: unsplit.
integer domain_nx
number of cells for each dimension in level-one mesh
integer, parameter unitpar
file handle for IO
character(len=std_len) filename_spectrum
Base file name for synthetic EUV spectrum output.
logical resume_previous_run
If true, restart a previous run from the latest snapshot.
double precision global_time
The global simulation time.
integer, dimension(nsavehi, nfile) itsave
Save output of type N on iterations itsave(:, N)
double precision time_max
End time for the simulation.
integer, dimension(3, 3) kr
Kronecker delta tensor.
logical, dimension(:), allocatable logflag
double precision time_init
Start time for the simulation.
logical stretch_uncentered
If true, adjust mod_geometry routines to account for grid stretching (but the flux computation will n...
logical firstprocess
If true, call initonegrid_usr upon restarting.
integer snapshotini
Resume from the snapshot with this index.
double precision small_temperature
error handling
double precision xprob
minimum and maximum domain boundaries for each dimension
integer it
Number of time steps taken.
character(len=std_len) filename_euv
Base file name for synthetic EUV emission output.
logical, dimension(:), allocatable loglimit
double precision, dimension(:), allocatable dg
extent of grid blocks in domain per dimension, in array over levels
integer it_init
initial iteration count
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
integer ditregrid
Reconstruct the AMR grid once every ditregrid iteration(s)
logical saveprim
If true, convert from conservative to primitive variables in output.
character(len=std_len) convert_type
Which format to use when converting.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer itfixgrid
Fix the AMR grid after this many time steps.
integer, parameter filecollapse_
Constant indicating collapsed output.
logical angmomfix
Enable to strictly conserve the angular momentum (works both in cylindrical and spherical coordinates...
double precision, dimension(:), allocatable amr_wavefilter
refinement: lohner estimate wavefilter setting
integer, parameter nlevelshi
The maximum number of levels in the grid refinement.
double precision location_slit
location of the slit
logical save_physical_boundary
True for save physical boundary cells in dat files.
logical stagger_grid
True for using stagger grid.
double precision time_convert_factor
Conversion factor for time unit.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer icomm
The MPI communicator.
integer, dimension(:), allocatable ng
number of grid blocks in domain per dimension, in array over levels
logical reset_time
If true, reset iteration count and global_time to original values, and start writing snapshots at ind...
integer, parameter nsavehi
Maximum number of saves that can be defined by tsave or itsave.
integer mype
The rank of the current MPI task.
double precision dtpar
If dtpar is positive, it sets the timestep dt, otherwise courantpar is used to limit the time step ba...
character(len=std_len) typediv
integer block_nx
number of cells for each dimension in grid block excluding ghostcells
integer type_block_io
MPI type for IO: block excluding ghost cells.
double precision, dimension(nfile) tsavestart
Start of read out (not counting specified read outs)
integer, dimension(nfile) ditsave
Repeatedly save output of type N when ditsave(N) time steps have passed.
integer, parameter plevel_
logical, dimension(ndim) collapse
If collapse(DIM) is true, generate output integrated over DIM.
integer, dimension(:), allocatable, parameter d
character(len=std_len) usr_filename
User parameter file.
logical ghost_copy
whether copy values instead of interpolation in ghost cells of finer blocks
double precision length_convert_factor
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision imex222_lambda
IMEX-222(lambda) one-parameter family of schemes.
double precision courantpar
The Courant (CFL) number used for the simulation.
double precision wall_time_max
Ending wall time (in hours) for the simulation.
integer ixm
the mesh range (within a block with ghost cells)
integer ierrmpi
A global MPI error return code.
logical autoconvert
If true, already convert to output format during the run.
integer, dimension(:), allocatable flux_method
Which flux scheme of spatial discretization to use (per grid level)
character(len=std_len) collapse_type
integer slowsteps
If > 1, then in the first slowsteps-1 time steps dt is reduced by a factor .
integer ssprk_order
SSPRK choice of methods (both threestep and fourstep, Shu-Osher 2N* implementation) also fivestep SSP...
double precision image_rotate
rotation of image
character(len=std_len) resolution_sxr
resolution of the SXR image
integer snapshotnext
IO: snapshot and collapsed views output numbers/labels.
integer, dimension(ndim) nstretchedblocks_baselevel
(even) number of (symmetrically) stretched blocks at AMR level 1, per dimension
integer npe
The number of MPI tasks.
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision, dimension(ndim) qstretch_baselevel
stretch factor between cells at AMR level 1, per dimension
double precision, dimension(ndim, 2) writespshift
integer imex_switch
IMEX_232 choice and parameters.
double precision time_between_print
to monitor timeintegration loop at given wall-clock time intervals
integer, parameter unitterm
Unit for standard output.
double precision, dimension(nfile) dtsave
Repeatedly save output of type N when dtsave(N) simulation time has passed.
logical, dimension(:), allocatable w_write
integer iprob
problem switch allowing different setups in same usr_mod.t
character(len=std_len) restart_from_file
If not 'unavailable', resume from snapshot with this base file name.
logical, dimension(ndim) periodb
True for dimensions with periodic boundaries.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter filelog_
Constant indicating log output.
logical final_dt_reduction
If true, allow final dt reduction for matching time_max on output.
logical, dimension(:), allocatable writelevel
integer, parameter fileout_
Constant indicating regular output.
double precision los_theta
direction of the line of sight (LOS)
character(len=std_len) typetvd
Which type of TVD method to use.
integer nbufferx
Number of cells as buffer zone.
double precision, dimension(:), allocatable entropycoef
double precision, dimension(:,:), allocatable dx
double precision tfixgrid
Fix the AMR grid after this time.
integer, parameter unitsnapshot
integer nghostcells
Number of ghost cells surrounding a grid.
double precision, dimension(:), allocatable w_refine_weight
Weights of variables used to calculate error for mesh refinement.
character(len=std_len) typedimsplit
character(len=std_len) typeaverage
character(len= *), parameter undefined
double precision, dimension(nsavehi, nfile) tsave
Save output of type N on times tsave(:, N)
double precision spectrum_window_max
logical convert
If true and restart_from_file is given, convert snapshots to other file formats.
character(len=std_len) resolution_euv
resolution of the EUV image
integer wavelength
wavelength for output
integer collapselevel
The level at which to produce line-integrated / collapsed output.
logical reset_grid
If true, rebuild the AMR grid upon restarting.
character(len=std_len) resolution_spectrum
resolution of the spectrum
double precision, dimension(:), allocatable refine_threshold
Error tolerance for refinement decision.
character(len=std_len) base_filename
Base file name for simulation output, which will be followed by a number.
double precision small_density
double precision spectrum_window_min
spectral window
double precision dtmin
Stop the simulation when the time step becomes smaller than this value.
integer refine_max_level
Maximal number of AMR levels.
integer, parameter fileslice_
Constant indicating slice output.
double precision, dimension(:), allocatable derefine_ratio
integer, parameter nfile
Number of output methods.
integer max_blocks
The maximum number of grid blocks in a processor.
integer, parameter fileanalysis_
Constant indicating analysis output (see Writing a custom analysis subroutine)
integer rk3_switch
RK3 Butcher table.
character(len=std_len) typefilelog
Which type of log to write: 'normal', 'special', 'regression_test'.
integer direction_slit
direction of the slit
logical pass_wall_time
If true, wall time is up, modify snapshotnext for later overwrite.
double precision, dimension(1:3) x_origin
where the is the origin (X=0,Y=0) of image
logical final_dt_exit
Force timeloop exit when final dt < dtmin.
integer, dimension(nfile) isaveit
integer, dimension(:,:), allocatable node
integer, dimension(nfile) isavet
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
integer log_fh
MPI file handle for logfile.
Module for reading input and writing output.
subroutine saveamrfile(ifile)
pure integer function, public count_ix(ixOL)
Compute number of elements in index range.
logical function snapshot_exists(ix)
integer function get_snapshot_index(filename)
subroutine read_par_files()
Read in the user-supplied parameter-file.
subroutine, public create_output_file(fh, ix, extension, suffix)
Standard method for creating a new output file.
character(len= *), parameter fmt_r
subroutine, public block_shape_io(igrid, n_ghost, ixOL, n_values)
Determine the shape of a block for output (whether to include ghost cells, and on which sides)
subroutine get_volume_average_func(func, f_avg, volume)
Compute the volume average of func(w) over the leaves of the grid.
subroutine get_global_minima(wmin)
Compute global minima of iw variables over the leaves of the grid.
logical save_now
whether a manually inserted snapshot is saved
subroutine get_global_maxima(wmax)
Compute global maxima of iw variables over the leaves of the grid.
subroutine printlog_regression_test()
Print a log that can be used to check whether the code still produces the same output (regression tes...
subroutine read_arguments()
Read the command line arguments passed to amrvac.
subroutine get_volume_average(power, mode, volume)
Compute mean(w**power) over the leaves of the grid. The first mode (power=1) corresponds to the mean,...
subroutine read_snapshot
Routine to read in snapshots (.dat files). When it cannot recognize the file version,...
integer nw_found
number of w found in dat files
character(len= *), parameter fmt_r2
integer, parameter version_number
Version number of the .dat file output.
subroutine get_fields_string(line, delims, n_max, fields, n_found, fully_read)
Routine to find entries in a string.
subroutine get_volume_coverage(vol_cov)
Compute how much of the domain is covered by each grid level. This routine does not take a non-Cartes...
character(len= *), parameter fmt_i
subroutine printlog_default
Write volume-averaged values and other information to the log file.
subroutine, public snapshot_write_header(fh, offset_tree, offset_block)
Write header for a snapshot.
subroutine, public snapshot_write_header1(fh, offset_tree, offset_block, dataset_names, nw_vars)
Write header for a snapshot, generalize cons_wnames and nw.
subroutine write_snapshot
subroutine snapshot_read_header(fh, offset_tree, offset_block)
Read header for a snapshot.
integer, dimension(3), parameter compatible_versions
List of compatible versions.
subroutine read_snapshot_old()
Module with slope/flux limiters.
Definition: mod_limiter.t:2
double precision cada3_radius
radius of the asymptotic region [0.001, 10], larger means more accurate in smooth region but more ove...
Definition: mod_limiter.t:13
double precision schmid_rad
Definition: mod_limiter.t:14
Module containing all the particle routines.
Definition: mod_particles.t:2
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_write_info), pointer phys_write_info
Definition: mod_physics.t:79
integer phys_wider_stencil
To use wider stencils in flux calculations. A value of 1 will extend it by one cell in both direction...
Definition: mod_physics.t:23
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition: mod_physics.t:19
logical phys_energy
Solve energy equation or not.
Definition: mod_physics.t:30
Writes D-1 slice, can do so in various formats, depending on slice_type.
Definition: mod_slice.t:2
integer slicenext
the file number of slices
Definition: mod_slice.t:13
subroutine write_slice
Definition: mod_slice.t:30
character(len=std_len) slice_type
choose data type of slice: vtu, vtuCC, dat, or csv
Definition: mod_slice.t:22
integer, dimension(nslicemax) slicedir
The slice direction for each slice.
Definition: mod_slice.t:19
integer nslices
Number of slices to output.
Definition: mod_slice.t:16
double precision, dimension(nslicemax) slicecoord
Slice coordinates, see Slice output.
Definition: mod_slice.t:10
Module for handling problematic values in simulations, such as negative pressures.
integer, public small_values_daverage
Average over this many cells in each direction.
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
character(len=20), public small_values_method
How to handle small values.
Module for handling split source terms (split from the fluxes)
Definition: mod_source.t:2
integer ittimelast
Definition: mod_timing.t:15
double precision timelast
Definition: mod_timing.t:16
Module with all the methods that users can customize in AMRVAC.
procedure(p_no_args), pointer usr_print_log
procedure(p_no_args), pointer usr_write_analysis
procedure(transform_w), pointer usr_transform_w
The data structure that contains information about a tree node/grid block.
Definition: mod_forest.t:11