Density diffusion flux for the MHD model.
In Euler or MHD fluid systems, the diffusion of density is driven indirectly through pressure gradients and inter- or intra- species collisional processes (viscosity, resisitivity, etc.).
In some cases, enhanced density diffusion may be desirable for stability purposes or to approximate processes not captured by the model.
Density diffusion flux
The density diffusion flux in conservative form
\[
\frac{\partial \rho_p}{\partial t} = - \nabla \cdot (-D \nabla \rho_p )
\]
where \(D\) is the diffusion coefficient (see options for the calculation of \(D\) below), and the flux vector in cartesian cooridnates
\[
\overline{\overline{\mathcal{F}}} = \begin{bmatrix}
- D \frac{\partial \rho_p}{\partial x} & - D \frac{\partial \rho_p}{\partial y}
& - D \frac{\partial \rho_p}{\partial z} \\
0 & 0 & 0 \\
\vdots & \vdots & \vdots
\end{bmatrix}
\]
Diffusion coefficient options
The diffusion coefficient can have the following forms:
_ddiff_type == "basic"
. Constant density diffusion coefficient
\[
D_{basic} = \text{_ddiff}.
\]
_ddiff_type == "nimrod"
. A NIMROD-style density diffusion, where diffusion is a funciton of the density gradient
\[
D_{nim} = \min(\text{_ddiff} \cdot dx^2 \; |v \cdot \nabla n / n|, \; \text{_ddiff_limit})
\]
The suggested maximum time step is given by
\[
\Delta t_{\max} = CFL_{diff} \frac{\Delta x_{eff}^2}{D}
\]
where
\[
\Delta x_{eff} = CFL \Delta x.
\]
Cylindrical coordinates
In the case of density diffusion in cylindrical coordinates, cylindrical source terms are necessary. This app must be used in conjuction with density_diffusion_cyl_source.
- See also
- Warpy constructor: warpy.apps.mhd.rmhd.rmhd.density_diffusion_flux_mhd
|
| DensityDiffusionFlux () |
|
| ~DensityDiffusionFlux () override |
|
void | setup (const WxCryptSet &wxc) override |
| Setup.
|
|
const std::vector< int > & | getInputVariableIndexes (int flag) const override |
| Grab input variable indices.
|
|
const std::vector< int > & | getAuxiliaryVariableIndexes (int flag) const override |
| Grab auxiliary variable indices.
|
|
const std::vector< int > & | getOutputVariableIndexes (int flag) const override |
| Grab output variable indices.
|
|
real | numerical_flux (const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const override |
| Redefinition of the numerical flux function.
|
|
real | internal_flux (const real *q, const real *aux, const solverVariables_t *pSV, std::vector< std::vector< real > > &internalFlux) const override |
| Redefinition of the interal flux function.
|
|
real | numerical_flux_ldg (const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const |
| The numerical flux function using the LDG gradient method.
|
|
| WmApplication () |
|
virtual | ~WmApplication ()=default |
|
virtual void | setup (const WxCryptSet &wxc) |
|
virtual const std::vector< int > & | getAuxiliaryVariableIndexes (int flag=WMAPPLICATIONFLAG_NONE) const |
|
virtual const std::vector< int > & | getCrossVariableIndexes (int flag=WMAPPLICATIONFLAG_NONE) const |
|
const std::vector< std::string > & | getBoundaryNames () const |
|
bool | isOnBoundary (const std::string &boundaryName) const |
|
bool | has (int flag) const |
|
virtual real | numerical_flux (const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const |
|
virtual real | internal_flux (const real *q, const real *aux, const solverVariables_t *pSV, std::vector< std::vector< real > > &internalFlux) const |
|
virtual real | source (const real *q, const real *aux, const elementGeometry_t *pEG, real *source) const |
|
virtual void | bc_q (const real *q_in, const real *aux_in, const real *aux_out, const solverVariables_t *pFV, real *q_out) const |
| Boundary Condition Application which sets the boundary condition on ghost nodes.
|
|
virtual real | bcNumericalFlux (const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, real *numericalFlux) const |
|
virtual void | conserved_to_primitive (const real *q, const real *aux, real *w) const |
|
virtual void | primitive_to_conserved (const real *w, const real *aux, real *q) const |
|
virtual void | evaluate_function (const real *q, const real *aux, const solverVariables_t *pSV, real *result) const |
|
virtual void | bc_q_kinetic (const real *q_in, const real *aux_in, const solverVariables_t *pFV, real *q_out) const |
|
std::shared_ptr< std::string > | app_name () |
|
virtual const std::vector< int > & | getInputVariableIndexes (int flag=0) const |
|
virtual const std::vector< int > & | getOutputVariableIndexes (int flag=0) const |
|