WARPXM v1.10.0
Loading...
Searching...
No Matches
wxm::apps::rmhd Namespace Reference

Namespaces

namespace  bc
 

Classes

class  CurrentDensity
 Class that calculates the current gradient J from deritives of B: That is, J = skin_depth_norm * curl(B) More...
 
class  CylSource
 Implements the source terms in the resistitive diffusion terms that are added to Ideal MHD to produce the resistive MHD equations arising from a cylindrical geometry. More...
 
class  DensityDiffusionCylSource
 Cylindrical source for density diffusion in the MHD model. More...
 
class  DensityDiffusionFlux
 Density diffusion flux for the MHD model. More...
 
class  ElectricField
 Class that calculates the electric field from Ohm's Law using ideal and resistive terms E = -u x B + (dp/L) * (nupt) * eta * J. More...
 
class  ohmic_source_t
 class that calculates a source associated with reduced Ohmic heating More...
 
class  WmApplication_Resistive_Diffusion
 

Functions

real current_density_cart (const real skin_depth_norm, const real *grad_B[3], real *J)
 
real current_density_cyl (const real skin_depth_norm, const real r, const real *B, const real *grad_B[3], const int ir, const int ith, const int iz, real *J)
 
real get_eta (const real *q, const real J[3], const bool const_resistivity, const real constant_resistivity_eta, const std::string resistivity_type, const real max_eta, const real min_eta, const real n_vac, const real eta_vac, const real gas_gamma, const real ompt, const real nupt, const real lnlam=10.0, const real Ai=1.0, const real Zi=1.0, const real Ae=1.0/1836.0, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real get_resistivity (const real *q, const real gas_gamma, const real lnlam=10, const real Ai=1.0, const real Zi=1.0, const real Ae=1.0/1836.0, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 Compute Spitzer-like resistivity.
 
real get_resistivity_chodura (const real *q, const real J[3], const real gas_gamma, const real ompt, const real nupt, const real lnlam=10, const real Ai=1.0, const real Zi=1.0, const real Ae=1.0/1836.0, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 Compute Chodura resistivity.
 
real rmhd_internal_diffusion_flux (const real nu_p_norm, const real skin_depth_norm, const solverVariables_t *pSV, const real *q, const real *J, const real eta, std::vector< std::vector< real > > &internal_diffusion_flux)
 
void noslipWallCondition_MHD (const real gas_gamma, const solverVariables_t *pSV, const real *QC, const real vx_w, const real vy_w, const real vz_w, const real bx_w, const real by_w, const real bz_w, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real set_time_step (const solverVariables_t *pSV, const real j_coefficient, const real diffusion_coefficient, const real eta, const real cfl_diff=0.5)
 Suggested time step calculation for RMHD resistive diffusion.
 

Function Documentation

◆ current_density_cart()

real wxm::apps::rmhd::current_density_cart ( const real  skin_depth_norm,
const real grad_B[3],
real J 
)

◆ current_density_cyl()

real wxm::apps::rmhd::current_density_cyl ( const real  skin_depth_norm,
const real  r,
const real B,
const real grad_B[3],
const int  ir,
const int  ith,
const int  iz,
real J 
)

◆ get_eta()

real wxm::apps::rmhd::get_eta ( const real q,
const real  J[3],
const bool  const_resistivity,
const real  constant_resistivity_eta,
const std::string  resistivity_type,
const real  max_eta,
const real  min_eta,
const real  n_vac,
const real  eta_vac,
const real  gas_gamma,
const real  ompt,
const real  nupt,
const real  lnlam = 10.0,
const real  Ai = 1.0,
const real  Zi = 1.0,
const real  Ae = 1.0/1836.0,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ get_resistivity()

real wxm::apps::rmhd::get_resistivity ( const real q,
const real  gas_gamma,
const real  lnlam = 10,
const real  Ai = 1.0,
const real  Zi = 1.0,
const real  Ae = 1.0/1836.0,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

Compute Spitzer-like resistivity.

Compute Spitzer-like resistivity [Goldston and Rutherford IoP (1995), Eq. 11.29] note that for parallel resistivity, this is ~2x too large; see [G&R Eq. 11.31 and associated discussion]

The normalized eta is (see Datta dissertation Appx. A (Eq. A.2.1) for normalization):

\[ \eta = \sqrt{2} Z_i \ln(\lambda) A_e^{0.5} T_e^{-1.5} \]

Warning
Users should be careful that their chosen normalization is consistent with the normalization in Datta dissertation Appx. A (Eq. A.2.1). For example, if the electron mass is not 1/1836, this resistivity calculation will not be consistent with that choice.
Parameters
qThe raw MHD variables.
gas_gammaRatio of specific heats.
lnlamCoulomb logarithm.
AiProton normalized ion mass.
ZiIon charge.
AeProton normalized electron mass.
rho_minNormalized minimum mass density.
press_minNormalized minimum mass density.
Returns
Spitzer resistivity eta.

◆ get_resistivity_chodura()

real wxm::apps::rmhd::get_resistivity_chodura ( const real q,
const real  J[3],
const real  gas_gamma,
const real  ompt,
const real  nupt,
const real  lnlam = 10,
const real  Ai = 1.0,
const real  Zi = 1.0,
const real  Ae = 1.0/1836.0,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

Compute Chodura resistivity.

Compute Chodura resistivity [Sgro and Nielson, Phys. Fluids 1976].

Warning
Users should be careful that their chosen normalization is consistent with the normalization in Datta dissertation Appx. A (Eq. A.2.1). For example, if the electron mass is not 1/1836, this resistivity calculation will not be consistent with that choice.
Parameters
qThe raw MHD variables.
JCurrent density.
gas_gammaRatio of specific heats.
omptNormalized proton plasma frequency.
nuptNormalized proton-proton collision frequency.
lnlamCoulomb logarithm.
AiProton normalized ion mass.
ZiIon charge.
AeProton normalized electron mass.
rho_minNormalized minimum mass density.
Returns
Pressure calculated assuming pressure and density floors.

◆ noslipWallCondition_MHD()

void wxm::apps::rmhd::noslipWallCondition_MHD ( const real  gas_gamma,
const solverVariables_t pSV,
const real QC,
const real  vx_w,
const real  vy_w,
const real  vz_w,
const real  bx_w,
const real  by_w,
const real  bz_w,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ rmhd_internal_diffusion_flux()

real wxm::apps::rmhd::rmhd_internal_diffusion_flux ( const real  nu_p_norm,
const real  skin_depth_norm,
const solverVariables_t pSV,
const real q,
const real J,
const real  eta,
std::vector< std::vector< real > > &  internal_diffusion_flux 
)

◆ set_time_step()

real wxm::apps::rmhd::set_time_step ( const solverVariables_t pSV,
const real  j_coefficient,
const real  diffusion_coefficient,
const real  eta,
const real  cfl_diff = 0.5 
)

Suggested time step calculation for RMHD resistive diffusion.

dB/dt + curl [ nuptau * (dp/L)^2 * eta * curl B ] = 0 j_coefficient = dp/L diffusion_coefficient = nuptau * (dp/L)

for diffusive time step limit:

  • dt_max = cfl_diff dx_eff^2 / D
  • let dx_eff = cfl*dx, where cfl is the warpxm-style cfl (which takes into account the spatial order), and dx is the element size
  • thus, dt_max = cfl^2*dx^2 / D
Parameters
pSVStruct of elemenent geometry related information.
j_coefficientdp/L
diffusion_coefficientnuptau * (dp/L)
etaResistivity.
cfl_diffCoefficient on the CFL for diffusion. Default 0.5.
Returns
Suggested maximum time step.