WARPXM v1.10.0
Loading...
Searching...
No Matches
imhd.h File Reference
#include <cmath>
#include "warpxm/warpxm_config.h"
#include "apps/general_functions/geometry.h"

Go to the source code of this file.

Functions

real max_wave_speed (const real gas_gamma, const real *q, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real max_wave_speed_magnitude (const real gas_gamma, const real *q, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real min_wave_speed (real gas_gamma, const real *q, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real calc_fast_magnetosonic (real gas_gamma, const real *q, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real calc_slow_magnetosonic (real gas_gamma, const real *q, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void calc_magnetosonic_speeds (const real gas_gamma, const real *q, real *c, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 Calculate magnetosonic speeds in local frame, normal to element boundary.
 
real rusanov_flux (const real gas_gamma, const real *ql, const real *qr, const real *F_L, const real *F_R, real *F_num, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real hll_flux (const real gas_gamma, const real *ql, const real *qr, const real *F_L, const real *F_R, real *F_num, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real hlld_flux (const real gas_gamma, const real *ql, const real *qr, const real *F_L, const real *F_R, real *F_num, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
real get_avg (real a, real b)
 
real get_delta_lambda_k (const real lambda_k_R, const real lambda_k_L)
 
real get_lambda_k_star (const real lambda_k_roe, const real lambda_k_R, const real lambda_k_L)
 
real roe_flux (const real gas_gamma, const real *ql, const real *qr, const real *F_L, const real *F_R, real *F_num, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void imhd_numerical_flux_local_frame (const real gas_gamma, const real *q, real *f, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 IMHD numerical flux vector calculated in frame rotated to element boundary.
 
real imhd_numerical_flux (const real gas_gamma, const real *q_l, const real *q_r, const real *aux_l, const real *aux_r, const solverVariables_t *pFV, const std::string &flux_type, real *numericalFlux, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 IMHD numerical flux vector calculated in global frame.
 
void imhd_internal_flux (const real gas_gamma, const solverVariables_t *pSV, const real *q, std::vector< std::vector< real > > &internalFlux, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 IMHD internal flux tensor calculated in global frame.
 
void CtoNC_MHD (const real gas_gamma, const real QC[8], real QNC[8], const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 Converts conserved variables to nonconserved (primitive) variables.
 
void NCtoC_MHD (const real gas_gamma, const real QNC[8], real QC[8], const real rho_min=std::numeric_limits< real >::epsilon())
 Converts nonconserved (primitive) variables to conserved variables.
 
void dqdw_MHD (const real gas_gamma, const real Q[8], real dqdw[][8], const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void dwdq_MHD (const real gas_gamma, const real Q[8], real dwdq[][8], const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void conducting_wall_noslip (const real gas_gamma, const solverVariables_t *pFV, const real *QC, const real vx_w, const real vy_w, const real vz_w, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void conducting_wall_noslip_old (const real gas_gamma, const solverVariables_t *pFV, const real *QC, const real vx_w, const real vy_w, const real vz_w, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void conducting_wall_noslip_gradients (const real gas_gamma, const solverVariables_t *pFV, const real *QC, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void conducting_wall_freeslip (const real gas_gamma, const solverVariables_t *pFV, const real *QC, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void conducting_wall_freeslip_gradients (const real gas_gamma, const solverVariables_t *pFV, const real *QC, const real *aux_in, const real *aux_out, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void twofluid_to_imhd (const real *QC, const solverVariables_t *pSV, real oc_sd, real *QC_w)
 
void copyout (const real *QC, real *QC_w)
 
void copyout_gradients (const real *grad_QC, real *grad_QC_w)
 
void freeslipWallCondition_imhd (const real gas_gamma, const solverVariables_t *pFV, const real *QC, real *QC_w, const real rho_min=std::numeric_limits< real >::epsilon(), const real press_min=std::numeric_limits< real >::epsilon())
 
void applyLdgPenaltyEta (const solverVariables_t *pFV, const real *q_l, const real *q_r, real *numericalFlux)
 Apply LDG eta penalty to numerical flux.
 

Function Documentation

◆ applyLdgPenaltyEta()

void applyLdgPenaltyEta ( const solverVariables_t pFV,
const real q_l,
const real q_r,
real numericalFlux 
)

Apply LDG eta penalty to numerical flux.

See Datta dissertation, Section 3.2 for construction of the LDG method as implemented here. The "eta" penalty-weighting parameter in the code is refered to as "alpha" in the dissertation, and the stabalizing penalty "tau" in Hesthaven "Nodal Discontinuous Galerkin Methods" 2008.

See Hesthaven, Chapter 7.2 for further reading on the topic, and the selection of values for eta/alpha/tau.

Parameters
constsolverVariables_t pFV - solver variable structure
constreal q_l - input variable, left
constreal q_r - input variable, right
realnumericalFlux - numerical fluxes

◆ calc_fast_magnetosonic()

real calc_fast_magnetosonic ( real  gas_gamma,
const real q,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ calc_magnetosonic_speeds()

void calc_magnetosonic_speeds ( const real  gas_gamma,
const real q,
real c,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

Calculate magnetosonic speeds in local frame, normal to element boundary.

@warn Intended for use with numerical flux calcs. Looks for speeds normal to element boundary in local frame (x direction).

Parameters
gas_gammaAdiabatic index.
qArray of MHD variables.
cArray of speeds.
rho_minMinimum density.
press_minMinimum pressure.

◆ calc_slow_magnetosonic()

real calc_slow_magnetosonic ( real  gas_gamma,
const real q,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ conducting_wall_freeslip()

void conducting_wall_freeslip ( const real  gas_gamma,
const solverVariables_t pFV,
const real QC,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ conducting_wall_freeslip_gradients()

void conducting_wall_freeslip_gradients ( const real  gas_gamma,
const solverVariables_t pFV,
const real QC,
const real aux_in,
const real aux_out,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ conducting_wall_noslip()

void conducting_wall_noslip ( const real  gas_gamma,
const solverVariables_t pFV,
const real QC,
const real  vx_w,
const real  vy_w,
const real  vz_w,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ conducting_wall_noslip_gradients()

void conducting_wall_noslip_gradients ( const real  gas_gamma,
const solverVariables_t pFV,
const real QC,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ conducting_wall_noslip_old()

void conducting_wall_noslip_old ( const real  gas_gamma,
const solverVariables_t pFV,
const real QC,
const real  vx_w,
const real  vy_w,
const real  vz_w,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ copyout()

void copyout ( const real QC,
real QC_w 
)

◆ copyout_gradients()

void copyout_gradients ( const real grad_QC,
real grad_QC_w 
)

◆ CtoNC_MHD()

void CtoNC_MHD ( const real  gas_gamma,
const real  QC[8],
real  QNC[8],
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

Converts conserved variables to nonconserved (primitive) variables.

Parameters
QCResultant conserved variables
QNCNonconserved variables
rho_minMinimum density.
press_minMinimum pressure.

◆ dqdw_MHD()

void dqdw_MHD ( const real  gas_gamma,
const real  Q[8],
real  dqdw[][8],
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ dwdq_MHD()

void dwdq_MHD ( const real  gas_gamma,
const real  Q[8],
real  dwdq[][8],
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ freeslipWallCondition_imhd()

void freeslipWallCondition_imhd ( const real  gas_gamma,
const solverVariables_t pFV,
const real QC,
real QC_w,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ get_avg()

real get_avg ( real  a,
real  b 
)

◆ get_delta_lambda_k()

real get_delta_lambda_k ( const real  lambda_k_R,
const real  lambda_k_L 
)

◆ get_lambda_k_star()

real get_lambda_k_star ( const real  lambda_k_roe,
const real  lambda_k_R,
const real  lambda_k_L 
)

◆ hll_flux()

real hll_flux ( const real  gas_gamma,
const real ql,
const real qr,
const real F_L,
const real F_R,
real F_num,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ hlld_flux()

real hlld_flux ( const real  gas_gamma,
const real ql,
const real qr,
const real F_L,
const real F_R,
real F_num,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ imhd_internal_flux()

void imhd_internal_flux ( const real  gas_gamma,
const solverVariables_t pSV,
const real q,
std::vector< std::vector< real > > &  internalFlux,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

IMHD internal flux tensor calculated in global frame.

See also
MHD writeup document
Parameters
pCIMHD constant struct.
qVector of MHD variables.
pFVSolver variable struct. Holds element and node geometry info(?).
internalFluxRank two tensor of IMHD internal fluxes being returned.
rho_minMinimum density.
press_minMinimum pressure.

◆ imhd_numerical_flux()

real imhd_numerical_flux ( const real  gas_gamma,
const real q_l,
const real q_r,
const real aux_l,
const real aux_r,
const solverVariables_t pFV,
const std::string &  flux_type,
real numericalFlux,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

IMHD numerical flux vector calculated in global frame.

See also
MHD writeup document

Procedure for calculating the numerical flux

  1. Rotate vector into element face frame.
  2. Calculate flux in face frame for left and right variables.
  3. Calculate numerical flux.
  4. Rotate back to global frame.
  5. Assign numerical fluxes to numericalFlux.
Parameters
pCIMHD constant struct.
q_lVector of MHD variables on left of element boundary.
q_rVector of MHD variable on right of element boundary.
aux_lNOT USED.
aux_rNOT USED.
pFVSolver variable struct. Holds element and node geometry info(?).
flux_typeString with name of numerical flux type.
numericalFluxVector of IMHD numerical fluxes being returned.
rho_minMinimum density.
press_minMinimum pressure.

◆ imhd_numerical_flux_local_frame()

void imhd_numerical_flux_local_frame ( const real  gas_gamma,
const real q,
real f,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

IMHD numerical flux vector calculated in frame rotated to element boundary.

In the rotated frame, x is normal to the boundary and y and z are tangential.

Ideal MHD rotated flux vector:

\begin{align} \overline{\mathcal{F}}_x=\begin{bmatrix} \rho u \\ \rho u u - (B_xB_x - \frac{1}{2}B^2) + P \\ \rho u v - B_xB_y \\ \rho u w - B_xB_z \\ u(\varepsilon+P + \frac{1}{2}B^2) - (\boldsymbol{u} \cdot \boldsymbol{B}) B_x \\ 0 \\ u B_y - B_x v \\ u B_z - B_x w \end{bmatrix} \end{align}

See also
MHD writeup document
Parameters
gas_gammaRatio of specific heats.
qIMHD fluid variables rotated into the element surface frame.
fVector of IMHD numerical fluxes.
rho_minMinimum density.
press_minMinimum pressure.

◆ max_wave_speed()

real max_wave_speed ( const real  gas_gamma,
const real q,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ max_wave_speed_magnitude()

real max_wave_speed_magnitude ( const real  gas_gamma,
const real q,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ min_wave_speed()

real min_wave_speed ( real  gas_gamma,
const real q,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ NCtoC_MHD()

void NCtoC_MHD ( const real  gas_gamma,
const real  QNC[8],
real  QC[8],
const real  rho_min = std::numeric_limits< real >::epsilon() 
)

Converts nonconserved (primitive) variables to conserved variables.

Parameters
QNCNonconserved variables
QCResultant conserved variables
rho_minMinimum density.

◆ roe_flux()

real roe_flux ( const real  gas_gamma,
const real ql,
const real qr,
const real F_L,
const real F_R,
real F_num,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ rusanov_flux()

real rusanov_flux ( const real  gas_gamma,
const real ql,
const real qr,
const real F_L,
const real F_R,
real F_num,
const real  rho_min = std::numeric_limits< real >::epsilon(),
const real  press_min = std::numeric_limits< real >::epsilon() 
)

◆ twofluid_to_imhd()

void twofluid_to_imhd ( const real QC,
const solverVariables_t pSV,
real  oc_sd,
real QC_w 
)