Discontinuous Galerkin Models#

Compressible Euler Equations#

Grudge operators modeling compressible, inviscid flows (Euler)

Model definitions#

class grudge.models.euler.EulerOperator(dcoll: DiscretizationCollection, bdry_conditions=None, flux_type='lf', gamma=1.4, quadrature_tag=None)[source]#

This operator discretizes the Euler equations:

\[\partial_t \mathbf{Q} + \nabla\cdot\mathbf{F} = 0,\]

where \(\mathbf{Q}\) is the state vector containing density, momentum, and total energy, and \(\mathbf{F}\) is the vector of inviscid fluxes (see euler_volume_flux())

Predefined initial conditions#

grudge.models.euler.vortex_initial_condition(x_vec, t=0, center=5, mach_number=0.5, epsilon=1, gamma=1.4)[source]#

Initial condition adapted from Section 2 (equation 2) of:

K. Mattsson, M. Svärd, M. Carpenter, and J. Nordström (2006). High-order accurate computations for unsteady aerodynamics. DOI.

Helper routines and array containers#

class grudge.models.euler.ConservedEulerField(mass: meshmode.dof_array.DOFArray, energy: meshmode.dof_array.DOFArray, momentum: numpy.ndarray)[source]#
grudge.models.euler.conservative_to_primitive_vars(cv_state: ConservedEulerField, gamma=1.4)[source]#

Converts from conserved variables (density, momentum, total energy) into primitive variables (density, velocity, pressure).

Parameters:
  • cv_state – A ConservedEulerField containing the conserved variables.

  • gamma – The isentropic expansion factor for a single-species gas (default set to 1.4).

Returns:

A Tuple containing the primitive variables: (density, velocity, pressure).

grudge.models.euler.compute_wavespeed(cv_state: ConservedEulerField, gamma=1.4)[source]#

Computes the total translational wavespeed.

Parameters:
  • cv_state – A ConservedEulerField containing the conserved variables.

  • gamma – The isentropic expansion factor for a single-species gas (default set to 1.4).

Returns:

A DOFArray containing local wavespeeds.

grudge.models.euler.euler_volume_flux(dcoll: DiscretizationCollection, cv_state: ConservedEulerField, gamma=1.4)[source]#

Computes the (non-linear) volume flux for the Euler operator.

Parameters:
  • cv_state – A ConservedEulerField containing the conserved variables.

  • gamma – The isentropic expansion factor for a single-species gas (default set to 1.4).

Returns:

A ConservedEulerField containing the volume fluxes.

grudge.models.euler.euler_numerical_flux(dcoll: DiscretizationCollection, tpair: TracePair, gamma=1.4, lf_stabilization=False)[source]#

Computes the interface numerical flux for the Euler operator.

Parameters:
  • tpair – A grudge.trace_pair.TracePair containing the conserved variables on the interior and exterior sides of element facets.

  • gamma – The isentropic expansion factor for a single-species gas (default set to 1.4).

  • lf_stabilization – A boolean denoting whether to apply Lax-Friedrichs dissipation.

Returns:

A ConservedEulerField containing the interface fluxes.