[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
PHiLiP::NumericalFlux::SymmetricInternalPenalty< dim, nstate, real > Class Template Reference

Symmetric interior penalty method. More...

#include <viscous_numerical_flux.hpp>

Inheritance diagram for PHiLiP::NumericalFlux::SymmetricInternalPenalty< dim, nstate, real >:
Collaboration diagram for PHiLiP::NumericalFlux::SymmetricInternalPenalty< dim, nstate, real >:

Public Member Functions

 SymmetricInternalPenalty (std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input, std::shared_ptr< ArtificialDissipationBase< dim, nstate >> artificial_dissipation_input)
 Constructor.
 
std::array< real, nstate > evaluate_solution_flux (const std::array< real, nstate > &soln_int, const std::array< real, nstate > &soln_ext, const dealii::Tensor< 1, dim, real > &normal_int) const override
 Evaluate solution flux at the interface. More...
 
std::array< real, nstate > evaluate_auxiliary_flux (const dealii::types::global_dof_index current_cell_index, const dealii::types::global_dof_index neighbor_cell_index, const real artificial_diss_coeff_int, const real artificial_diss_coeff_ext, const std::array< real, nstate > &soln_int, const std::array< real, nstate > &soln_ext, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_int, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_ext, const dealii::Tensor< 1, dim, real > &normal_int, const real &penalty, const bool on_boundary=false) const override
 Evaluate auxiliary flux at the interface. More...
 
- Public Member Functions inherited from PHiLiP::NumericalFlux::NumericalFluxDissipative< dim, nstate, real >
 NumericalFluxDissipative (std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input, std::shared_ptr< ArtificialDissipationBase< dim, nstate >> artificial_dissipation_input)
 Constructor.
 
virtual ~NumericalFluxDissipative ()=default
 Abstract class must have a virtual destructor and an implementation.
 

Additional Inherited Members

- Protected Attributes inherited from PHiLiP::NumericalFlux::NumericalFluxDissipative< dim, nstate, real >
const std::shared_ptr< Physics::PhysicsBase< dim, nstate, real > > pde_physics
 Associated physics.
 
const std::shared_ptr< ArtificialDissipationBase< dim, nstate > > artificial_dissip
 Link to artificial dissipation.
 

Detailed Description

template<int dim, int nstate, typename real>
class PHiLiP::NumericalFlux::SymmetricInternalPenalty< dim, nstate, real >

Symmetric interior penalty method.

Definition at line 94 of file viscous_numerical_flux.hpp.

Member Function Documentation

◆ evaluate_auxiliary_flux()

template<int dim, int nstate, typename real >
std::array< real, nstate > PHiLiP::NumericalFlux::SymmetricInternalPenalty< dim, nstate, real >::evaluate_auxiliary_flux ( const dealii::types::global_dof_index  current_cell_index,
const dealii::types::global_dof_index  neighbor_cell_index,
const real  artificial_diss_coeff_int,
const real  artificial_diss_coeff_ext,
const std::array< real, nstate > &  soln_int,
const std::array< real, nstate > &  soln_ext,
const std::array< dealii::Tensor< 1, dim, real >, nstate > &  soln_grad_int,
const std::array< dealii::Tensor< 1, dim, real >, nstate > &  soln_grad_ext,
const dealii::Tensor< 1, dim, real > &  normal_int,
const real &  penalty,
const bool  on_boundary = false 
) const
overridevirtual

Evaluate auxiliary flux at the interface.

\[ \hat{A} = {{ A \nabla u_h }} - \mu {{ A }} [[ u_h ]] \]

Note that \(\mu\) must be chosen to have a stable scheme.

Implements PHiLiP::NumericalFlux::NumericalFluxDissipative< dim, nstate, real >.

Definition at line 159 of file viscous_numerical_flux.cpp.

◆ evaluate_solution_flux()

template<int dim, int nstate, typename real >
std::array< real, nstate > PHiLiP::NumericalFlux::SymmetricInternalPenalty< dim, nstate, real >::evaluate_solution_flux ( const std::array< real, nstate > &  soln_int,
const std::array< real, nstate > &  soln_ext,
const dealii::Tensor< 1, dim, real > &  normal_int 
) const
overridevirtual

Evaluate solution flux at the interface.

\[\hat{u} = {u_h} \]

Implements PHiLiP::NumericalFlux::NumericalFluxDissipative< dim, nstate, real >.

Definition at line 147 of file viscous_numerical_flux.cpp.


The documentation for this class was generated from the following files: