MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
|
Magneto-hydrodynamics module. More...
Functions/Subroutines | |
subroutine, public | mhd_phys_init () |
subroutine, public | mhd_ei_to_e (ixil, ixol, w, x) |
Transform internal energy to total energy. | |
subroutine, public | mhd_e_to_ei (ixil, ixol, w, x) |
Transform total energy to internal energy. | |
subroutine, public | mhd_get_v (w, x, ixil, ixol, v) |
Calculate v vector. | |
subroutine, public | multiplyambicoef (ixil, ixol, res, w, x) |
multiply res by the ambipolar coefficient The ambipolar coefficient is calculated as -mhd_eta_ambi/rho^2 The user may mask its value in the user file by implemneting usr_mask_ambipolar subroutine | |
subroutine, public | mhd_get_rho (w, x, ixil, ixol, rho) |
subroutine, public | get_normalized_divb (w, ixil, ixol, divb) |
get dimensionless div B = |divB| * volume / area / |B| | |
subroutine, public | get_current (w, ixil, ixol, idirmin, current) |
Calculate idirmin and the idirmin:3 components of the common current array make sure that dxlevel(^D) is set correctly. | |
double precision function, dimension(ixo^s), public | mhd_mag_en_all (w, ixil, ixol) |
Compute 2 times total magnetic energy. | |
subroutine, public | mhd_clean_divb_multigrid (qdt, qt, active) |
subroutine, public | mhd_face_to_center (ixol, s) |
calculate cell-center values from face-center values | |
subroutine, public | b_from_vector_potential (ixisl, ixil, ixol, ws, x) |
calculate magnetic field from vector potential | |
Variables | |
double precision, public | mhd_gamma = 5.d0/3.0d0 |
The adiabatic index. | |
double precision, public | mhd_adiab = 1.0d0 |
The adiabatic constant. | |
double precision, public | mhd_eta = 0.0d0 |
The MHD resistivity. | |
double precision, public | mhd_eta_hyper = 0.0d0 |
The MHD hyper-resistivity. | |
double precision, public | mhd_etah = 0.0d0 |
Hall resistivity. | |
double precision, public | mhd_eta_ambi = 0.0d0 |
The MHD ambipolar coefficient. | |
double precision, public, protected | mhd_trac_mask = 0.d0 |
Height of the mask used in the TRAC method. | |
double precision, public | mhd_glm_alpha = 0.5d0 |
GLM-MHD parameter: ratio of the diffusive and advective time scales for div b taking values within [0, 1]. | |
double precision, public, protected | mhd_reduced_c = 0.02d0*const_c |
Reduced speed of light for semirelativistic MHD: 2% of light speed. | |
double precision, public | hypertc_kappa |
The thermal conductivity kappa in hyperbolic thermal conduction. | |
double precision, public, protected | he_abundance =0.1d0 |
Helium abundance over Hydrogen. | |
double precision, public, protected | h_ion_fr =1d0 |
Ionization fraction of H H_ion_fr = H+/(H+ + H) | |
double precision, public, protected | he_ion_fr =1d0 |
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He) | |
double precision, public, protected | he_ion_fr2 =1d0 |
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+) | |
double precision, public, protected | rr =1d0 |
integer, public | equi_rho0_ = -1 |
equi vars indices in the stateequi_vars array | |
integer, public | equi_pe0_ = -1 |
integer, public, protected | mhd_n_tracer = 0 |
Number of tracer species. | |
integer, public, protected | rho_ |
Index of the density (in the w array) | |
integer, dimension(:), allocatable, public, protected | mom |
Indices of the momentum density. | |
integer, public, protected | c |
Indices of the momentum density for the form of better vectorization. | |
integer, public, protected | m |
integer, public, protected | c_ |
integer, public, protected | e_ |
Index of the energy density (-1 if not present) | |
integer, public, protected | b |
integer, public, protected | p_ |
Index of the gas pressure (-1 if not present) should equal e_. | |
integer, public, protected | q_ |
Index of the heat flux q. | |
integer, public, protected | psi_ |
Indices of the GLM psi. | |
integer, public, protected | te_ |
Indices of temperature. | |
integer, public, protected | tcoff_ |
Index of the cutoff temperature for the TRAC method. | |
integer, public, protected | tweight_ |
integer, dimension(:), allocatable, public, protected | tracer |
Indices of the tracers. | |
integer, dimension(2 *^nd), public, protected | boundary_divbfix_skip =0 |
To skip * layer of ghost cells during divB=0 fix for boundary. | |
logical, public, protected | mhd_energy = .true. |
Whether an energy equation is used. | |
logical, public, protected | mhd_thermal_conduction = .false. |
Whether thermal conduction is used. | |
logical, public, protected | mhd_radiative_cooling = .false. |
Whether radiative cooling is added. | |
logical, public, protected | mhd_hyperbolic_thermal_conduction = .false. |
Whether thermal conduction is used. | |
logical, public, protected | mhd_viscosity = .false. |
Whether viscosity is added. | |
logical, public, protected | mhd_gravity = .false. |
Whether gravity is added. | |
logical, public, protected | mhd_rotating_frame = .false. |
Whether rotating frame is activated. | |
logical, public, protected | mhd_hall = .false. |
Whether Hall-MHD is used. | |
logical, public, protected | mhd_ambipolar = .false. |
Whether Ambipolar term is used. | |
logical, public, protected | mhd_ambipolar_sts = .false. |
Whether Ambipolar term is implemented using supertimestepping. | |
logical, public, protected | mhd_ambipolar_exp = .false. |
Whether Ambipolar term is implemented explicitly. | |
logical, public, protected | mhd_particles = .false. |
Whether particles module is added. | |
logical, public, protected | mhd_magnetofriction = .false. |
Whether magnetofriction is added. | |
logical, public, protected | mhd_glm = .false. |
Whether GLM-MHD is used to control div B. | |
logical, public, protected | mhd_glm_extended = .true. |
Whether extended GLM-MHD is used with additional sources. | |
logical, public, protected | mhd_trac = .false. |
Whether TRAC method is used. | |
integer, public, protected | mhd_trac_type =1 |
Which TRAC method is used. | |
integer, public, protected | mhd_trac_finegrid =4 |
Distance between two adjacent traced magnetic field lines (in finest cell size) | |
logical, public, protected | mhd_internal_e = .false. |
Whether internal energy is solved instead of total energy. | |
logical, public, protected | mhd_hydrodynamic_e = .false. |
Whether hydrodynamic energy is solved instead of total energy. | |
logical, public, protected | source_split_divb = .false. |
Whether divB cleaning sources are added splitting from fluid solver. | |
logical, public, protected | mhd_semirelativistic = .false. |
Whether semirelativistic MHD equations (Gombosi 2002 JCP) are solved. | |
logical, public, protected | mhd_partial_ionization = .false. |
Whether plasma is partially ionized. | |
logical, public, protected | mhd_cak_force = .false. |
Whether CAK radiation line force is activated. | |
logical, public, protected | mhd_4th_order = .false. |
MHD fourth order. | |
logical, public | has_equi_rho0 = .false. |
whether split off equilibrium density | |
logical, public | has_equi_pe0 = .false. |
whether split off equilibrium thermal pressure | |
logical, public | mhd_equi_thermal = .false. |
logical, public, protected | mhd_dump_full_vars = .false. |
whether dump full variables (when splitting is used) in a separate dat file | |
integer, public, protected | mhd_divb_nth = 1 |
Whether divB is computed with a fourth order approximation. | |
logical, public | divbwave = .true. |
Add divB wave in Roe solver. | |
logical, public | clean_initial_divb = .false. |
clean initial divB | |
logical, public, protected | eq_state_units = .true. |
logical, dimension(2 *^nd), public, protected | boundary_divbfix =.true. |
To control divB=0 fix for boundary. | |
logical, public, protected | b0field_forcefree =.true. |
B0 field is force-free. | |
logical, public | partial_energy = .false. |
Whether an internal or hydrodynamic energy equation is used. | |
character(len=std_len), public, protected | typedivbfix = 'linde' |
Method type to clean divergence of B. | |
character(len=std_len), public, protected | type_ct = 'uct_contact' |
Method type of constrained transport. | |
type(tc_fluid), allocatable, public | tc_fl |
type of fluid for thermal conduction | |
type(te_fluid), allocatable, public | te_fl_mhd |
type of fluid for thermal emission synthesis | |
type(rc_fluid), allocatable, public | rc_fl |
type of fluid for radiative cooling | |
procedure(mask_subroutine), pointer, public | usr_mask_ambipolar => null() |
procedure(sub_convert), pointer, public | mhd_to_primitive => null() |
procedure(sub_convert), pointer, public | mhd_to_conserved => null() |
procedure(sub_get_pthermal), pointer, public | mhd_get_pthermal => null() |
procedure(sub_get_pthermal), pointer, public | mhd_get_temperature => null() |
Magneto-hydrodynamics module.
subroutine, public mod_mhd_phys::b_from_vector_potential | ( | integer, intent(in) | ixis, |
integer, intent(in) | l, | ||
integer, intent(in) | ixi, | ||
l, | |||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixis^s,1:nws), intent(inout) | ws, | ||
double precision, dimension(ixi^s,1:ndim), intent(in) | x | ||
) |
calculate magnetic field from vector potential
Definition at line 7437 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::get_current | ( | double precision, dimension(ixi^s,1:nw), intent(in) | w, |
integer, intent(in) | ixi, | ||
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
integer, intent(out) | idirmin, | ||
double precision, dimension(ixi^s,7-2*ndir:3) | current | ||
) |
Calculate idirmin and the idirmin:3 components of the common current array make sure that dxlevel(^D) is set correctly.
Definition at line 5467 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::get_normalized_divb | ( | double precision, dimension(ixi^s,1:nw), intent(in) | w, |
integer, intent(in) | ixi, | ||
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixi^s) | divb | ||
) |
get dimensionless div B = |divB| * volume / area / |B|
Definition at line 5433 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::mhd_clean_divb_multigrid | ( | double precision, intent(in) | qdt, |
double precision, intent(in) | qt, | ||
logical, intent(inout) | active | ||
) |
[in] | qdt | Current time step |
[in] | qt | Current time |
[in,out] | active | Output if the source is active |
Definition at line 6563 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::mhd_e_to_ei | ( | integer, intent(in) | ixi, |
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixi^s, nw), intent(inout) | w, | ||
double precision, dimension(ixi^s, 1:ndim), intent(in) | x | ||
) |
Transform total energy to internal energy.
Definition at line 1944 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::mhd_ei_to_e | ( | integer, intent(in) | ixi, |
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixi^s, nw), intent(inout) | w, | ||
double precision, dimension(ixi^s, 1:ndim), intent(in) | x | ||
) |
Transform internal energy to total energy.
Definition at line 1878 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::mhd_face_to_center | ( | integer, intent(in) | ixo, |
integer, intent(in) | l, | ||
type(state) | s | ||
) |
calculate cell-center values from face-center values
Definition at line 7372 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::mhd_get_rho | ( | double precision, dimension(ixi^s,1:nw), intent(in) | w, |
double precision, dimension(ixi^s,1:ndim), intent(in) | x, | ||
integer, intent(in) | ixi, | ||
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixi^s), intent(out) | rho | ||
) |
subroutine, public mod_mhd_phys::mhd_get_v | ( | double precision, dimension(ixi^s,nw), intent(in) | w, |
double precision, dimension(ixi^s,1:ndim), intent(in) | x, | ||
integer, intent(in) | ixi, | ||
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixi^s,ndir), intent(out) | v | ||
) |
Calculate v vector.
Definition at line 2388 of file mod_mhd_phys.t.
double precision function, dimension(ixo^s), public mod_mhd_phys::mhd_mag_en_all | ( | double precision, dimension(ixi^s, nw), intent(in) | w, |
integer, intent(in) | ixi, | ||
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l | |||
) |
Compute 2 times total magnetic energy.
Definition at line 6032 of file mod_mhd_phys.t.
subroutine, public mod_mhd_phys::mhd_phys_init |
subroutine, public mod_mhd_phys::multiplyambicoef | ( | integer, intent(in) | ixi, |
integer, intent(in) | l, | ||
integer, intent(in) | ixo, | ||
l, | |||
double precision, dimension(ixi^s), intent(inout) | res, | ||
double precision, dimension(ixi^s,1:nw), intent(in) | w, | ||
double precision, dimension(ixi^s,1:ndim), intent(in) | x | ||
) |
multiply res by the ambipolar coefficient The ambipolar coefficient is calculated as -mhd_eta_ambi/rho^2 The user may mask its value in the user file by implemneting usr_mask_ambipolar subroutine
Definition at line 4339 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::b |
Definition at line 75 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::b0field_forcefree =.true. |
B0 field is force-free.
Definition at line 180 of file mod_mhd_phys.t.
logical, dimension(2*^nd), public, protected mod_mhd_phys::boundary_divbfix =.true. |
To control divB=0 fix for boundary.
Definition at line 178 of file mod_mhd_phys.t.
integer, dimension(2*^nd), public, protected mod_mhd_phys::boundary_divbfix_skip =0 |
To skip * layer of ghost cells during divB=0 fix for boundary.
Definition at line 94 of file mod_mhd_phys.t.
integer public protected mod_mhd_phys::c |
Indices of the momentum density for the form of better vectorization.
Definition at line 71 of file mod_mhd_phys.t.
integer public protected mod_mhd_phys::c_ |
Definition at line 71 of file mod_mhd_phys.t.
logical, public mod_mhd_phys::clean_initial_divb = .false. |
clean initial divB
Definition at line 172 of file mod_mhd_phys.t.
logical, public mod_mhd_phys::divbwave = .true. |
Add divB wave in Roe solver.
Definition at line 170 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::e_ |
Index of the energy density (-1 if not present)
Definition at line 73 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::eq_state_units = .true. |
Definition at line 176 of file mod_mhd_phys.t.
integer, public mod_mhd_phys::equi_pe0_ = -1 |
Definition at line 63 of file mod_mhd_phys.t.
integer, public mod_mhd_phys::equi_rho0_ = -1 |
equi vars indices in the stateequi_vars array
Definition at line 62 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::h_ion_fr =1d0 |
Ionization fraction of H H_ion_fr = H+/(H+ + H)
Definition at line 46 of file mod_mhd_phys.t.
logical, public mod_mhd_phys::has_equi_pe0 = .false. |
whether split off equilibrium thermal pressure
Definition at line 161 of file mod_mhd_phys.t.
logical, public mod_mhd_phys::has_equi_rho0 = .false. |
whether split off equilibrium density
Definition at line 159 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::he_abundance =0.1d0 |
Helium abundance over Hydrogen.
Definition at line 43 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::he_ion_fr =1d0 |
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
Definition at line 49 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::he_ion_fr2 =1d0 |
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
Definition at line 52 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::hypertc_kappa |
The thermal conductivity kappa in hyperbolic thermal conduction.
Definition at line 39 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::m |
Definition at line 71 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_4th_order = .false. |
MHD fourth order.
Definition at line 157 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_adiab = 1.0d0 |
The adiabatic constant.
Definition at line 20 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_ambipolar = .false. |
Whether Ambipolar term is used.
Definition at line 123 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_ambipolar_exp = .false. |
Whether Ambipolar term is implemented explicitly.
Definition at line 127 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_ambipolar_sts = .false. |
Whether Ambipolar term is implemented using supertimestepping.
Definition at line 125 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_cak_force = .false. |
Whether CAK radiation line force is activated.
Definition at line 155 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::mhd_divb_nth = 1 |
Whether divB is computed with a fourth order approximation.
Definition at line 166 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_dump_full_vars = .false. |
whether dump full variables (when splitting is used) in a separate dat file
Definition at line 164 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_energy = .true. |
Whether an energy equation is used.
Definition at line 107 of file mod_mhd_phys.t.
logical, public mod_mhd_phys::mhd_equi_thermal = .false. |
Definition at line 162 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_eta = 0.0d0 |
The MHD resistivity.
Definition at line 22 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_eta_ambi = 0.0d0 |
The MHD ambipolar coefficient.
Definition at line 28 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_eta_hyper = 0.0d0 |
The MHD hyper-resistivity.
Definition at line 24 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_etah = 0.0d0 |
Hall resistivity.
Definition at line 26 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_gamma = 5.d0/3.0d0 |
The adiabatic index.
Definition at line 18 of file mod_mhd_phys.t.
procedure(sub_get_pthermal), pointer, public mod_mhd_phys::mhd_get_pthermal => null() |
Definition at line 219 of file mod_mhd_phys.t.
procedure(sub_get_pthermal), pointer, public mod_mhd_phys::mhd_get_temperature => null() |
Definition at line 221 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_glm = .false. |
Whether GLM-MHD is used to control div B.
Definition at line 133 of file mod_mhd_phys.t.
double precision, public mod_mhd_phys::mhd_glm_alpha = 0.5d0 |
GLM-MHD parameter: ratio of the diffusive and advective time scales for div b taking values within [0, 1].
Definition at line 35 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_glm_extended = .true. |
Whether extended GLM-MHD is used with additional sources.
Definition at line 135 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_gravity = .false. |
Whether gravity is added.
Definition at line 117 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_hall = .false. |
Whether Hall-MHD is used.
Definition at line 121 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_hydrodynamic_e = .false. |
Whether hydrodynamic energy is solved instead of total energy.
Definition at line 146 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_hyperbolic_thermal_conduction = .false. |
Whether thermal conduction is used.
Definition at line 113 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_internal_e = .false. |
Whether internal energy is solved instead of total energy.
Definition at line 143 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_magnetofriction = .false. |
Whether magnetofriction is added.
Definition at line 131 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::mhd_n_tracer = 0 |
Number of tracer species.
Definition at line 65 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_partial_ionization = .false. |
Whether plasma is partially ionized.
Definition at line 153 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_particles = .false. |
Whether particles module is added.
Definition at line 129 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_radiative_cooling = .false. |
Whether radiative cooling is added.
Definition at line 111 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::mhd_reduced_c = 0.02d0*const_c |
Reduced speed of light for semirelativistic MHD: 2% of light speed.
Definition at line 37 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_rotating_frame = .false. |
Whether rotating frame is activated.
Definition at line 119 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_semirelativistic = .false. |
Whether semirelativistic MHD equations (Gombosi 2002 JCP) are solved.
Definition at line 151 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_thermal_conduction = .false. |
Whether thermal conduction is used.
Definition at line 109 of file mod_mhd_phys.t.
procedure(sub_convert), pointer, public mod_mhd_phys::mhd_to_conserved => null() |
Definition at line 217 of file mod_mhd_phys.t.
procedure(sub_convert), pointer, public mod_mhd_phys::mhd_to_primitive => null() |
Definition at line 216 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_trac = .false. |
Whether TRAC method is used.
Definition at line 137 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::mhd_trac_finegrid =4 |
Distance between two adjacent traced magnetic field lines (in finest cell size)
Definition at line 141 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::mhd_trac_mask = 0.d0 |
Height of the mask used in the TRAC method.
Definition at line 32 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::mhd_trac_type =1 |
Which TRAC method is used.
Definition at line 139 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::mhd_viscosity = .false. |
Whether viscosity is added.
Definition at line 115 of file mod_mhd_phys.t.
integer, dimension(:), allocatable, public, protected mod_mhd_phys::mom |
Indices of the momentum density.
Definition at line 69 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::p_ |
Index of the gas pressure (-1 if not present) should equal e_.
Definition at line 77 of file mod_mhd_phys.t.
logical, public mod_mhd_phys::partial_energy = .false. |
Whether an internal or hydrodynamic energy equation is used.
Definition at line 184 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::psi_ |
Indices of the GLM psi.
Definition at line 81 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::q_ |
Index of the heat flux q.
Definition at line 79 of file mod_mhd_phys.t.
type(rc_fluid), allocatable, public mod_mhd_phys::rc_fl |
type of fluid for radiative cooling
Definition at line 200 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::rho_ |
Index of the density (in the w array)
Definition at line 67 of file mod_mhd_phys.t.
double precision, public, protected mod_mhd_phys::rr =1d0 |
Definition at line 56 of file mod_mhd_phys.t.
logical, public, protected mod_mhd_phys::source_split_divb = .false. |
Whether divB cleaning sources are added splitting from fluid solver.
Definition at line 148 of file mod_mhd_phys.t.
type(tc_fluid), allocatable, public mod_mhd_phys::tc_fl |
type of fluid for thermal conduction
Definition at line 196 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::tcoff_ |
Index of the cutoff temperature for the TRAC method.
Definition at line 85 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::te_ |
Indices of temperature.
Definition at line 83 of file mod_mhd_phys.t.
type(te_fluid), allocatable, public mod_mhd_phys::te_fl_mhd |
type of fluid for thermal emission synthesis
Definition at line 198 of file mod_mhd_phys.t.
integer, dimension(:), allocatable, public, protected mod_mhd_phys::tracer |
Indices of the tracers.
Definition at line 88 of file mod_mhd_phys.t.
integer, public, protected mod_mhd_phys::tweight_ |
Definition at line 86 of file mod_mhd_phys.t.
character(len=std_len), public, protected mod_mhd_phys::type_ct = 'uct_contact' |
Method type of constrained transport.
Definition at line 192 of file mod_mhd_phys.t.
character(len=std_len), public, protected mod_mhd_phys::typedivbfix = 'linde' |
Method type to clean divergence of B.
Definition at line 190 of file mod_mhd_phys.t.
procedure(mask_subroutine), pointer, public mod_mhd_phys::usr_mask_ambipolar => null() |
Definition at line 215 of file mod_mhd_phys.t.