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