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