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

Detailed Description

Plasma parameters and evaluation routines.

Authors
Chris Hansen
Date
April 2013

Functions/Subroutines

pure real(r8) function, dimension(2) brag_comb_transport (n, t, mu, bmag, lambda, me_factor)
 Evalute Braginskii transport coefficients for a combined temperature eq.
 
pure real(r8) function, dimension(2) brag_elec_transport (n, t, bmag, lambda, me_factor)
 Evalute Braginskii transport coefficients for the electron fluid.
 
pure real(r8) function, dimension(2) brag_ion_transport (n, t, mu, bmag, lambda)
 Evalute Braginskii transport coefficients for the ion fluid.
 
pure real(r8) function elec_ion_therm_rate (te, n, mu)
 Evalute electron-ion thermalization rate.
 
pure real(r8) function electron_drift_speed (j, n)
 Evalute drift velocity for Electrons.
 
pure real(r8) function eta_chodura (j, t, n, m, m_factor)
 Evalute Chodura resistivity.
 
pure real(r8) function eta_spitzer (t, lam)
 Evalute Spitzer resistivity.
 
pure real(r8) function ion_visc (n, t, mu, lam)
 Evalute simple ion viscosity.
 
pure real(r8) function log_lambda (t, n)
 Evalute ln( Lambda )
 
pure real(r8) function plasma_freq (n, m)
 Evalute plasma frequency for singly charged species.
 
pure real(r8) function sound_speed (t, m)
 Evalute sound speed.
 

Variables

real(r8), parameter elec_charge = 1.602176565E-19_r8
 \( q_e \) - Electron charge [C]
 
real(r8), parameter elec_mass = 9.10938291E-31_r8
 \( m_e \) - Electron mass [Kg]
 
real(r8), parameter ep0 = 8.854187817E-12_r8
 \( \epsilon_0 \) - Permitivity of free-space
 
real(r8), parameter mu0 = pi*(4.E-7_r8)
 \( \mu_0 \) - Permeability of free-space
 
real(r8), parameter proton_mass = 1.672621777E-27_r8
 \( m_p \) - Proton mass [Kg]
 

Function/Subroutine Documentation

◆ brag_comb_transport()

pure real(r8) function, dimension(2) brag_comb_transport ( real(r8), intent(in)  n,
real(r8), intent(in)  t,
real(r8), intent(in)  mu,
real(r8), intent(in)  bmag,
real(r8), intent(in), optional  lambda,
real(r8), intent(in), optional  me_factor 
)

Evalute Braginskii transport coefficients for a combined temperature eq.

Note
The highest physical conduction coefficient is chosen for each direction (usually parallel->electron, perpendicular->ion)
Parameters
[in]nDensity [m^-3]
[in]tTemperature [eV]
[in]muReduced ion mass
[in]bmagMagnetic field strength [T]
[in]lambdaCoulomb logarithm {11} (optional)
[in]me_factorElectron mass enhancement {1} (optional)
Returns
\( \left( \chi_{\parallel,e}, min(\chi_{\perp,i},\chi_{\parallel,i}) \right) \)

◆ brag_elec_transport()

pure real(r8) function, dimension(2) brag_elec_transport ( real(r8), intent(in)  n,
real(r8), intent(in)  t,
real(r8), intent(in)  bmag,
real(r8), intent(in), optional  lambda,
real(r8), intent(in), optional  me_factor 
)

Evalute Braginskii transport coefficients for the electron fluid.

Parameters
[in]nDensity [m^-3]
[in]tTemperature [eV]
[in]bmagMagnetic field strength [T]
[in]lambdaCoulomb logarithm (optional: default = 11)
[in]me_factorElectron mass enhancement {1} (optional)
Returns
\( \left( \chi_{\parallel,e}, min(\chi_{\perp,e},\chi_{\parallel,e}) \right) \)

◆ brag_ion_transport()

pure real(r8) function, dimension(2) brag_ion_transport ( real(r8), intent(in)  n,
real(r8), intent(in)  t,
real(r8), intent(in)  mu,
real(r8), intent(in)  bmag,
real(r8), intent(in), optional  lambda 
)

Evalute Braginskii transport coefficients for the ion fluid.

Note
The highest physical conduction coefficient is chosen for each direction (usually parallel->electron, perpendicular->ion)
Parameters
[in]nDensity [m^-3]
[in]tTemperature [eV]
[in]muReduced ion mass
[in]bmagMagnetic field strength [T]
[in]lambdaCoulomb logarithm (optional: default = 11)
Returns
\( \left( \chi_{\parallel,i}, min(\chi_{\perp,i},\chi_{\parallel,i}) \right) \)

