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