Numerical Flux Functions
This page documents the numerical flux functions available in the code
bc_solvers.jl should be renamed to this
EulerEquationMod.RoeSolver
— Method.EulerEquationMod.RoeSolver
This calculates the Roe flux for boundary conditions at a node. The inputs must be in conservative variables.
Inputs
params : ParamType
q : conservative variables of the fluid
qg : conservative variables of the boundary
aux_vars : vector of all auxiliary variables at this node
nrm : scaled face normal vector (x-y space, outward normal of the face q lives on)
Outputs
flux : vector to populate with solution
Aliasing restrictions: none of the inputs can alias params.res_vals1, params.res_vals2, params.q_vals, params.flux_vals1, or params.sat
EulerEquationMod.RoeSolver
— Method.The main Roe solver. Populates flux
with the computed flux.
EulerEquationMod.RoeSolver_revm
— Method.###EulerEquationMod.RoeSolver_revm
Reverse mode of EulerEquationMod.RoeSolver
. This function computes the reverse mode of the Roe flux w.r.t the mesh metrics
Inputs
params
: Parameter objectq
: Conservative variable of the fluidqg
: Conservative variable of the boundary or the adjacent elementaux_vars
: Auxiliary variablesnrm
: Element face normal vector in the physical spaceflux_bar
: Flux value which needs to get differentiated
Output
nrm_bar
: derivaitve of the flux_bar w.r.t the mesh metrics
Aliasing Restrictions: Same as the forward function
EulerEquationMod.calcEulerFlux_Ducros
— Method.Calculates the numerical flux function associated with the Ducros flux splitting. Methods are available for 2D and 3D.
Inputs:
params:
qL: the left state
qR: the right state
aux_vars: the aux vars for the left state
dir: the direction vector
Inputs/Outputs:
F: vector to be populated with the flux
EulerEquationMod.calcEulerFlux_IR
— Method.This function calculates the Ismail-Roe numerical flux at a node in a specified direction
Inputs:
params: ParamType
qL: left state vector
qR: right state vector
aux_vars: auxiliary variable vector for qL
dir: a direction vector of length Tdim
Inputs/Outputs:
F: a numDofPerNode x Tdim matrix where each column will be populated with the flux in the direction specified by the corresponding column of nrm
Aliasing Restrictions: none
EulerEquationMod.calcEulerFlux_IR
— Method.This function computes the Ismail-Roe flux in Tdim directions at once for a given state. This is more efficient Than calling the single direction method Tdim times. Methods are available for 2 and 3 dimensions.
Inputs:
params: ParamType
qL: left state vector
qR: right state vector
aux_vars: auxiliary variable vector for qL
dir: a Tdim x Tdim matrix with each column containing a normal vector
Inputs/Outputs:
F: a numDofPerNode x Tdim matrix where each column will be populated with the flux in the direction specified by the corresponding column of nrm
Aliasing restrictions: none
EulerEquationMod.calcEulerFlux_IRSLF
— Method.This is the second method that takes in a normal vector directly. See the first method for a description of what this function does.
Inputs qL, qR: vectors conservative variables at left and right states aux_vars: aux_vars for qL nrm: a normal vector in xy space
Inputs/Outputs F: vector to be updated with the result
Alising restrictions: See getEntropyLFStab
EulerEquationMod.calcEulerFlux_IRSLW
— Method.This function is similar to calcEulerFlux_IRSLF, but uses Lax-Wendroff dissipation rather than Lax-Friedrich.
Aliasing restrictions: see getEntropyLWStab
EulerEquationMod.calcEulerFlux_standard
— Method.This function computes the split form flux corresponding to the standard flux in Tdim directions at once for a given state. This is more efficient Than calling the single direction method Tdim times. Methods are available for 2 and 3 dimensions.
Inputs: params: ParamType qL: left state vector qR: right state vector aux_vars: auxiliary variable vector for qL dir: a Tdim x Tdim matrix with each column containing a normal vector
Inputs/Outputs: F: a numDofPerNode x Tdim matrix where each column will be populated with the flux in the direction specified by the corresponding column of nrm
Aliasing restrictions: none
EulerEquationMod.calcLFFlux
— Method.Calculates the Lax-Friedrich flux function on the conservative variables
EulerEquationMod.calcSAT
— Method.###EulerEquationMod.calcSAT
Computes the simultaneous approximation term for use in computing the numerical flux, namely:
0.5*(|A_hat| - A)(qL - qR)
Arguments
params
: Parameter object of type ParamType`roe_vars: [u, v, H] at roe average state
dq
: difference between left and right statesnrm
: Normal to face in the physical spacesat
: Simultaneous approximation Term
EulerEquationMod.calcSAT_revm
— Method.###EulerEquationMod.calcSAT_revm
Reverse mode of calcSAT
Inputs
params
: Parameter object of type ParamTypenrm
: Normal to face in the physical spacedq
: Boundary condition penalty variablevel
: Velocities along X & Y directions in the physical spaceH
: Total enthalpysat_bar
: Inpute seed for sat flux whose derivative needs to be computed
Output
nrm_bar
: derivative ofsat_bar
w.r.t physical normal vector