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