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