Face Element Integrals
This page describes the functions that evaluate the face element integrals for a single interface. The functions that loop over all the interfaces are located on the face integrals page. These integrals require data from all the nodes of the elements rather than the face nodes as with regular face integrals.
These integrals are used by the entropy stable scheme, and some of them internally use a numerical flux function. This flux function must satisfy an entropy property for the resulting scheme to be entropy stable! The IR flux function is typically used.
Functions
EulerEquationMod.calcECFaceIntegral
— Method.Calculate the face integrals in an entropy conservative manner for a given interface. Unlike standard face integrals, this requires data from the entirety of both elements, not just data interpolated to the face
resL and resR are updated with the results of the computation for the left and right elements, respectively.
Note that nrm_xy must contains the normal vector in x-y space at the face nodes.
The flux function must be symmetric!
Aliasing restrictions: none, although its unclear what the meaning of this function would be if resL and resR alias
Performance note: the version in the tests is the same speed as this one for p=1 Omega elements and about 10% faster for p=4 elements, but would not be able to take advantage of the sparsity of R for SBP Gamma elements
EulerEquationMod.calcECFaceIntegral
— Method.Method for sparse faces. See other method for details
Aliasing restrictions: params.flux_vals1 must not be in use
EulerEquationMod.calcESLFFaceIntegral
— Method.Calculate the face integral in an entropy stable manner using Lax-Friedrich type dissipation. This uses calcECFaceIntegral and calcLFEntropyPenaltyIntegral internally, see those functions for details.
EulerEquationMod.calcESLW2FaceIntegral
— Method.Calculate the face integral in an entropy stable manner using Lax-Wendroff type dissipation. This uses calcECFaceIntegral and calcLW2EntropyPenaltyIntegral internally, see those functions for details.
EulerEquationMod.calcESLWFaceIntegral
— Method.Calculate the face integral in an entropy stable manner using approximate Lax-Wendroff type dissipation. This uses calcECFaceIntegral and calcLWEntropyPenaltyIntegral internally, see those functions for details.
EulerEquationMod.calcEntropyFix
— Method.This function modifies the eigenvalues of the euler flux jacobian such that if any value is zero, a little dissipation is still added. The absolute values of the eigenvalues modified eigenvalues are calculated.
Methods are available for 2 and 3 dimensions
This function depends on the ordering of the eigenvalues produced by calcEvals.
Inputs: params: ParamType, used to dispatch to 2 or 3D method
Inputs/Outputs: Lambda: vector of eigenvalues to be modified
Aliasing restrictions: none
Calculate a term that provably dissipates (mathematical) entropy using a Lax-Friedrich type of dissipation. This requires data from the left and right element volume nodes, rather than face nodes for a regular face integral.
Note that nrm_face must contain the scaled face normal vector in x-y space at the face nodes, and qL, qR, resL, and resR are the arrays for the entire element, not just the face.
Aliasing restrictions: params.nrm2, params.A0, w_vals_stencil, w_vals2_stencil
Method for sparse faces. See other method for details
Aliasing restrictions: params: v_vals, v_vals2, q_vals, A0, res_vals1, res_vals2
Calculate a term that provably dissipates (mathematical) entropy using a Lax-Wendroff type of dissipation. This requires data from the left and right element volume nodes, rather than face nodes for a regular face integral.
Note nrm_face must contain the scaled normal vector in x-y space at the face nodes, and qL, qR, resL, and resR are the arrays for the entire element, not just the face.
Implementation Detail: Because the scaling does not exist in arbitrary directions for 3D, the function projects q into n-t coordinates, computes the eigendecomposition there, and then rotates back
Aliasing restrictions: from params the following fields are used: Y, S2, Lambda, res_vals1, res_vals2, w_vals_stencil, w_vals2_stencil, v_vals, v_vals2, q_vals, q_vals2, nrm2, P
Method for sparse faces. See other method for details
Aliasing restrictions: params: v_vals, v_vals2, q_vals, q_vals2, A0, res_vals1, res_vals2 A0, S2, Lambda, nrm2, P
Calculate a term that provably dissipates (mathematical) entropy using a an approximation to Lax-Wendroff type of dissipation. This requires data from the left and right element volume nodes, rather than face nodes for a regular face integral.
Note that nrm_face must contain the scaled normal vector in x-y space at the face nodes, and qL, qR, resL, and resR are the arrays for the entire element, not just the face.
The approximation to Lax-Wendroff is the computation of
for i=1:Tdim abs(niY_iS2_iLambda_iY_i.') end
rather than computing the flux jacobian in the normal direction.
Aliasing restrictions: from params the following fields are used: Y, S2, Lambda, res_vals1, res_vals2, res_vals3, w_vals_stencil, w_vals2_stencil, v_vals, v_vals2, q_vals, q_vals2
EulerEquationMod.getFaceElementFunctors
— Method.Populates the field(s) of the EulerData object with FaceElementIntegralType
functors as specified by the options dictionary
Inputs
mesh: an AbstractMesh
sbp: an SBP operator
opts: the options dictionary
Inputs/Outputs
eqn: the EulerData object
EulerEquationMod.getEntropyLFStab
— Method.This function computes the entropy dissipation term using Lax-Friedrich type dissipation. The term is evaluated using simple averaging of qL and qR. The term is subtracted off of F.
This function is dimension agnostic, but only works for conservative variables.
Aliasing restrictions: params.q_vals3, see also getEntropyLFStab_inner
EulerEquationMod.getEntropyLFStab_inner
— Method.Updates the vector F with the stabilization term from Carpenter, Fisher, Nielsen, Frankel, Entrpoy stable spectral collocation schemes for the Navier-Stokes equatiosn: Discontinuous interfaces. The term is subtracted off from F.
The q_avg vector should some average of qL and qR, but the type of averaging is left up to the user.
This function is agnostic to dimension, but only works for conservative variables.
Aliasing: from params the following arrays are used: A0, v_vals v_vals2.
EulerEquationMod.getIRA0
— Method.Computes dq/dv, where q are the conservative variables and v are the IR entropy variables. This is equiavlent to calcA0 scaled by gamma_1, but computed from the conservative variables, which is much less expensive.
Methods are available for 2 and 3 dimensions. A0 is overwritten with the result
EulerEquationMod.getLambdaMaxSimple
— Method.Calculates the maximum magnitude eigenvalue of the Euler flux jacobian at the arithmatic average of two states.
This functions works in both 2D and 3D Inputs: params: ParamType, conservative variable qL: left state qR: right state dir: direction vector (does not have to be unit vector)
Outputs: lambda_max: eigenvalue of maximum magnitude
Aliasing restrictions: params.q_vals3 must be unused
Flux Functors
EulerEquationMod.ECFaceIntegral
— Type.Entropy conservative term only
Lax-Friedrich entropy penalty term only
Lax-Wendroff entropy penalty term only
Approximate Lax-Wendroff entropy penalty term only
Entropy conservative integral + Lax-Friedrich penalty
Entropy conservative integral + Lax-Wendroff penalty
Entropy conservative integral + approximate Lax-Wendroff penalty