Module Solvers

The Solvers module implements the algorithms which evaluate the semi-discrete residual corresponding to the discretization of an AbstractConservationLaw in space using the operators and geometric information contained within the SpatialDiscretization to create an ODEProblem, which can be solved using your choice of time-marching method from OrdinaryDiffEq.jl.

Overview

StableSpectralElements.jl is based on a semi-discrete or "method-of-lines" approach to the numerical solution of partial differential equations, in which the spatial discretization is performed first in order to obtain a system of ordinary differential equations of the form $\underline{u}'(t) = \underline{R}(\underline{u}(t),t)$, where $\underline{u}(t) \in \mathbb{R}^{N_p \cdot N_c \cdot N_e}$ is the global solution array containing the $N_p$ coefficients for each of the $N_c$ solution variables and each of the $N_e$ mesh elements. These systems can be solved using standard time-marching methods from OrdinaryDiffEq.jl by computing the semi-discrete residual $\underline{R}(\underline{u}(t),t)$ using the in-place function semi_discrete_residual!.

Reference

StableSpectralElements.Solvers.semi_discrete_residual!Function
semi_discrete_residual!(dudt::AbstractArray{Float64,3},
                       u::AbstractArray{Float64, 3},
                       solver::Solver,
                       t::Float64)

In-place function for evaluating the right-hand side of $\underline{u}'(t) = \underline{R}(\underline{u}(t),t)$ within the ordinary differential equation solver.

Arguments

  • dudt::AbstractArray{Float64,3}: Time derivative (first index is nodal or modal coefficient, second index is solution variable, third index is mesh element)
  • u::AbstractArray{Float64,3}: Current global solution state (first index is nodal or modal coefficient, second index is solution variable, third index is mesh element)
  • solver::Solver: Parameter of type Solver containing all the information defining the spatial discretization, the algorithms used to evaluate the semi-discrete residual, and preallocated arrays used for temporary storage; type parameters of Solver are used to dispatch the resulting algorithm based on the equation type, operator type, mass matrix solver (i.e. CholeskySolver, DiagonalSolver, or WeightAdjustedSolver), residual form (i.e. StandardForm or FluxDifferencingForm), and parallelism (i.e. Serial or Threaded).
  • t::Float64: Time variable
source
StableSpectralElements.Solvers.FluxDifferencingFormType
FluxDifferencingForm(; mapping_form = SkewSymmetricMapping(),
               inviscid_numerical_flux = LaxFriedrichsNumericalFlux(),
               viscous_numerical_flux = BR1(),
               two_point_flux = EntropyConservativeFlux())

Type used for dispatch indicating the use of a flux-differencing (e.g. entropy-stable or kinetic-energy-preserving) discontinuous spectral-element formulation based on two-point flux functions.

source
StableSpectralElements.Solvers.FluxDifferencingOperatorsType
FluxDifferencingOperators{S_type, C_type, V_type, Vt_type, R_type, Rt_type} <:
   AbstractDiscretizationOperators

Set of operators and nodal values of geometric factors used to evaluate the semi-discrete residual for the FluxDifferencingForm in reference space. Note that no physical-operator formulation is implemented for the flux-differencing form.

source
StableSpectralElements.Solvers.SolverType
Solver{ConservationLaw <: AbstractConservationLaw, 
       Operators <: AbstractDiscretizationOperators,
       MassSolver <: AbstractMassMatrixSolver,
       ResidualForm <: AbstractMassMatrixSolver,
       Parallelism <: AbstractParallelism,
       PreAllocatedArrays <: AbstractPreAllocatedArrays}

Composite type defining the spatial discretization which is passed into semi_discrete_residual! function to compute the time derivative within an ordinary differential equation solver from OrdinaryDiffEq.jl. The algorithm used to compute the semi-discrete residual is dispatched based on the type parameters.

source
StableSpectralElements.Solvers.StandardFormType
StandardForm(; mapping_form = SkewSymmetricMapping(),
               inviscid_numerical_flux = LaxFriedrichsNumericalFlux(),
               viscous_numerical_flux = BR1())

Type used for dispatch indicating the use of a standard (i.e. not flux-differencing) weak-form discontinuous spectral-element method. The mapping_form argument can be set to StandardMapping() to recover a standard discontinuous Galerkin formulation.

source