The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
|
Module for modeling extended MHD with a mixed conforming/lagrange basis set.
The equations currently solved are the equations of extended xMHD with scalar temperature and pressure.
\[ \frac{\partial n}{\partial t} + \nabla \cdot \left( n u \right) = D \nabla^2 n \]
\[ \rho \left[ \frac{\partial u}{\partial t} + u \cdot \nabla u \right] = J \times B - \nabla nk(T_e + T_i) - \nabla \cdot \Pi \]
\[ \frac{\partial B}{\partial t} = - \nabla \times \left( -u \times B + \eta J + \frac{1}{ne} \left( J \times B - \nabla nkT_e \right) + \frac{m_e}{ne^2} \frac{\partial J}{\partial t} \right) \]
\[ \frac{n}{\gamma-1} \left[ \frac{\partial T_i}{\partial t} + u \cdot \nabla T_i \right] = -nkT_i \nabla \cdot u - \nabla \cdot q_i + Q_i \]
For two temperature
\[ \frac{n}{\gamma-1} \left[ \frac{\partial T_e}{\partial t} + u_e \cdot \nabla T_e \right] = -nkT_e \nabla \cdot u_e - \nabla \cdot q_e + Q_e \]
\[ u_e = u - \frac{J}{ne} \]
Thermal transport
\[ q_s = - n \left[ \chi_{\parallel,s} \hat{b} \hat{b} + \chi_{\perp,s} \left( I - \hat{b} \hat{b} \right) \right] \cdot \nabla T_s \]
With xmhd_brag=T
and single temperature ( \( T_e = T_i \))
\[ \chi_{\parallel,i} = \chi_{\parallel,e}, \chi_{\perp,i} = min(\chi_{\perp,i},\chi_{\parallel,i}) \]
Heat sources
For single temperature ( \( T_e = T_i \))
\[ Q_i = \left[ \eta J^2 - \left( \nabla u \right)^T : \Pi \right] / 2 \]
For two temperature
\[ Q_i = - \left( \nabla u \right)^T : \Pi, Q_e = \eta J^2 \]
Viscosity
\[ W = \left( \nabla u + (\nabla u)^T - \frac{2}{3} I \nabla \cdot u \right) \]
With ‘visc_type='kin’`
\[ \Pi = - \nu \nabla u \]
With ‘visc_type='iso’`
\[ \Pi = - \nu W \]
With ‘visc_type='ani’`
\[ \Pi = - \left[ \nu_{\parallel} \hat{b} \hat{b} + \nu_{\perp} \left( I - \hat{b} \hat{b} \right) \right] \cdot W \]
Runtime options are available when using the driver (xmhd_run) and post-processing (xmhd_plot) routines, see their documentation for available options.
Data Types | |
type | ml_xmhd_vecspace |
Needs docs. More... | |
type | oft_xmhd_driver |
Forcing object class for xMHD. More... | |
type | oft_xmhd_errmatrix |
NL metric matrix for xMHD. More... | |
type | oft_xmhd_massmatrix |
Mass metric matrix for Reduced MHD. More... | |
type | oft_xmhd_probe |
Probe object class for xMHD. More... | |
interface | xmhd_driver_apply |
Modify solution vectors to apply forcing. More... | |
type | xmhd_interp |
Interpolate a xMHD operator fields. More... | |
type | xmhd_interp_cache |
Interpolate a xMHD operator fields. More... | |
type | xmhd_loc_values |
Unpack object for local fields evaluated by xmhd_interp. More... | |
type | xmhd_ops |
xMHD operator container More... | |
interface | xmhd_probe_apply |
Sample probe signals from solution. More... | |
type | xmhd_sub_fields |
Object for sub-fields used by oft_xmhd_push and oft_xmhd_pop. More... | |
Functions/Subroutines | |
subroutine | oft_xmhd_create (new, level, cache, native) |
Create new xMHD solution vector. | |
subroutine | oft_xmhd_inject (self, afine, acors) |
Inject a fine level xMHD solution to the next coarsest level. | |
subroutine | oft_xmhd_interp (self, acors, afine) |
Interpolate a coarse level xMHD solution to the next finest level. | |
subroutine, public | oft_xmhd_pop (u, sub_fields) |
Extract subfields from xMHD solution vector. | |
subroutine, public | oft_xmhd_push (sub_fields, u) |
Replace subfields in xMHD solution vector with specified fields. | |
subroutine | oft_xmhd_rst_load (u, filename, path, t, dt, rst_version_out) |
Load xMHD solution state from a restart file. | |
subroutine | oft_xmhd_rst_save (u, t, dt, filename, path) |
Save xMHD solution state to a restart file. | |
subroutine | xmhd_alloc_ops |
Setup and allocate operators used in xMHD advance. | |
subroutine | xmhd_build_ops (fullin) |
Update Jacobian matrix with new fields. | |
subroutine | xmhd_diag (u, diag_vals, ueq) |
Compute diagnostic values from fields. | |
subroutine | xmhd_driver_rst_dummy (self, rst_file) |
Save or load driver info from restart file. | |
subroutine | xmhd_errmatrix_apply (self, a, b) |
Compute the NL metric for solution field. | |
subroutine | xmhd_interp_apply (self, cell, f, gop, val) |
Reconstruct xMHD solution fields. | |
subroutine | xmhd_interp_delete (self) |
Destroy interpolator for xMHD solution fields. | |
subroutine | xmhd_interp_setup (self, mesh) |
Setup interpolator for xMHD solution fields. | |
subroutine | xmhd_interp_unpack (vals, vloc) |
Unpack interpolation fields into xMHD field structure. | |
subroutine, public | xmhd_lin_run (equil_fields, pert_fields, escale) |
Main driver subroutine for extended MHD time advance. | |
subroutine | xmhd_massmatrix_apply (self, a, b) |
Compute the mass matrix for solution field. | |
subroutine | xmhd_mfnk_update (uin) |
Update Jacobian matrices on all levels with new fields. | |
subroutine, public | xmhd_plot (probes) |
Plotting subroutine for xMHD post-processing. | |
subroutine | xmhd_probe_flush (self) |
Flush internal write buffers for probe object. | |
subroutine | xmhd_profile (u) |
Simple subroutine to compute timing of different solve phases. | |
subroutine | xmhd_read_settings (dt, lin_tol, nl_tol, rst_ind, nsteps, rst_freq, nclean, maxextrap, ittarget, nl_update) |
Read settings for extended MHD model from input files. | |
subroutine, public | xmhd_run (initial_fields, driver, probes, profile_only) |
Main driver subroutine for extended MHD time advance. | |
subroutine | xmhd_set_bc |
Set BC flags. | |
subroutine | xmhd_set_level (level) |
Set the current level for xMHD model. | |
subroutine | xmhd_set_ops (uin) |
Update Jacobian matrices on all levels with new solution. | |
subroutine | xmhd_setup_regions () |
Setup material regions from XML input file. | |
subroutine | xmhd_setup_rep |
Setup composite FE representation and ML environment. | |
pure real(r8) function | xmhd_upwind_weight (v_mag, he, nu) |
Evalute upwinding parameter. | |
Variables | |
character(len=2) | bbc = 'bc' |
Magnetic field BC ('bc','ic') | |
real(r8) | d2_dens = -1.d0 |
Particle hyper-diffusivity. | |
real(r8) | d_dens = 1.E0_r8 |
Particle diffusivity. | |
real(r8), public | den_floor = -1.d0 |
Density floor. | |
real(r8), public | den_scale = -1.d0 |
Density scale factor. | |
real(r8) | dt_scale = 1._r8 |
Time step scaling factor. | |
real(r8) | eta = 1.E-2_r8 |
Resistivity \( \eta / \mu_0 \). | |
real(r8) | eta_hyper = -1.d0 |
Hyper resistivity. | |
real(r8), dimension(:), allocatable | eta_reg |
Resistivity scale factors. | |
real(r8) | eta_temp = -1._r8 |
Reference temperature for resistivity. | |
real(r8), public | g_accel = 0.d0 |
Gravity (oriented in the -Z direction) | |
integer(i4) | j2_ind = -1 |
Index of hyper-res aux field. | |
real(r8), public | j2_scale = -1.d0 |
Hyper-resistivity aux variable scale factor. | |
real(r8) | jac_dt = 1.E-3_r8 |
Time step. | |
real(r8) | k_boltz = elec_charge |
Boltzmann constant. | |
real(r8) | kappa_par = 1.d0 |
Parallel thermal conductivity factor. | |
real(r8) | kappa_perp = 1.d0 |
Perpendicular thermal conductivity factor. | |
real(r8) | m_ion = -1.d0 |
Ion mass in SI units. | |
real(r8) | me_factor = 100.E0_r8 |
Artificial electron mass factor. | |
class(oft_mesh), pointer | mesh |
type(oft_mf_matrix), target | mfmat |
Matrix free operator. | |
class(multigrid_mesh), pointer | mg_mesh |
type(oft_matrix_ptr), dimension(:), pointer | ml_int => NULL() |
MG interpolation operators. | |
type(oft_matrix_ptr), dimension(:), pointer | ml_j => NULL() |
MG Jacobian operators. | |
type(oft_ml_fem_comp_type), target | ml_xmhd_rep |
ML container for field representation. | |
real(r8) | mu_ion = 2.d0 |
Ion mass in atomic units. | |
integer(i4) | n2_ind = -1 |
Index of hyper-diff aux field. | |
real(r8), public | n2_scale = -1.d0 |
Hyper-diffusivity aux variable scale factor. | |
character(len=1) | nbc = 'd' |
Density BC ('d','n') | |
real(r8), dimension(:,:), allocatable | neg_flag |
real(r8), dimension(:,:), allocatable | neg_source |
real(r8) | nu_par = 0.d0 |
Fluid (parallel) viscosity (note: dynamic viscosity is nu=nu_par*den_scale/n0) | |
real(r8) | nu_perp = 0.d0 |
Fluid perpendicular viscosity if anisotropic, as nu_par. | |
integer(i4), dimension(fem_max_levels) | nu_xmhd = 1 |
Number of smoother iterations. | |
class(oft_hcurl_fem), pointer | oft_hcurl => NULL() |
class(oft_h1_fem), pointer | oft_hgrad => NULL() |
class(oft_scalar_fem), pointer | oft_lagrange => NULL() |
type(xmhd_ops), pointer | oft_xmhd_ops => NULL() |
Operator container. | |
type(xmhd_ops), dimension(:), pointer | oft_xmhd_ops_ml => NULL() |
MG operator container. | |
procedure(oft_1d_func), pointer, public | res_profile => NULL() |
Resistivity profile. | |
logical, dimension(:), allocatable | solid_cell |
Solid region flag. | |
character(len=1) | tbc = 'd' |
Temperature BC ('d','n') | |
real(r8) | te_factor = 1.d0 |
Electron temperature factor (single temperature) | |
integer(i4) | te_ind = -1 |
Index of electron temperature field. | |
real(r8), public | temp_floor = -1.d0 |
Temperature floor. | |
real(r8) | temp_gamma = 5.d0/3.d0 |
Ratio of specific heats. | |
character(len=4) | vbc = 'all' |
Velocity BC ('all','norm','tang') | |
integer(i4) | vbc_type = 1 |
Velocity BC type. | |
real(r8), public | vel_scale = 1.d3 |
Velocity scale factor. | |
integer(i4) | visc_itype = 1 |
Viscosity type flag (1->kin,2->iso,3->ani) | |
character(len=3) | visc_type = 'kin' |
Viscosity type ('kin','iso','ani') | |
logical, public | xmhd_adv_b = .TRUE. |
Advance magnetic field. | |
logical | xmhd_adv_den = .TRUE. |
Advance density. | |
logical | xmhd_adv_temp = .TRUE. |
Advance temperature. | |
logical | xmhd_advec = .TRUE. |
Include fluid advection. | |
integer(i4) | xmhd_blevel = 0 |
Highest level on base meshes. | |
logical, public | xmhd_bnorm_force = .TRUE. |
Force B-norm to zero in cleaning. | |
logical | xmhd_brag = .FALSE. |
Braginskii thermal conduction. | |
logical | xmhd_diss_centered = .FALSE. |
Time-center dissipation terms. | |
real(r8) | xmhd_eps = 1.d-10 |
Epsilon for magnetic field normalization. | |
logical | xmhd_hall = .FALSE. |
Include Hall terms. | |
real(r8), dimension(:,:), pointer, contiguous | xmhd_hcurl_cop => NULL() |
real(r8), dimension(:,:), pointer, contiguous | xmhd_hcurl_rop => NULL() |
real(r8), dimension(:,:), pointer, contiguous | xmhd_hgrad_rop => NULL() |
logical | xmhd_jcb = .TRUE. |
Include JxB force on fluid. | |
real(r8), dimension(:,:), pointer, contiguous | xmhd_lag_gop => NULL() |
real(r8), dimension(:), pointer, contiguous | xmhd_lag_rop => NULL() |
integer(i4) | xmhd_lev = 1 |
Active FE level. | |
integer(i4) | xmhd_level = 1 |
Active FE level. | |
integer(i4) | xmhd_lin_level = 1 |
Highest linear element level. | |
logical | xmhd_linear = .FALSE. |
Compute linearized matrix. | |
integer(i4), dimension(:,:), allocatable | xmhd_mat_mask |
Matrix block mask. | |
logical | xmhd_mfnk = .FALSE. |
Use Jacobian free non-linear advance. | |
integer(i4), public | xmhd_minlev = -1 |
Lowest MG level. | |
type(oft_ml_fem_type), pointer, public | xmhd_ml_h1 => NULL() |
type(oft_ml_fem_type), pointer, public | xmhd_ml_h1grad => NULL() |
type(oft_ml_fem_type), pointer, public | xmhd_ml_hcurl => NULL() |
type(oft_ml_fem_comp_type), pointer, public | xmhd_ml_hcurl_grad => NULL() |
type(oft_ml_fem_type), pointer, public | xmhd_ml_lagrange => NULL() |
type(ml_xmhd_vecspace), target | xmhd_ml_vecspace |
type(oft_ml_fem_comp_type), pointer, public | xmhd_ml_vlagrange => NULL() |
logical | xmhd_monitor_div = .FALSE. |
Monitor div error? | |
integer(i4) | xmhd_nlevels = 1 |
Number of total levels. | |
integer(i4), dimension(fem_max_levels) | xmhd_nparts = 1 |
Number of local partitions for preconditioning. | |
logical | xmhd_ohmic = .FALSE. |
Resistive heating. | |
integer(i4) | xmhd_opcount = 1 |
Number of time steps since preconditioner update. | |
real(r8) | xmhd_opdt = 1.d0 |
Time step for current preconditioner object. | |
class(oft_solver), pointer | xmhd_pre => NULL() |
Preconditioner object. | |
type(xml_node), pointer | xmhd_pre_node => NULL() |
preconditioner XML node | |
integer(i4) | xmhd_prefreq = 30 |
Desired update frequency for preconditioner. | |
type(oft_fem_comp_type), pointer | xmhd_rep => NULL() |
Active field representation. | |
type(xml_node), pointer | xmhd_root_node => NULL() |
xMHD XML node | |
integer(i4), parameter | xmhd_rst_version = 3 |
Restart file version number. | |
logical | xmhd_rw = .FALSE. |
Resistive wall present. | |
logical | xmhd_skip_update = .FALSE. |
Internal flag for skipping Jac update. | |
integer(i4), public | xmhd_taxis = 3 |
Axis for toroidal flux and current. | |
logical | xmhd_therm_equil = .FALSE. |
Include cross-species thermalization. | |
logical | xmhd_two_temp = .FALSE. |
Run 2 Temperature model. | |
logical | xmhd_upwind = .TRUE. |
Use upwinding. | |
logical | xmhd_vbcdir = .FALSE. |
Velocity BC type flag. | |
logical | xmhd_visc_heat = .FALSE. |
Viscous heating. | |
|
private |
Create new xMHD solution vector.
[out] | new | Vector to create |
[in] | level | FE level (optional) |
[in] | cache | Allow caching? (optional) |
[in] | native | Force native representation? (optional) |
|
private |
Inject a fine level xMHD solution to the next coarsest level.
[in,out] | afine | Fine level vector to inject |
[in,out] | acors | Injected solution on coarse level |
|
private |
Interpolate a coarse level xMHD solution to the next finest level.
[in,out] | acors | Coarse level vector to interpolate |
[in,out] | afine | Interpolated solution on fine level |
subroutine, public oft_xmhd_pop | ( | class(oft_vector), intent(inout) | u, |
type(xmhd_sub_fields), intent(inout) | sub_fields | ||
) |
Extract subfields from xMHD solution vector.
[in,out] | u | Source xMHD vector |
[in,out] | sub_fields | Destination fields |
subroutine, public oft_xmhd_push | ( | type(xmhd_sub_fields), intent(inout) | sub_fields, |
class(oft_vector), intent(inout) | u | ||
) |
Replace subfields in xMHD solution vector with specified fields.
[in,out] | sub_fields | Source fields |
[in,out] | u | Destination xMHD vector |
|
private |
Load xMHD solution state from a restart file.
[in,out] | u | Solution to load |
[in] | filename | Name of restart file |
[in] | path | Path to store solution vector in file |
[out] | t | Time of loaded solution |
[out] | dt | Timestep at loaded time |
[out] | rst_version_out | Version number |
|
private |
Save xMHD solution state to a restart file.
[in,out] | u | Solution to save |
[in] | t | Current solution time |
[in] | dt | Current timestep |
[in] | filename | Name of restart file |
[in] | path | Path to store solution vector in file |
|
private |
Setup and allocate operators used in xMHD advance.
|
private |
Update Jacobian matrix with new fields.
[in,out] | fullin | Interpolator to evaluate fields |
|
private |
Compute diagnostic values from fields.
[in,out] | u | Field for diagnostic evaluation |
[out] | diag_vals | Resulting diagnostic values [7] |
[in,out] | ueq | Equilibrium field (linear calculation) |
|
private |
Save or load driver info from restart file.
[in,out] | self | Forcing object |
[in] | rst_file | Name of restart file |
|
private |
Compute the NL metric for solution field.
b = F(a)
[in,out] | self | Error matrix object |
[in,out] | a | Source field |
[in,out] | b | Result of metric function |
|
private |
Reconstruct xMHD solution fields.
[in,out] | self | Interpolation object |
[in] | cell | Cell for interpolation |
[in] | f | Position in cell in logical coord [4] |
[in] | gop | Logical gradient vectors at f [3,4] |
[out] | val | Reconstructed field at f [30] |
|
private |
Destroy interpolator for xMHD solution fields.
[in,out] | self | Interpolation object |
|
private |
Setup interpolator for xMHD solution fields.
[in,out] | self | Interpolation object |
|
private |
Unpack interpolation fields into xMHD field structure.
[in] | vals | Packed fields [40] |
[out] | vloc | Unpacked fields |
subroutine, public xmhd_lin_run | ( | type(xmhd_sub_fields), intent(inout) | equil_fields, |
type(xmhd_sub_fields), intent(inout) | pert_fields, | ||
real(r8), intent(in), optional | escale | ||
) |
Main driver subroutine for extended MHD time advance.
Runtime options are set in the main input file using the group xmhd_options
group, see xmhd_run.
[in,out] | equil_fields | Equilibrium fields |
[in,out] | pert_fields | Perurbed fields |
[in] | escale | Desired total energy to be used for rescaling (optional) |
|
private |
Compute the mass matrix for solution field.
b = M(a)
[in,out] | self | Mass matrix object |
[in,out] | a | Source field |
[in,out] | b | Result of metric function |
|
private |
Update Jacobian matrices on all levels with new fields.
[in,out] | uin | Current field |
subroutine, public xmhd_plot | ( | class(oft_xmhd_probe), intent(inout), optional | probes | ) |
Plotting subroutine for xMHD post-processing.
Runtime options are set in the main input file using the group xmhd_plot_options
. Additionally, options from the xmhd_options
group are also read to get run parameters, see xmhd_run for a list of options.
Option group: xmhd_plot_options
Option | Description | Type [dim] |
---|---|---|
probe_only=F | Sample probes but do not create plots | bool |
plot_div=F | Plot \( \nabla \cdot B \) field | bool |
file_list="none" | List of history files to read during post-processing (w/o file ext) | str(40) |
t0=0 | Starting time for sampling | real |
t1=1 | End time for sampling | real |
dt=1E-6 | Time between sampling points | real |
rst_start=0 | First restart file to read | int |
rst_end=2000 | Last restart file to read | int |
[in,out] | probes | Probe object (optional) |
|
private |
Flush internal write buffers for probe object.
[in,out] | self | Probe object |
|
private |
Simple subroutine to compute timing of different solve phases.
|
private |
Read settings for extended MHD model from input files.
Runtime options are set in the main input file using the xmhd_options
group.
Option group: xmhd_options
Option | Description | Type [dim] |
---|---|---|
xmhd_jcb=T | Apply magnetic force to fluid | bool |
xmhd_advec=T | Advect fluid | bool |
xmhd_adv_den=T | Advance density | bool |
xmhd_adv_temp=T | Advance temperature | bool |
xmhd_hall=F | Include Hall terms in Ohm's law | bool |
xmhd_ohmic=F | Include Ohmic heating | bool |
xmhd_visc_heat=F | Include Viscous heating | bool |
xmhd_brag=F | Use Braginskii thermal conduction coefficients | bool |
xmhd_upwind=T | Use upwinding for advection terms | bool |
xmhd_therm_equil=F | Include thermal equilibration between electrons/ions | bool |
xmhd_diss_centered=F | Use time-centered dissipation terms for (V, Ti, Te) | bool |
bbc="bc" | Magnetic field boundary condition ("bc", "ic") | str(2) |
vbc="all" | Velocity boundary condition ("all") | str(3) |
nbc="d" | Density boundary condition ("d", "n") | str(1) |
tbc="d" | Temperature boundary condition ("d", "n") | str(1) |
visc_type="kin" | Fluid viscosity type ("kin", "iso", "ani") | str(3) |
dt=1 | Maximum allowed time step | real |
eta=1 | Plasma resistivity, in diffusivity units | real |
eta_temp=-1 | Base resistivity temperature (neg to disable) | real |
eta_hyper=-1 | Plasma hyper-resistivity, in diffusivity units (neg to disable) | real |
nu_par=0 | Plasma (parallel) viscosity, in diffusivity units | real |
nu_perp=0 | Plasma perpendicular viscosity, in diffusivity units | real |
d_dens=1 | Plasma diffusivity, in diffusivity units | real |
d2_dens=-1 | Plasma hyper-diffusivity, in diffusivity units (neg to disable) | real |
kappa_par=1 | Plasma thermal conduction, in diffusivity units | real |
kappa_perp=1 | Plasma thermal conduction, in diffusivity units | real |
nsteps=2000 | Number of time steps to take | int |
rst_freq=10 | Frequency of dump files | int |
lin_tol=1E-9 | Linear solver tolerance | real |
nl_tol=1E-5 | Non-Linear solver tolerance | real |
nu_xmhd=1 | Number of smoother iterations for MG preconditioner | int [nlevels] |
xmhd_nparts=1 | Number of local partitions for factor based preconditioning | int |
nclean=500 | Frequency of divergence cleaning | int |
rst_ind=0 | Index of restart file to initialize from | int |
maxextrap=2 | Extrapolation order for NL initial guess | int |
ittarget=60 | Target number of linear iterations per time step | int |
mu_ion=2 | Ion mass in atomic units | real |
me_factor=100 | Electron mass factor | real |
te_factor=1 | Electron temperature factor (single temperature) | real |
xmhd_prefreq=30 | Frequency of full preconditioner update | int |
xmhd_monitor_div=F | Monitor divergence error | bool |
xmhd_mfnk=F | Use Matrix-Free NL solve | bool |
[out] | dt | Maximum timestep |
[out] | lin_tol | Linear solver tolerance |
[out] | nl_tol | Nonlinear solver tolerance |
[out] | rst_ind | Index of file for restart (rst_ind>0) |
[out] | nsteps | Number of timesteps |
[out] | rst_freq | Frequency to save restart files |
[out] | nclean | Frequency to clean divergence |
[out] | maxextrap | Extrapolation order for initial guess |
[out] | ittarget | Maximum number of linear iterations |
[out] | nl_update | Maximum number of linear iterations |
subroutine, public xmhd_run | ( | type(xmhd_sub_fields), intent(inout) | initial_fields, |
class(oft_xmhd_driver), intent(inout), optional | driver, | ||
class(oft_xmhd_probe), intent(inout), optional | probes, | ||
logical, intent(in), optional | profile_only | ||
) |
Main driver subroutine for extended MHD time advance.
Runtime options are set in the main input file using the group xmhd_options
group, see xmhd_read_settings.
[in,out] | initial_fields | Initial conditions |
[in,out] | driver | Forcing object |
[in,out] | probes | Probe object |
[in] | profile_only | Profile operator timing and stop? |
|
private |
Set BC flags.
|
private |
Set the current level for xMHD model.
[in] | level | Desired level |
|
private |
Update Jacobian matrices on all levels with new solution.
[in,out] | uin | Current solution |
|
private |
Setup material regions from XML input file.
|
private |
Setup composite FE representation and ML environment.
|
private |
Evalute upwinding parameter.
Upwinding parameter for inconsistent streamline upwinded Petrov-Galerkin method (distorted bases applied to advective terms only).
[in] | v_mag | Advection velocity [m/s] |
[in] | he | Element size [m] |
[in] | nu | Dissipation amplitude [m^2/s^2] |
|
private |
Magnetic field BC ('bc','ic')
|
private |
Particle hyper-diffusivity.
|
private |
Particle diffusivity.
real(r8), public den_floor = -1.d0 |
Density floor.
real(r8), public den_scale = -1.d0 |
Density scale factor.
|
private |
Time step scaling factor.
|
private |
Resistivity \( \eta / \mu_0 \).
|
private |
Hyper resistivity.
|
private |
Resistivity scale factors.
|
private |
Reference temperature for resistivity.
real(r8), public g_accel = 0.d0 |
Gravity (oriented in the -Z direction)
|
private |
Index of hyper-res aux field.
real(r8), public j2_scale = -1.d0 |
Hyper-resistivity aux variable scale factor.
|
private |
Time step.
|
private |
Boltzmann constant.
|
private |
Parallel thermal conductivity factor.
|
private |
Perpendicular thermal conductivity factor.
|
private |
Ion mass in SI units.
|
private |
Artificial electron mass factor.
|
private |
|
private |
Matrix free operator.
|
private |
|
private |
MG interpolation operators.
|
private |
MG Jacobian operators.
|
private |
ML container for field representation.
|
private |
Ion mass in atomic units.
|
private |
Index of hyper-diff aux field.
real(r8), public n2_scale = -1.d0 |
Hyper-diffusivity aux variable scale factor.
|
private |
Density BC ('d','n')
|
private |
|
private |
|
private |
Fluid (parallel) viscosity (note: dynamic viscosity is nu=nu_par*den_scale/n0)
|
private |
Fluid perpendicular viscosity if anisotropic, as nu_par.
|
private |
Number of smoother iterations.
|
private |
|
private |
|
private |
|
private |
Operator container.
|
private |
MG operator container.
procedure(oft_1d_func), pointer, public res_profile => NULL() |
Resistivity profile.
|
private |
Solid region flag.
|
private |
Temperature BC ('d','n')
|
private |
Electron temperature factor (single temperature)
|
private |
Index of electron temperature field.
real(r8), public temp_floor = -1.d0 |
Temperature floor.
|
private |
Ratio of specific heats.
|
private |
Velocity BC ('all','norm','tang')
|
private |
Velocity BC type.
real(r8), public vel_scale = 1.d3 |
Velocity scale factor.
|
private |
Viscosity type flag (1->kin,2->iso,3->ani)
|
private |
Viscosity type ('kin','iso','ani')
logical, public xmhd_adv_b = .TRUE. |
Advance magnetic field.
|
private |
Advance density.
|
private |
Advance temperature.
|
private |
Include fluid advection.
|
private |
Highest level on base meshes.
logical, public xmhd_bnorm_force = .TRUE. |
Force B-norm to zero in cleaning.
|
private |
Braginskii thermal conduction.
|
private |
Time-center dissipation terms.
|
private |
Epsilon for magnetic field normalization.
|
private |
Include Hall terms.
|
private |
|
private |
|
private |
|
private |
Include JxB force on fluid.
|
private |
|
private |
|
private |
Active FE level.
|
private |
Active FE level.
|
private |
Highest linear element level.
|
private |
Compute linearized matrix.
|
private |
Matrix block mask.
|
private |
Use Jacobian free non-linear advance.
integer(i4), public xmhd_minlev = -1 |
Lowest MG level.
type(oft_ml_fem_type), pointer, public xmhd_ml_h1 => NULL() |
type(oft_ml_fem_type), pointer, public xmhd_ml_h1grad => NULL() |
type(oft_ml_fem_type), pointer, public xmhd_ml_hcurl => NULL() |
type(oft_ml_fem_comp_type), pointer, public xmhd_ml_hcurl_grad => NULL() |
type(oft_ml_fem_type), pointer, public xmhd_ml_lagrange => NULL() |
|
private |
type(oft_ml_fem_comp_type), pointer, public xmhd_ml_vlagrange => NULL() |
|
private |
Monitor div error?
|
private |
Number of total levels.
|
private |
Number of local partitions for preconditioning.
|
private |
Resistive heating.
|
private |
Number of time steps since preconditioner update.
|
private |
Time step for current preconditioner object.
|
private |
Preconditioner object.
|
private |
preconditioner XML node
|
private |
Desired update frequency for preconditioner.
|
private |
Active field representation.
|
private |
xMHD XML node
integer(i4), parameter xmhd_rst_version = 3 |
Restart file version number.
|
private |
Resistive wall present.
|
private |
Internal flag for skipping Jac update.
integer(i4), public xmhd_taxis = 3 |
Axis for toroidal flux and current.
|
private |
Include cross-species thermalization.
|
private |
Run 2 Temperature model.
|
private |
Use upwinding.
|
private |
Velocity BC type flag.
|
private |
Viscous heating.