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