MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_global_parameters.t
Go to the documentation of this file.
1!> This module contains definitions of global parameters and variables and some
2!> generic functions/subroutines used in AMRVAC.
3!>
5 use mpi
11
12 implicit none
13 public
14
15 !> minimum and maximum domain boundaries for each dimension
16 double precision :: xprob^l
17 !> store unstretched cell size of current level
18 double precision :: dxlevel(^nd)
19 !> stretch factor between cells at AMR level 1, per dimension
20 double precision :: qstretch_baselevel(^nd)
21
22 !> physical extent of stretched border in symmetric stretching
23 double precision :: xstretch^d
24
25 !> to monitor timeintegration loop at given wall-clock time intervals
26 double precision :: time_between_print
27
28 !> accumulated wall-clock time spent on boundary conditions
29 double precision :: time_bc
30
31 !> global time step
32 double precision :: dt
33
34 !> The Courant (CFL) number used for the simulation
35 double precision :: courantpar
36
37
38 !> If dtpar is positive, it sets the timestep dt, otherwise courantpar is used
39 !> to limit the time step based on the Courant condition.
40 double precision :: dtpar
41
42 !> For resistive MHD, the time step is also limited by the diffusion time:
43 !> \f$ dt < dtdiffpar \times dx^2/eta \f$
44 double precision :: dtdiffpar
45
46 !> The global simulation time
47 double precision :: global_time
48
49 !> Start time for the simulation
50 double precision :: time_init
51
52 !> End time for the simulation
53 double precision :: time_max
54
55 !> Ending wall time (in hours) for the simulation
56 double precision :: wall_time_max
57
58 !> Stop the simulation when the time step becomes smaller than this value
59 double precision :: dtmin
60
61 !> Conversion factor for length unit
62 double precision :: length_convert_factor
63
64 !> Conversion factor for time unit
65 double precision :: time_convert_factor
66
67 !> Fix the AMR grid after this time
68 double precision :: tfixgrid
69
70 ! Physics factors
71 !> Physical scaling factor for length
72 double precision :: unit_length=1.d0
73
74 !> Physical scaling factor for time
75 double precision :: unit_time=1.d0
76
77 !> Physical scaling factor for density
78 double precision :: unit_density=1.d0
79
80 !> Physical scaling factor for velocity
81 double precision :: unit_velocity=1.d0
82
83 !> Physical scaling factor for temperature
84 double precision :: unit_temperature=1.d0
85
86 !> Physical scaling factor for pressure
87 double precision :: unit_pressure=1.d0
88
89 !> Physical scaling factor for magnetic field
90 double precision :: unit_magneticfield=1.d0
91
92 !> Physical scaling factor for number density
93 double precision :: unit_numberdensity=1.d0
94
95 !> Physical scaling factor for charge
96 double precision :: unit_charge=1.d0
97
98 !> Physical scaling factor for mass
99 double precision :: unit_mass=1.d0
100
101 !> Normalised speed of light
102 double precision :: c_norm=1.d0
103
104 !> Normalised radiation constant
105 double precision :: arad_norm=1.d0
106
107 !> Physical scaling factor for Opacity
108 double precision :: unit_opacity=1.d0
109
110 !> Physical scaling factor for radiation flux
111 double precision :: unit_radflux=1.d0
112
113 !> Physical scaling factor for radiation energy density
114 double precision :: unit_erad=1.d0
115
116 !> Physical factors useful for radiation fld
118
119 !> error handling
121 double precision :: phys_trac_mask
122
123 !> amplitude of background dipolar, quadrupolar, octupolar, user's field
124 double precision :: bdip=0.d0
125 double precision :: bquad=0.d0
126 double precision :: boct=0.d0
127 double precision :: busr=0.d0
128 !> Stores the memory and load imbalance, used in printlog
129 double precision :: xload, xmemory
130 double precision :: tvdlfeps
131
132 !> RK2(alfa) method parameters from Butcher tableau
133 double precision :: rk_a21,rk_b1,rk_b2
134 !> IMEX-222(lambda) one-parameter family of schemes
135 double precision :: imex222_lambda
140 !> IMEX_ARS3 parameter ars_gamma
141 double precision :: ars_gamma
143 double precision :: imex_b3,imex_c2,imex_c3
144 !> IMEX_CB3a extra parameters
145 double precision :: imex_a22, imex_a33, imex_ha32
146
147 !> global fastest wave speed needed in fd scheme and glm method
148 double precision :: cmax_global
149
150 !> global fastest flow speed needed in glm method
151 double precision :: vmax_global
152
153 !> times for enhancing spatial resolution for EUV image/spectra
155 !> the white light emission below it (unit=Rsun) is not visible
156 double precision :: r_occultor
157 !> direction of the line of sight (LOS)
158 double precision :: los_theta,los_phi
159 !> rotation of image
160 double precision :: image_rotate
161 !> where the is the origin (X=0,Y=0) of image
162 double precision :: x_origin(1:3)
163 !> Base file name for synthetic EUV spectrum output
164 character(len=std_len) :: filename_spectrum
165 !> spectral window
167 !> location of the slit
168 double precision :: location_slit
169 !> for spherical coordinate, region below it (unit=Rsun) is treated as not transparent
170 double precision :: r_opt_thick
171
172 !> domain percentage cut off shifted from each boundary when converting data
173 double precision :: writespshift(^nd,2)
174
175 double precision, allocatable :: entropycoef(:)
176
177 !> extent of grid blocks in domain per dimension, in array over levels
178 double precision, dimension(:), allocatable :: dg^d
179 !> Corner coordinates
180 double precision, allocatable :: rnode(:,:)
181 double precision, allocatable :: rnode_sub(:,:)
182
183 !> spatial steps for all dimensions at all levels
184 double precision, allocatable :: dx(:,:)
185
186 !> Error tolerance for refinement decision
187 double precision, allocatable :: refine_threshold(:)
188 !> Error tolerance ratio for derefinement decision
189 double precision, allocatable :: derefine_ratio(:)
190
191 !> Stretching factors and first cell size for each AMR level and dimension
192 double precision, allocatable :: qstretch(:,:), dxfirst(:,:), &
193 dxfirst_1mq(:,:), dxmid(:,:)
194
195 !> Conversion factors the primitive variables
196 double precision, allocatable :: w_convert_factor(:)
197
198 !> Weights of variables used to calculate error for mesh refinement
199 double precision, allocatable :: w_refine_weight(:)
200
201 !> refinement: lohner estimate wavefilter setting
202 double precision, allocatable :: amr_wavefilter(:)
203
204 !> Maximum number of saves that can be defined by tsave or itsave
205 integer, parameter :: nsavehi=100
206 !> Number of output methods
207 integer, parameter :: nfile = 5
208
209 !> Save output of type N on times tsave(:, N)
210 double precision :: tsave(nsavehi,nfile)
211
212 double precision :: tsavelast(nfile)
213
214 !> Repeatedly save output of type N when dtsave(N) simulation time has passed
215 double precision :: dtsave(nfile)
216
217 !> Start of read out (not counting specified read outs)
218 double precision :: tsavestart(nfile)
219
220 !> Number of spatial dimensions for grid variables
221 integer, parameter :: ndim=^nd
222
223 !> The number of MPI tasks
224 integer :: npe
225
226 !> The rank of the current MPI task
227 integer :: mype
228
229 !> The MPI communicator
230 integer :: icomm
231
232 !> A global MPI error return code
233 integer :: ierrmpi
234
235 !> MPI file handle for logfile
236 integer :: log_fh
237 !> MPI type for block including ghost cells and its size
239 !> MPI type for block coarsened by 2, and for its children blocks
241 !> MPI type for staggered block coarsened by 2, and for its children blocks
242 integer :: type_coarse_block_stg(^nd,2^d&), type_sub_block_stg(^ND,2^D&)
243 !> MPI type for IO: block excluding ghost cells
245 !> MPI type for IO of staggered variables
247 !> MPI type for IO: cell corner (xc) or cell center (xcc) coordinates
249 !> MPI type for IO: cell corner (wc) or cell center (wcc) variables
251
252 ! geometry and domain setups
253 !> the mesh range of a physical block without ghost cells
254 integer :: ixm^ll
255
256 !> Indices for cylindrical coordinates FOR TESTS, negative value when not used:
257 integer :: r_ = -1
258 integer :: phi_ = -1
259 integer :: z_ = -1
260
261 !> Number of spatial dimensions (components) for vector variables
262 integer :: ndir=ndim
263
264 !> starting dimension for electric field
265 {^ifoned
266 integer, parameter :: sdim=3
267 }
268 {^iftwod
269 integer, parameter :: sdim=3
270 }
271 {^ifthreed
272 integer, parameter :: sdim=1
273 }
274
275 !> number of cells for each dimension in level-one mesh
276 integer :: domain_nx^d
277
278 !> number of cells for each dimension in grid block excluding ghostcells
279 integer :: block_nx^d
280
281 !> Lower index of grid block arrays (always 1)
282 integer, parameter :: {ixglo^d = 1|, }
283
284 !> Upper index of grid block arrays
285 integer :: ixghi^d
286
287 !> Lower index of stagger grid block arrays (always 0)
288 integer, parameter :: {ixgslo^d = 0|, }
289
290 !> Upper index of stagger grid block arrays
291 integer :: ixgshi^d
292
293 !> Number of ghost cells surrounding a grid
294 integer :: nghostcells = 2
295
296 integer, parameter :: stretch_none = 0 !< No stretching
297 integer, parameter :: stretch_uni = 1 !< Unidirectional stretching from a side
298 integer, parameter :: stretch_symm = 2 !< Symmetric stretching around the center
299
300 !> What kind of stretching is used per dimension
301 integer :: stretch_type(ndim)
302 !> (even) number of (symmetrically) stretched
303 !> blocks at AMR level 1, per dimension
305 !> (even) number of (symmetrically) stretched blocks per level and dimension
306 integer, allocatable :: nstretchedblocks(:,:)
307
308 !> grid hierarchy info (level and grid indices)
309 integer, parameter :: nodehi=^nd+1
310 integer, parameter :: plevel_=1
311 integer, parameter :: pig^d_=plevel_+^d
312
313 integer, allocatable :: node(:,:)
314 integer, allocatable :: node_sub(:,:)
315
316 !> grid location info (corner coordinates and grid spacing)
317 integer, parameter :: rnodehi=3*^nd
318 integer, parameter :: rpxmin0_=0
319 integer, parameter :: rpxmin^d_=rpxmin0_+^d
320 integer, parameter :: rpxmax0_=^nd
321 integer, parameter :: rpxmax^d_=rpxmax0_+^d
322 integer, parameter :: rpdx^d_=2*^nd+^d
323
324 !> index number of the latest existing data file
326
327 !> Save output of type N on iterations itsave(:, N)
328 integer :: itsave(nsavehi,nfile)
329
330 integer :: itsavelast(nfile)
331
332 !> Repeatedly save output of type N when ditsave(N) time steps have passed
333 integer :: ditsave(nfile)
334
335 integer :: isavet(nfile)
336
337 integer :: isaveit(nfile)
338
339 !> The level at which to produce line-integrated / collapsed output
340 integer :: collapselevel
341
342 !> Number of saved files of each type
343 integer :: n_saves(1:nfile)
344
345 !> IO: snapshot and collapsed views output numbers/labels
347
348 !> IO: slice output number/label
349 integer :: slicenext
350
351 !> Constant indicating log output
352 integer, parameter :: filelog_ = 1
353
354 !> Constant indicating regular output
355 integer, parameter :: fileout_ = 2
356
357 !> Constant indicating slice output
358 integer, parameter :: fileslice_ = 3
359
360 !> Constant indicating collapsed output
361 integer, parameter :: filecollapse_ = 4
362
363 !> Constant indicating analysis output (see @ref analysis.md)
364 integer, parameter :: fileanalysis_ = 5
365
366 !> Unit for standard input
367 integer, parameter :: unitstdin=5
368
369 !> Unit for standard output
370 integer, parameter :: unitterm=6
371
372 !> Unit for error messages
373 integer, parameter :: uniterr=6
374
375 !> file handle for IO
376 integer, parameter :: unitpar=9
377 integer, parameter :: unitconvert=10
378 integer, parameter :: unitslice=11
379 integer, parameter :: unitsnapshot=12
380 integer, parameter :: unitcollapse=13
381 integer, parameter :: unitanalysis=14
382
383 !> Number of auxiliary variables that are only included in output
384 integer :: nwauxio
385
387
388 !> Resume from the snapshot with this index
389 integer :: snapshotini
390
391 !> number of equilibrium set variables, besides the mag field
392 integer :: number_equi_vars = 0
393
394 integer :: phys_trac_type=1
396
397 !> integer switchers for type courant
398 integer, parameter :: type_maxsum=1
399 integer, parameter :: type_summax=2
400 integer, parameter :: type_minimum=3
401
402 ! AMR switches
403 !> The maximum number of grid blocks in a processor
404 integer :: max_blocks
405
406 !> The maximum number of levels in the grid refinement
407 integer, parameter :: nlevelshi = 20
408
409 !> Maximal number of AMR levels
411
412 !> Fix the AMR grid after this many time steps
413 integer :: itfixgrid
414
415 !> Reconstruct the AMR grid once every ditregrid iteration(s)
416 integer :: ditregrid
417
418 !> select types of refine criterion
420
421 !> Number of cells as buffer zone
422 integer :: nbufferx^d
423
424 integer :: levmin
425 integer :: levmax
426 integer :: levmax_sub
427
428 ! Miscellaneous
429 !> problem switch allowing different setups in same usr_mod.t
430 integer :: iprob
431
432 !> Kronecker delta tensor
433 integer :: kr(3,3)
434
435 !> Levi-Civita tensor
436 integer :: lvc(3,3,3)
437
438 !> How to compute the CFL-limited time step
439 integer :: type_courant=1
440
441 !> Number of time steps taken
442 integer :: it
443
444 !> Stop the simulation after this many time steps have been taken
445 integer :: it_max
446
447 !> initial iteration count
448 integer :: it_init
449
450 !> If > 1, then in the first slowsteps-1 time steps dt is reduced
451 !> by a factor \f$ 1 - (1- step/slowsteps)^2 \f$
452 integer :: slowsteps
453
454 ! Method switches
455
456 !> Index of the sub-step in a multi-step time integrator
457 integer :: istep
458
459 !> How many sub-steps the time integrator takes
460 integer :: nstep
461
462 !> flux schemes
463 integer, parameter :: fs_hll=1
464 integer, parameter :: fs_hllc=2
465 integer, parameter :: fs_hlld=3
466 integer, parameter :: fs_hllcd=4
467 integer, parameter :: fs_tvdlf=5
468 integer, parameter :: fs_tvdmu=6
469 integer, parameter :: fs_tvd=7
470 integer, parameter :: fs_hancock=8
471 integer, parameter :: fs_cd=9
472 integer, parameter :: fs_cd4=10
473 integer, parameter :: fs_fd=11
474 integer, parameter :: fs_source=12
475 integer, parameter :: fs_nul=13
476
477 !> time stepper type
478 integer :: t_stepper=0
479 integer, parameter :: onestep=1
480 integer, parameter :: twostep=2
481 integer, parameter :: threestep=3
482 integer, parameter :: fourstep=4
483 integer, parameter :: fivestep=5
484
485 !> time integrator method
486 integer :: t_integrator=0
487 integer, parameter :: forward_euler=1
488 integer, parameter :: predictor_corrector=2
489 integer, parameter :: ssprk3=3
490 integer, parameter :: ssprk4=4
491 integer, parameter :: ssprk5=5
492
493 integer, parameter :: imex_euler=6
494 integer, parameter :: imex_sp=7
495 integer, parameter :: rk2_alf=8
496 integer, parameter :: ssprk2=9
497 integer, parameter :: imex_midpoint=10
498 integer, parameter :: imex_trapezoidal=11
499 integer, parameter :: imex_222=12
500
501 integer, parameter :: rk3_bt=13
502 integer, parameter :: imex_ars3=14
503 integer, parameter :: imex_232=15
504 integer, parameter :: imex_cb3a=16
505
506 integer, parameter :: rk4=17
507
508 !> number of grid blocks in domain per dimension, in array over levels
509 integer, dimension(:), allocatable :: ng^d
510
511 !> Which flux scheme of spatial discretization to use (per grid level)
512 integer, allocatable :: flux_method(:)
513
514 !> The spatial discretization for the predictor step when using a two
515 !> step PC method
516 integer, allocatable :: typepred1(:)
517
518 !> Type of slope limiter used for reconstructing variables on cell edges
519 integer, allocatable :: type_limiter(:)
520
521 !> Type of slope limiter used for computing gradients or divergences, when
522 !> typegrad or typediv are set to 'limited'
523 integer, allocatable :: type_gradient_limiter(:)
524
525 !> background magnetic field location indicator
526 integer :: b0i=0
527
528 !> Limiter used for prolongation to refined grids and ghost cells
529 integer :: prolong_limiter=0
530
531 !> bound (left/min and right.max) speed of Riemann fan
532 integer :: boundspeed
533
534 integer :: nxdiffusehllc
535 !> SSPRK choice of methods (both threestep and fourstep, Shu-Osher 2N* implementation)
536 !> also fivestep SSPRK54
537 integer :: ssprk_order
538 !> RK3 Butcher table
539 integer :: rk3_switch
540 !> IMEX_232 choice and parameters
541 integer :: imex_switch
542
543 !> Array indicating the type of boundary condition per variable and per
544 !> physical boundary
545 integer, allocatable :: typeboundary(:, :)
546 !> boundary condition types
547 integer, parameter :: bc_special=1
548 integer, parameter :: bc_cont=2
549 integer, parameter :: bc_symm=3
550 integer, parameter :: bc_asymm=4
551 integer, parameter :: bc_periodic=5
552 integer, parameter :: bc_aperiodic=6
553 integer, parameter :: bc_noinflow=7
554 integer, parameter :: bc_data=8
555 integer, parameter :: bc_character=9
556 integer, parameter :: bc_icarus=10
557
558 !> wavelength for output
559 integer :: wavelength
560 ! minimum and maximum energy of SXR (keV)
562 !> wave length for spectrum
563 integer :: spectrum_wl
564 !> direction of the slit (for dat resolution only)
565 integer :: direction_slit
566
567
568 !> Cartesian geometry or not
569 logical :: slab
570
571 !> uniform Cartesian geometry or not (stretched Cartesian)
572 logical :: slab_uniform
573
574 !> each cell has its own timestep or not
575 logical :: local_timestep = .false.
576
577 !> whether or not to save an output file
578 logical :: save_file(nfile)
579
580 !> If true, adjust mod_geometry routines to account for grid stretching (but
581 !> the flux computation will not)
583 !> True if a dimension is stretched
584 logical :: stretched_dim(ndim)
585
586 !> If true, restart a previous run from the latest snapshot
588
589 !> If true and restart_from_file is given, convert snapshots to
590 !> other file formats
591 logical :: convert
592
593 !> If true, already convert to output format during the run
594 logical :: autoconvert
595
596 !> If true, convert from conservative to primitive variables in output
597 logical :: saveprim
598
599 !> do time evolving
600 logical :: time_advance
601
602 !> Force timeloop exit when final dt < dtmin
603 logical :: final_dt_exit
604
605 !> If true, reset iteration count and global_time to original values, and
606 !> start writing snapshots at index 0
607 logical :: reset_time
608
609 !> If true, reset iteration count to 0
610 logical :: reset_it
611
612 !> If true, allow final dt reduction for matching time_max on output
614
615 !> If true, call initonegrid_usr upon restarting
616 logical :: firstprocess
617
618 !> If true, wall time is up, modify snapshotnext for later overwrite
619 logical :: pass_wall_time
620
621 !> If true, do H-correction to fix the carbuncle problem at grid-aligned shocks
622 logical :: h_correction=.false.
623
624 !> If true, rebuild the AMR grid upon restarting
625 logical :: reset_grid
626 !> True for using stagger grid
627 logical :: stagger_grid=.false.
628 !> True for record electric field
629 logical :: record_electric_field=.false.
630
631 !> resolution of the images
632 logical :: dat_resolution
633
634 !> If collapse(DIM) is true, generate output integrated over DIM
635 logical :: collapse(ndim)
636 !> IO switches for conversion
637 logical :: nocartesian
638
639 !> Use particles module or not
640 logical :: use_particles=.false.
641
642 !> Use multigrid (only available in 2D and 3D)
643 logical :: use_multigrid = .false.
644
645 !> prolongate primitive variables in level-jump ghost cells
646 logical :: prolongprimitive=.false.
647
648 !> coarsen primitive variables in level-jump ghost cells
649 logical :: coarsenprimitive=.false.
650
651 !> Save a snapshot before crash a run met unphysical values
652 logical :: crash=.false.
653
654 !> check and optionally fix unphysical small values (density, gas pressure)
655 logical :: check_small_values=.true.
656
657 !> fix small values with average or replace methods
658 logical :: fix_small_values=.false.
659
660 !> split magnetic field as background B0 field
661 logical :: b0field=.false.
662 logical :: b0fieldalloccoarse=.false.
663 !> Use SI units (.true.) or use cgs units (.false.)
664 logical :: si_unit=.false.
665
666 !> Use TRAC for MHD or 1D HD
667 logical :: phys_trac=.false.
668
669 !> Whether to apply flux conservation at refinement boundaries
670 logical :: fix_conserve_global = .true.
673 !> Use split or unsplit way to add user's source terms, default: unsplit
675 !> if any normal source term is added in split fasion
676 logical :: any_source_split=.false.
677 logical :: dimsplit
678 !> whether IMEX in use or not
680
681 !> need global maximal wave speed
682 logical :: need_global_cmax=.false.
683
684 ! Boundary region parameters
685
686 !> True for dimensions with periodic boundaries
687 logical :: periodb(ndim)
688
689 !> Indicates whether there is a pole at a boundary
690 logical :: poleb(2,ndim)
691
692 !> True for dimensions with aperiodic boundaries
693 logical :: aperiodb(ndim)
694
695 !> True for save physical boundary cells in dat files
697
698 !> whether copy values instead of interpolation in ghost cells of finer blocks
699 logical :: ghost_copy=.false.
700
701 !> if there is an internal boundary
703
704 !> use arcsec as length unit of images/spectra
706 !> big image
707 logical :: big_image
708
709 !> True if a block has any physical boundary
710 logical, allocatable :: phyboundblock(:)
711
712 !> if true write the w variable in output
713 logical, allocatable :: w_write(:)
714
715 logical, allocatable :: writelevel(:)
716
717 logical, allocatable :: loglimit(:), logflag(:)
718
719 ! Parameters
720 character(len=*), parameter :: undefined = 'undefined'
721
722 !> Names of the output methods
723 character(len=40), parameter :: output_names(nfile) = &
724 ['log ', 'normal ', 'slice ', 'collapsed', 'analysis ']
725 !> Which format to use when converting
726 !>
727 !> Options are: tecplot, tecplotCC, vtu, vtuCC, vtuB, vtuBCC,
728 !> tecplotmpi, tecplotCCmpi, vtumpi, vtuCCmpi, vtuBmpi, vtuBCCmpi, pvtumpi, pvtuCCmpi,
729 !> pvtuBmpi, pvtuBCCmpi, tecline, teclinempi, onegrid
730 character(len=std_len) :: convert_type
731
732 character(len=std_len) :: collapse_type
733
734 !> User parameter file
735 character(len=std_len) :: usr_filename
736
737 !> Base file name for simulation output, which will be followed by a number
738 character(len=std_len) :: base_filename
739
740 !> If not 'unavailable', resume from snapshot with this base file name
741 character(len=std_len) :: restart_from_file
742
743 !> Which type of log to write: 'normal', 'special', 'regression_test'
744 character(len=std_len) :: typefilelog
745
746 character(len=std_len) :: typeaverage
747 character(len=std_len) :: typedimsplit
748 character(len=std_len) :: geometry_name='default'
749 character(len=std_len) :: typepoly
750 !> Base file name for synthetic EUV emission output
751 character(len=std_len) :: filename_euv
752 !> Base file name for synthetic SXR emission output
753 character(len=std_len) :: filename_sxr
754 !> Base file name for synthetic white light
755 character(len=std_len) :: filename_whitelight
756 !> white light observation instrument
757 character(len=std_len) :: whitelight_instrument
758
759 !> Which type of TVD method to use
760 character(len=std_len) :: typetvd
761
762 character(len=std_len) :: typediv,typegrad
763
764 !> Which par files are used as input
765 character(len=std_len), allocatable :: par_files(:)
766
767 !> Which type of entropy fix to use with Riemann-type solvers
768 character(len=std_len), allocatable :: typeentropy(:)
769
770 !> Block pointer for using one block and its previous state
771 type(state), pointer :: block
772
773 !$OMP THREADPRIVATE(block,dxlevel,b0i)
774
775contains
776
777 !> Cross product of two vectors
778 pure subroutine cross_product(ixI^L,ixO^L,a,b,axb)
779 integer, intent(in) :: ixi^l, ixo^l
780 double precision, intent(in) :: a(ixi^s,3), b(ixi^s,3)
781 double precision, intent(out) :: axb(ixi^s,3)
782
783 axb(ixo^s,1)=a(ixo^s,2)*b(ixo^s,3)-a(ixo^s,3)*b(ixo^s,2)
784 axb(ixo^s,2)=a(ixo^s,3)*b(ixo^s,1)-a(ixo^s,1)*b(ixo^s,3)
785 axb(ixo^s,3)=a(ixo^s,1)*b(ixo^s,2)-a(ixo^s,2)*b(ixo^s,1)
786 end subroutine cross_product
787
788end module mod_global_parameters
Module with basic data types used in amrvac.
This module contains variables that describe the connectivity of the mesh and also data structures fo...
Module for physical and numeric constants.
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(nfile) tsavelast
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, dimension(ndim) aperiodb
True for dimensions with aperiodic boundaries.
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.
double precision arad_norm
Normalised radiation constant.
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
integer ixgshi
Upper index of stagger grid block arrays.
logical reset_it
If true, reset iteration count to 0.
double precision unit_charge
Physical scaling factor for charge.
integer, parameter bc_noinflow
integer type_coarse_block
MPI type for block coarsened by 2, and for its children blocks.
integer ixghi
Upper index of grid block arrays.
integer, parameter stretch_uni
Unidirectional stretching from a side.
pure subroutine cross_product(ixil, ixol, a, b, axb)
Cross product of two vectors.
character(len=std_len) geometry_name
integer, dimension(3, 3, 3) lvc
Levi-Civita tensor.
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
logical activate_unit_arcsec
use arcsec as length unit of images/spectra
character(len=std_len) typepoly
logical source_split_usr
Use split or unsplit way to add user's source terms, default: unsplit.
integer, parameter imex_euler
double precision unit_opacity
Physical scaling factor for Opacity.
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, dimension(nfile) save_file
whether or not to save an output file
logical any_source_split
if any normal source term is added in split fasion
logical resume_previous_run
If true, restart a previous run from the latest snapshot.
double precision global_time
The global simulation time.
integer type_block_xc_io
MPI type for IO: cell corner (xc) or cell center (xcc) coordinates.
integer, dimension(nsavehi, nfile) itsave
Save output of type N on iterations itsave(:, N)
double precision unit_mass
Physical scaling factor for mass.
logical use_imex_scheme
whether IMEX in use or not
integer istep
Index of the sub-step in a multi-step time integrator.
double precision time_max
End time for the simulation.
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision xstretch
physical extent of stretched border in symmetric stretching
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, parameter bc_character
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 ars_gamma
IMEX_ARS3 parameter ars_gamma.
double precision unit_numberdensity
Physical scaling factor for number density.
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.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter type_maxsum
integer switchers for type courant
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.
integer prolong_limiter
Limiter used for prolongation to refined grids and ghost cells.
double precision, dimension(:), allocatable amr_wavefilter
refinement: lohner estimate wavefilter setting
double precision unit_length
Physical scaling factor for length.
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.
double precision vmax_global
global fastest flow speed needed in glm method
logical stagger_grid
True for using stagger grid.
double precision const_rad_a
Physical factors useful for radiation fld.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
double precision time_convert_factor
Conversion factor for time unit.
logical, dimension(:), allocatable phyboundblock
True if a block has any physical boundary.
integer, dimension(:,:), allocatable nstretchedblocks
(even) number of (symmetrically) stretched blocks per level and dimension
integer, dimension(^nd, 2^d &) type_coarse_block_stg
MPI type for staggered block coarsened by 2, and for its children blocks.
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...
double precision bdip
amplitude of background dipolar, quadrupolar, octupolar, user's field
integer b0i
background magnetic field location indicator
integer, parameter imex_trapezoidal
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...
integer, dimension(1:nfile) n_saves
Number of saved files of each type.
character(len=std_len) typediv
integer block_nx
number of cells for each dimension in grid block excluding ghostcells
integer type_block_io
MPI type for IO: block excluding ghost cells.
double precision, dimension(nfile) tsavestart
Start of read out (not counting specified read outs)
integer, dimension(nfile) ditsave
Repeatedly save output of type N when ditsave(N) time steps have passed.
integer, dimension(2^d &) type_sub_block
logical, dimension(ndim) collapse
If collapse(DIM) is true, generate output integrated over DIM.
integer, parameter unitstdin
Unit for standard input.
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, parameter nodehi
grid hierarchy info (level and grid indices)
integer ndir
Number of spatial dimensions (components) for vector variables.
integer, parameter uniterr
Unit for error messages.
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
logical slab
Cartesian geometry or not.
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
double precision, dimension(:,:), allocatable qstretch
Stretching factors and first cell size for each AMR level and dimension.
integer, parameter bc_periodic
integer type_block_wc_io
MPI type for IO: cell corner (wc) or cell center (wcc) variables.
integer type_courant
How to compute the CFL-limited time step.
integer, parameter bc_special
boundary condition types
integer snapshotnext
IO: snapshot and collapsed views output numbers/labels.
double precision unit_magneticfield
Physical scaling factor for magnetic field.
integer, parameter unitanalysis
logical, dimension(ndim) stretched_dim
True if a dimension is stretched.
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.
double precision unit_velocity
Physical scaling factor for velocity.
integer index_latest_data
index number of the latest existing data file
integer, dimension(nfile) itsavelast
integer, parameter stretch_none
No stretching.
integer imex_switch
IMEX_232 choice and parameters.
integer, parameter fs_hll
flux schemes
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 time_advance
do time evolving
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 c_norm
Normalised speed of light.
logical b0field
split magnetic field as background B0 field
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter filelog_
Constant indicating log output.
character(len=40), dimension(nfile), parameter output_names
Names of the output methods.
double precision unit_temperature
Physical scaling factor for temperature.
double precision unit_radflux
Physical scaling factor for radiation flux.
integer, dimension(:), allocatable type_limiter
Type of slope limiter used for reconstructing variables on cell edges.
logical final_dt_reduction
If true, allow final dt reduction for matching time_max on output.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
logical, dimension(:), allocatable writelevel
integer, parameter unitcollapse
integer, parameter fileout_
Constant indicating regular output.
double precision los_theta
direction of the line of sight (LOS)
integer, parameter type_summax
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 time_bc
accumulated wall-clock time spent on boundary conditions
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer type_block
MPI type for block including ghost cells and its size.
double precision rk_a21
RK2(alfa) method parameters from Butcher tableau.
integer, parameter predictor_corrector
double precision tfixgrid
Fix the AMR grid after this time.
integer, parameter unitsnapshot
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.
integer, parameter sdim
starting dimension for electric field
integer, parameter forward_euler
logical phys_trac
Use TRAC for MHD or 1D HD.
logical fix_conserve_global
Whether to apply flux conservation at refinement boundaries.
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 need_global_cmax
need global maximal wave speed
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.
integer, dimension(ndim) stretch_type
What kind of stretching is used per dimension.
logical crash
Save a snapshot before crash a run met unphysical values.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
integer t_stepper
time stepper type
double precision, dimension(:,:), allocatable rnode_sub
logical use_multigrid
Use multigrid (only available in 2D and 3D)
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
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.
integer, parameter rnodehi
grid location info (corner coordinates and grid spacing)
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 stretch_symm
Symmetric stretching around the center.
integer, parameter fileanalysis_
Constant indicating analysis output (see Writing a custom analysis subroutine)
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer rk3_switch
RK3 Butcher table.
integer, parameter bc_aperiodic
character(len=std_len) typefilelog
Which type of log to write: 'normal', 'special', 'regression_test'.
double precision imex_a22
IMEX_CB3a extra parameters.
integer direction_slit
direction of the slit (for dat resolution only)
integer type_block_io_stg
MPI type for IO of staggered variables.
logical, dimension(2, ndim) poleb
Indicates whether there is a pole at a boundary.
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
logical pass_wall_time
If true, wall time is up, modify snapshotnext for later overwrite.
integer, parameter ixgslo
Lower index of stagger grid block arrays (always 0)
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 record_electric_field
True for record electric field.
integer, parameter unitconvert
double precision unit_erad
Physical scaling factor for radiation energy density.
double precision, dimension(:,:), allocatable dxfirst
integer t_integrator
time integrator method
integer, parameter fs_hancock
logical final_dt_exit
Force timeloop exit when final dt < dtmin.
integer, dimension(nfile) isaveit
integer, dimension(:), allocatable type_gradient_limiter
Type of slope limiter used for computing gradients or divergences, when typegrad or typediv are set t...
integer number_equi_vars
number of equilibrium set variables, besides the mag field
integer, dimension(:,:), allocatable node
integer, dimension(:,:), allocatable node_sub
integer, parameter type_minimum
integer, dimension(nfile) isavet
integer, parameter imex_midpoint
double precision, dimension(:,:), allocatable dxfirst_1mq
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
double precision, dimension(:,:), allocatable dxmid
integer, parameter ixglo
Lower index of grid block arrays (always 1)
integer log_fh
MPI file handle for logfile.