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