◆ elec_ion_therm_rate()

pure real(r8) function elec_ion_therm_rate ( real(r8), intent(in)  te,
real(r8), intent(in)  n,
real(r8), intent(in)  mu 
)

Evalute electron-ion thermalization rate.

Definition from page 34 of NRL plasma formulary

Parameters
[in]teElectron temperature [eV]
[in]nDensity [m^-3]
[in]muReduced ion mass
Returns
\( \tau = 1 / (3.2E-15 * Z^2 \lambda / (\mu * T^{3/2})) \)

◆ electron_drift_speed()

pure real(r8) function electron_drift_speed ( real(r8), intent(in)  j,
real(r8), intent(in)  n 
)

Evalute drift velocity for Electrons.

Parameters
[in]jCurrent density [A/m^2]
[in]nNumber density [m^-3]
Returns
\( j / ne \)

◆ eta_chodura()

pure real(r8) function eta_chodura ( real(r8), intent(in)  j,
real(r8), intent(in)  t,
real(r8), intent(in)  n,
real(r8), intent(in)  m,
real(r8), intent(in), optional  m_factor 
)

Evalute Chodura resistivity.

Parameters
[in]jCurrent density [A/m^2]
[in]tTemperature [eV]
[in]nNumber density [m^-3]
[in]mIon mass [kg]
[in]m_factorElectron mass enhancement (optional)
Returns
\( \frac{C m_e \omega_{p,i}}{n e^2} ( 1 - e^{-v_{d,e} / f v_{s,i}} ) \)

◆ eta_spitzer()

pure real(r8) function eta_spitzer ( real(r8), intent(in)  t,
real(r8), intent(in), optional  lam 
)

Evalute Spitzer resistivity.

Parameters
[in]tTemperature [eV]
[in]lam!< Coulomb logarithm (optional: default = 11)
Returns
\( 5.253e-05 \frac{\lambda}{T^{3/2} \mu_0} \)

◆ ion_visc()

pure real(r8) function ion_visc ( real(r8), intent(in)  n,
real(r8), intent(in)  t,
real(r8), intent(in), optional  mu,
real(r8), intent(in), optional  lam 
)

Evalute simple ion viscosity.

\( \tau_i = 2.09 \times 10^7 T_i^{3/2} \mu^{1/2} / (n/10^6) \lambda \)

Parameters
[in]nDensity [m^-3]
[in]tTemperature [eV]
[in]muReduced ion mass (optional: default = 1.0)
[in]lamCoulomb logarithm (optional: default = 11)
Returns
\( 0.96 n k T_i \tau_i \)

◆ log_lambda()

pure real(r8) function log_lambda ( real(r8), intent(in)  t,
real(r8), intent(in)  n 
)

Evalute ln( Lambda )

Definition for electron-ion collisions in regime \( T_i m_e / m_i < 10 Z^2 < T_e \), see page 34 of NRL plasma formulary

Parameters
[in]tTemperature [eV]
[in]nDensity [m^-3]
Returns
\( \lambda = 24 - ln( \sqrt{n/10^6} / T ) \)

◆ plasma_freq()

pure real(r8) function plasma_freq ( real(r8), intent(in)  n,
real(r8), intent(in)  m 
)

Evalute plasma frequency for singly charged species.

Parameters
[in]nNumber density [m^-3]
[in]mParticle mass [kg]
Returns
\( \left( 4 \pi n e^2 / m \right)^{1/2} \)

◆ sound_speed()

pure real(r8) function sound_speed ( real(r8), intent(in)  t,
real(r8), intent(in)  m 
)

Evalute sound speed.

Parameters
[in]tFluid temperature [eV]
[in]mParticle mass [kg]
Returns
\( \left( 3 k T / m \right)^{1/2} \)

Variable Documentation

◆ elec_charge

real(r8), parameter elec_charge = 1.602176565E-19_r8

\( q_e \) - Electron charge [C]

◆ elec_mass

real(r8), parameter elec_mass = 9.10938291E-31_r8

\( m_e \) - Electron mass [Kg]

◆ ep0

real(r8), parameter ep0 = 8.854187817E-12_r8

\( \epsilon_0 \) - Permitivity of free-space

◆ mu0

real(r8), parameter mu0 = pi*(4.E-7_r8)

\( \mu_0 \) - Permeability of free-space

◆ proton_mass

real(r8), parameter proton_mass = 1.672621777E-27_r8

\( m_p \) - Proton mass [Kg]