The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines | Variables
xmhd Module Reference

Detailed Description

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.

Note
Electron mass is artificially enhanced by a factor using me_factor
Authors
Chris Hansen
Date
September 2013

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(r8d2_dens = -1.d0
 Particle hyper-diffusivity.
 
real(r8d_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(r8dt_scale = 1._r8
 Time step scaling factor.
 
real(r8eta = 1.E-2_r8
 Resistivity \( \eta / \mu_0 \).
 
real(r8eta_hyper = -1.d0
 Hyper resistivity.
 
real(r8), dimension(:), allocatable eta_reg
 Resistivity scale factors.
 
real(r8eta_temp = -1._r8
 Reference temperature for resistivity.
 
real(r8), public g_accel = 0.d0
 Gravity (oriented in the -Z direction)
 
integer(i4j2_ind = -1
 Index of hyper-res aux field.
 
real(r8), public j2_scale = -1.d0
 Hyper-resistivity aux variable scale factor.
 
real(r8jac_dt = 1.E-3_r8
 Time step.
 
real(r8k_boltz = elec_charge
 Boltzmann constant.
 
real(r8kappa_par = 1.d0
 Parallel thermal conductivity factor.
 
real(r8kappa_perp = 1.d0
 Perpendicular thermal conductivity factor.
 
real(r8m_ion = -1.d0
 Ion mass in SI units.
 
real(r8me_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(r8mu_ion = 2.d0
 Ion mass in atomic units.
 
integer(i4n2_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(r8nu_par = 0.d0
 Fluid (parallel) viscosity (note: dynamic viscosity is nu=nu_par*den_scale/n0)
 
real(r8nu_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(r8te_factor = 1.d0
 Electron temperature factor (single temperature)
 
integer(i4te_ind = -1
 Index of electron temperature field.
 
real(r8), public temp_floor = -1.d0
 Temperature floor.
 
real(r8temp_gamma = 5.d0/3.d0
 Ratio of specific heats.
 
character(len=4) vbc = 'all'
 Velocity BC ('all','norm','tang')
 
integer(i4vbc_type = 1
 Velocity BC type.
 
real(r8), public vel_scale = 1.d3
 Velocity scale factor.
 
integer(i4visc_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(i4xmhd_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(r8xmhd_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(i4xmhd_lev = 1
 Active FE level.
 
integer(i4xmhd_level = 1
 Active FE level.
 
integer(i4xmhd_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(i4xmhd_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(i4xmhd_opcount = 1
 Number of time steps since preconditioner update.
 
real(r8xmhd_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(i4xmhd_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.
 

Function/Subroutine Documentation

◆ oft_xmhd_create()

subroutine oft_xmhd_create ( class(oft_vector), intent(out), pointer  new,
integer(i4), intent(in), optional  level,
logical, intent(in), optional  cache,
logical, intent(in), optional  native 
)
private

Create new xMHD solution vector.

Parameters
[out]newVector to create
[in]levelFE level (optional)
[in]cacheAllow caching? (optional)
[in]nativeForce native representation? (optional)

◆ oft_xmhd_inject()

subroutine oft_xmhd_inject ( class(ml_xmhd_vecspace), intent(inout)  self,
class(oft_vector), intent(inout)  afine,
class(oft_vector), intent(inout)  acors 
)
private

Inject a fine level xMHD solution to the next coarsest level.

Note
The global xmhd_level is decremented by one in this subroutine
Parameters
[in,out]afineFine level vector to inject
[in,out]acorsInjected solution on coarse level

◆ oft_xmhd_interp()

subroutine oft_xmhd_interp ( class(ml_xmhd_vecspace), intent(inout)  self,
class(oft_vector), intent(inout)  acors,
class(oft_vector), intent(inout)  afine 
)
private

Interpolate a coarse level xMHD solution to the next finest level.

Note
The global xmhd_level is incremented by one in this subroutine
Parameters
[in,out]acorsCoarse level vector to interpolate
[in,out]afineInterpolated solution on fine level

◆ oft_xmhd_pop()

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.

Parameters
[in,out]uSource xMHD vector
[in,out]sub_fieldsDestination fields

◆ oft_xmhd_push()

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.

Parameters
[in,out]sub_fieldsSource fields
[in,out]uDestination xMHD vector

◆ oft_xmhd_rst_load()

subroutine oft_xmhd_rst_load ( class(oft_vector), intent(inout), pointer  u,
character(len=*), intent(in)  filename,
character(len=*), intent(in)  path,
real(r8), intent(out), optional  t,
real(r8), intent(out), optional  dt,
integer(i4), intent(out), optional  rst_version_out 
)
private

Load xMHD solution state from a restart file.

Parameters
[in,out]uSolution to load
[in]filenameName of restart file
[in]pathPath to store solution vector in file
[out]tTime of loaded solution
[out]dtTimestep at loaded time
[out]rst_version_outVersion number

◆ oft_xmhd_rst_save()

subroutine oft_xmhd_rst_save ( class(oft_vector), intent(inout), pointer  u,
real(r8), intent(in)  t,
real(r8), intent(in)  dt,
character(len=*), intent(in)  filename,
character(len=*), intent(in)  path 
)
private

Save xMHD solution state to a restart file.

Parameters
[in,out]uSolution to save
[in]tCurrent solution time
[in]dtCurrent timestep
[in]filenameName of restart file
[in]pathPath to store solution vector in file

◆ xmhd_alloc_ops()

subroutine xmhd_alloc_ops
private

Setup and allocate operators used in xMHD advance.

◆ xmhd_build_ops()

subroutine xmhd_build_ops ( class(fem_interp), intent(inout)  fullin)
private

Update Jacobian matrix with new fields.

Parameters
[in,out]fullinInterpolator to evaluate fields

◆ xmhd_diag()

subroutine xmhd_diag ( class(oft_vector), intent(inout), target  u,
real(r8), dimension(7), intent(out)  diag_vals,
class(oft_vector), intent(inout), optional  ueq 
)
private

Compute diagnostic values from fields.

Parameters
[in,out]uField for diagnostic evaluation
[out]diag_valsResulting diagnostic values [7]
[in,out]ueqEquilibrium field (linear calculation)

◆ xmhd_driver_rst_dummy()

subroutine xmhd_driver_rst_dummy ( class(oft_xmhd_driver), intent(inout)  self,
character(*), intent(in)  rst_file 
)
private

Save or load driver info from restart file.

Parameters
[in,out]selfForcing object
[in]rst_fileName of restart file

◆ xmhd_errmatrix_apply()

subroutine xmhd_errmatrix_apply ( class(oft_xmhd_errmatrix), intent(inout)  self,
class(oft_vector), intent(inout), target  a,
class(oft_vector), intent(inout)  b 
)
private

Compute the NL metric for solution field.

b = F(a)

Parameters
[in,out]selfError matrix object
[in,out]aSource field
[in,out]bResult of metric function

◆ xmhd_interp_apply()

subroutine xmhd_interp_apply ( class(xmhd_interp), intent(inout)  self,
integer(i4), intent(in)  cell,
real(r8), dimension(:), intent(in)  f,
real(r8), dimension(3,4), intent(in)  gop,
real(r8), dimension(:), intent(out)  val 
)
private

Reconstruct xMHD solution fields.

Parameters
[in,out]selfInterpolation object
[in]cellCell for interpolation
[in]fPosition in cell in logical coord [4]
[in]gopLogical gradient vectors at f [3,4]
[out]valReconstructed field at f [30]

◆ xmhd_interp_delete()

subroutine xmhd_interp_delete ( class(xmhd_interp), intent(inout)  self)
private

Destroy interpolator for xMHD solution fields.

Parameters
[in,out]selfInterpolation object

◆ xmhd_interp_setup()

subroutine xmhd_interp_setup ( class(xmhd_interp), intent(inout)  self,
class(oft_mesh), intent(inout), target  mesh 
)
private

Setup interpolator for xMHD solution fields.

  • Fetches local vector values for interpolation
    Parameters
    [in,out]selfInterpolation object

◆ xmhd_interp_unpack()

subroutine xmhd_interp_unpack ( real(r8), dimension(40), intent(in)  vals,
type(xmhd_loc_values), intent(out)  vloc 
)
private

Unpack interpolation fields into xMHD field structure.

Parameters
[in]valsPacked fields [40]
[out]vlocUnpacked fields

◆ xmhd_lin_run()

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.

Note
This method assumes that [B0, V0, Ne0, Ti0, Te0] constitute an equilibrium.
Parameters
[in,out]equil_fieldsEquilibrium fields
[in,out]pert_fieldsPerurbed fields
[in]escaleDesired total energy to be used for rescaling (optional)

◆ xmhd_massmatrix_apply()

subroutine xmhd_massmatrix_apply ( class(oft_xmhd_massmatrix), intent(inout)  self,
class(oft_vector), intent(inout), target  a,
class(oft_vector), intent(inout)  b 
)
private

Compute the mass matrix for solution field.

b = M(a)

Parameters
[in,out]selfMass matrix object
[in,out]aSource field
[in,out]bResult of metric function

◆ xmhd_mfnk_update()

subroutine xmhd_mfnk_update ( class(oft_vector), intent(inout), target  uin)
private

Update Jacobian matrices on all levels with new fields.

Parameters
[in,out]uinCurrent field

◆ xmhd_plot()

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
Parameters
[in,out]probesProbe object (optional)

◆ xmhd_probe_flush()

subroutine xmhd_probe_flush ( class(oft_xmhd_probe), intent(inout)  self)
private

Flush internal write buffers for probe object.

Parameters
[in,out]selfProbe object

◆ xmhd_profile()

subroutine xmhd_profile ( class(oft_vector), intent(inout)  u)
private

Simple subroutine to compute timing of different solve phases.

◆ xmhd_read_settings()

subroutine xmhd_read_settings ( real(r8), intent(out)  dt,
real(r8), intent(out)  lin_tol,
real(r8), intent(out)  nl_tol,
integer(i4), intent(out)  rst_ind,
integer(i4), intent(out)  nsteps,
integer(i4), intent(out)  rst_freq,
integer(i4), intent(out)  nclean,
integer(i4), intent(out)  maxextrap,
integer(i4), intent(out)  ittarget,
integer(i4), intent(out)  nl_update 
)
private

Read settings for extended MHD model from input files.

Runtime options are set in the main input file using the xmhd_options group.

Warning
Most default physical values in the namelist below will not be appropriate for new simulations. They are merely placeholders to indicate the type of the input.

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
Parameters
[out]dtMaximum timestep
[out]lin_tolLinear solver tolerance
[out]nl_tolNonlinear solver tolerance
[out]rst_indIndex of file for restart (rst_ind>0)
[out]nstepsNumber of timesteps
[out]rst_freqFrequency to save restart files
[out]ncleanFrequency to clean divergence
[out]maxextrapExtrapolation order for initial guess
[out]ittargetMaximum number of linear iterations
[out]nl_updateMaximum number of linear iterations

◆ xmhd_run()

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.

Parameters
[in,out]initial_fieldsInitial conditions
[in,out]driverForcing object
[in,out]probesProbe object
[in]profile_onlyProfile operator timing and stop?

◆ xmhd_set_bc()

subroutine xmhd_set_bc
private

Set BC flags.

◆ xmhd_set_level()

subroutine xmhd_set_level ( integer(i4), intent(in)  level)
private

Set the current level for xMHD model.

Parameters
[in]levelDesired level

◆ xmhd_set_ops()

subroutine xmhd_set_ops ( class(oft_vector), intent(inout), target  uin)
private

Update Jacobian matrices on all levels with new solution.

Parameters
[in,out]uinCurrent solution

◆ xmhd_setup_regions()

subroutine xmhd_setup_regions
private

Setup material regions from XML input file.

◆ xmhd_setup_rep()

subroutine xmhd_setup_rep
private

Setup composite FE representation and ML environment.

◆ xmhd_upwind_weight()

pure real(r8) function xmhd_upwind_weight ( real(r8), intent(in)  v_mag,
real(r8), intent(in)  he,
real(r8), intent(in)  nu 
)
private

Evalute upwinding parameter.

Upwinding parameter for inconsistent streamline upwinded Petrov-Galerkin method (distorted bases applied to advective terms only).

Returns
Upwinding weight
Parameters
[in]v_magAdvection velocity [m/s]
[in]heElement size [m]
[in]nuDissipation amplitude [m^2/s^2]

Variable Documentation

◆ bbc

character(len=2) bbc = 'bc'
private

Magnetic field BC ('bc','ic')

◆ d2_dens

real(r8) d2_dens = -1.d0
private

Particle hyper-diffusivity.

◆ d_dens

real(r8) d_dens = 1.E0_r8
private

Particle diffusivity.

◆ den_floor

real(r8), public den_floor = -1.d0

Density floor.

◆ den_scale

real(r8), public den_scale = -1.d0

Density scale factor.

◆ dt_scale

real(r8) dt_scale = 1._r8
private

Time step scaling factor.

◆ eta

real(r8) eta = 1.E-2_r8
private

Resistivity \( \eta / \mu_0 \).

◆ eta_hyper

real(r8) eta_hyper = -1.d0
private

Hyper resistivity.

◆ eta_reg

real(r8), dimension(:), allocatable eta_reg
private

Resistivity scale factors.

◆ eta_temp

real(r8) eta_temp = -1._r8
private

Reference temperature for resistivity.

◆ g_accel

real(r8), public g_accel = 0.d0

Gravity (oriented in the -Z direction)

◆ j2_ind

integer(i4) j2_ind = -1
private

Index of hyper-res aux field.

◆ j2_scale

real(r8), public j2_scale = -1.d0

Hyper-resistivity aux variable scale factor.

◆ jac_dt

real(r8) jac_dt = 1.E-3_r8
private

Time step.

◆ k_boltz

real(r8) k_boltz = elec_charge
private

Boltzmann constant.

◆ kappa_par

real(r8) kappa_par = 1.d0
private

Parallel thermal conductivity factor.

◆ kappa_perp

real(r8) kappa_perp = 1.d0
private

Perpendicular thermal conductivity factor.

◆ m_ion

real(r8) m_ion = -1.d0
private

Ion mass in SI units.

◆ me_factor

real(r8) me_factor = 100.E0_r8
private

Artificial electron mass factor.

◆ mesh

class(oft_mesh), pointer mesh
private

◆ mfmat

type(oft_mf_matrix), target mfmat
private

Matrix free operator.

◆ mg_mesh

class(multigrid_mesh), pointer mg_mesh
private

◆ ml_int

type(oft_matrix_ptr), dimension(:), pointer ml_int => NULL()
private

MG interpolation operators.

◆ ml_j

type(oft_matrix_ptr), dimension(:), pointer ml_j => NULL()
private

MG Jacobian operators.

◆ ml_xmhd_rep

type(oft_ml_fem_comp_type), target ml_xmhd_rep
private

ML container for field representation.

◆ mu_ion

real(r8) mu_ion = 2.d0
private

Ion mass in atomic units.

◆ n2_ind

integer(i4) n2_ind = -1
private

Index of hyper-diff aux field.

◆ n2_scale

real(r8), public n2_scale = -1.d0

Hyper-diffusivity aux variable scale factor.

◆ nbc

character(len=1) nbc = 'd'
private

Density BC ('d','n')

◆ neg_flag

real(r8), dimension(:,:), allocatable neg_flag
private

◆ neg_source

real(r8), dimension(:,:), allocatable neg_source
private

◆ nu_par

real(r8) nu_par = 0.d0
private

Fluid (parallel) viscosity (note: dynamic viscosity is nu=nu_par*den_scale/n0)

◆ nu_perp

real(r8) nu_perp = 0.d0
private

Fluid perpendicular viscosity if anisotropic, as nu_par.

◆ nu_xmhd

integer(i4), dimension(fem_max_levels) nu_xmhd = 1
private

Number of smoother iterations.

◆ oft_hcurl

class(oft_hcurl_fem), pointer oft_hcurl => NULL()
private

◆ oft_hgrad

class(oft_h1_fem), pointer oft_hgrad => NULL()
private

◆ oft_lagrange

class(oft_scalar_fem), pointer oft_lagrange => NULL()
private

◆ oft_xmhd_ops

type(xmhd_ops), pointer oft_xmhd_ops => NULL()
private

Operator container.

◆ oft_xmhd_ops_ml

type(xmhd_ops), dimension(:), pointer oft_xmhd_ops_ml => NULL()
private

MG operator container.

◆ res_profile

procedure(oft_1d_func), pointer, public res_profile => NULL()

Resistivity profile.

◆ solid_cell

logical, dimension(:), allocatable solid_cell
private

Solid region flag.

◆ tbc

character(len=1) tbc = 'd'
private

Temperature BC ('d','n')

◆ te_factor

real(r8) te_factor = 1.d0
private

Electron temperature factor (single temperature)

◆ te_ind

integer(i4) te_ind = -1
private

Index of electron temperature field.

◆ temp_floor

real(r8), public temp_floor = -1.d0

Temperature floor.

◆ temp_gamma

real(r8) temp_gamma = 5.d0/3.d0
private

Ratio of specific heats.

◆ vbc

character(len=4) vbc = 'all'
private

Velocity BC ('all','norm','tang')

◆ vbc_type

integer(i4) vbc_type = 1
private

Velocity BC type.

◆ vel_scale

real(r8), public vel_scale = 1.d3

Velocity scale factor.

◆ visc_itype

integer(i4) visc_itype = 1
private

Viscosity type flag (1->kin,2->iso,3->ani)

◆ visc_type

character(len=3) visc_type = 'kin'
private

Viscosity type ('kin','iso','ani')

◆ xmhd_adv_b

logical, public xmhd_adv_b = .TRUE.

Advance magnetic field.

◆ xmhd_adv_den

logical xmhd_adv_den = .TRUE.
private

Advance density.

◆ xmhd_adv_temp

logical xmhd_adv_temp = .TRUE.
private

Advance temperature.

◆ xmhd_advec

logical xmhd_advec = .TRUE.
private

Include fluid advection.

◆ xmhd_blevel

integer(i4) xmhd_blevel = 0
private

Highest level on base meshes.

◆ xmhd_bnorm_force

logical, public xmhd_bnorm_force = .TRUE.

Force B-norm to zero in cleaning.

◆ xmhd_brag

logical xmhd_brag = .FALSE.
private

Braginskii thermal conduction.

◆ xmhd_diss_centered

logical xmhd_diss_centered = .FALSE.
private

Time-center dissipation terms.

◆ xmhd_eps

real(r8) xmhd_eps = 1.d-10
private

Epsilon for magnetic field normalization.

◆ xmhd_hall

logical xmhd_hall = .FALSE.
private

Include Hall terms.

◆ xmhd_hcurl_cop

real(r8), dimension(:,:), pointer, contiguous xmhd_hcurl_cop => NULL()
private

◆ xmhd_hcurl_rop

real(r8), dimension(:,:), pointer, contiguous xmhd_hcurl_rop => NULL()
private

◆ xmhd_hgrad_rop

real(r8), dimension(:,:), pointer, contiguous xmhd_hgrad_rop => NULL()
private

◆ xmhd_jcb

logical xmhd_jcb = .TRUE.
private

Include JxB force on fluid.

◆ xmhd_lag_gop

real(r8), dimension(:,:), pointer, contiguous xmhd_lag_gop => NULL()
private

◆ xmhd_lag_rop

real(r8), dimension(:), pointer, contiguous xmhd_lag_rop => NULL()
private

◆ xmhd_lev

integer(i4) xmhd_lev = 1
private

Active FE level.

◆ xmhd_level

integer(i4) xmhd_level = 1
private

Active FE level.

◆ xmhd_lin_level

integer(i4) xmhd_lin_level = 1
private

Highest linear element level.

◆ xmhd_linear

logical xmhd_linear = .FALSE.
private

Compute linearized matrix.

◆ xmhd_mat_mask

integer(i4), dimension(:,:), allocatable xmhd_mat_mask
private

Matrix block mask.

◆ xmhd_mfnk

logical xmhd_mfnk = .FALSE.
private

Use Jacobian free non-linear advance.

◆ xmhd_minlev

integer(i4), public xmhd_minlev = -1

Lowest MG level.

◆ xmhd_ml_h1

type(oft_ml_fem_type), pointer, public xmhd_ml_h1 => NULL()

◆ xmhd_ml_h1grad

type(oft_ml_fem_type), pointer, public xmhd_ml_h1grad => NULL()

◆ xmhd_ml_hcurl

type(oft_ml_fem_type), pointer, public xmhd_ml_hcurl => NULL()

◆ xmhd_ml_hcurl_grad

type(oft_ml_fem_comp_type), pointer, public xmhd_ml_hcurl_grad => NULL()

◆ xmhd_ml_lagrange

type(oft_ml_fem_type), pointer, public xmhd_ml_lagrange => NULL()

◆ xmhd_ml_vecspace

type(ml_xmhd_vecspace), target xmhd_ml_vecspace
private

◆ xmhd_ml_vlagrange

type(oft_ml_fem_comp_type), pointer, public xmhd_ml_vlagrange => NULL()

◆ xmhd_monitor_div

logical xmhd_monitor_div = .FALSE.
private

Monitor div error?

◆ xmhd_nlevels

integer(i4) xmhd_nlevels = 1
private

Number of total levels.

◆ xmhd_nparts

integer(i4), dimension(fem_max_levels) xmhd_nparts = 1
private

Number of local partitions for preconditioning.

◆ xmhd_ohmic

logical xmhd_ohmic = .FALSE.
private

Resistive heating.

◆ xmhd_opcount

integer(i4) xmhd_opcount = 1
private

Number of time steps since preconditioner update.

◆ xmhd_opdt

real(r8) xmhd_opdt = 1.d0
private

Time step for current preconditioner object.

◆ xmhd_pre

class(oft_solver), pointer xmhd_pre => NULL()
private

Preconditioner object.

◆ xmhd_pre_node

type(xml_node), pointer xmhd_pre_node => NULL()
private

preconditioner XML node

◆ xmhd_prefreq

integer(i4) xmhd_prefreq = 30
private

Desired update frequency for preconditioner.

◆ xmhd_rep

type(oft_fem_comp_type), pointer xmhd_rep => NULL()
private

Active field representation.

◆ xmhd_root_node

type(xml_node), pointer xmhd_root_node => NULL()
private

xMHD XML node

◆ xmhd_rst_version

integer(i4), parameter xmhd_rst_version = 3

Restart file version number.

◆ xmhd_rw

logical xmhd_rw = .FALSE.
private

Resistive wall present.

◆ xmhd_skip_update

logical xmhd_skip_update = .FALSE.
private

Internal flag for skipping Jac update.

◆ xmhd_taxis

integer(i4), public xmhd_taxis = 3

Axis for toroidal flux and current.

◆ xmhd_therm_equil

logical xmhd_therm_equil = .FALSE.
private

Include cross-species thermalization.

◆ xmhd_two_temp

logical xmhd_two_temp = .FALSE.
private

Run 2 Temperature model.

◆ xmhd_upwind

logical xmhd_upwind = .TRUE.
private

Use upwinding.

◆ xmhd_vbcdir

logical xmhd_vbcdir = .FALSE.
private

Velocity BC type flag.

◆ xmhd_visc_heat

logical xmhd_visc_heat = .FALSE.
private

Viscous heating.