1 #ifndef __VISCOUS_NUMERICAL_FLUX__ 2 #define __VISCOUS_NUMERICAL_FLUX__ 4 #include <deal.II/base/tensor.h> 5 #include <deal.II/base/types.h> 6 #include "physics/physics.h" 7 #include "dg/artificial_dissipation.h" 10 namespace NumericalFlux {
13 template<
int dim,
int nstate,
typename real>
27 const std::array<real, nstate> &soln_int,
28 const std::array<real, nstate> &soln_ext,
29 const dealii::Tensor<1,dim,real> &normal_int)
const = 0;
33 const dealii::types::global_dof_index current_cell_index,
34 const dealii::types::global_dof_index neighbor_cell_index,
35 const real artificial_diss_coeff_int,
36 const real artificial_diss_coeff_ext,
37 const std::array<real, nstate> &soln_int,
38 const std::array<real, nstate> &soln_ext,
39 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_int,
40 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_ext,
41 const dealii::Tensor<1,dim,real> &normal_int,
43 const bool on_boundary =
false)
const = 0;
46 const std::shared_ptr < Physics::PhysicsBase<dim, nstate, real> >
pde_physics;
51 template<
int dim,
int nstate,
typename real>
66 const std::array<real, nstate> &soln_int,
67 const std::array<real, nstate> &soln_ext,
68 const dealii::Tensor<1,dim,real> &normal_int)
const override;
78 const dealii::types::global_dof_index current_cell_index,
79 const dealii::types::global_dof_index neighbor_cell_index,
80 const real artificial_diss_coeff_int,
81 const real artificial_diss_coeff_ext,
82 const std::array<real, nstate> &soln_int,
83 const std::array<real, nstate> &soln_ext,
84 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_int,
85 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_ext,
86 const dealii::Tensor<1,dim,real> &normal_int,
88 const bool on_boundary =
false)
const override;
93 template<
int dim,
int nstate,
typename real>
108 const std::array<real, nstate> &soln_int,
109 const std::array<real, nstate> &soln_ext,
110 const dealii::Tensor<1,dim,real> &normal_int)
const override;
120 const dealii::types::global_dof_index current_cell_index,
121 const dealii::types::global_dof_index neighbor_cell_index,
122 const real artificial_diss_coeff_int,
123 const real artificial_diss_coeff_ext,
124 const std::array<real, nstate> &soln_int,
125 const std::array<real, nstate> &soln_ext,
126 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_int,
127 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_ext,
128 const dealii::Tensor<1,dim,real> &normal_int,
130 const bool on_boundary =
false)
const override;
134 template<
int dim,
int nstate,
typename real>
149 const std::array<real, nstate> &soln_int,
150 const std::array<real, nstate> &soln_ext,
151 const dealii::Tensor<1,dim,real> &normal_int)
const override;
161 const dealii::types::global_dof_index current_cell_index,
162 const dealii::types::global_dof_index neighbor_cell_index,
163 const real artificial_diss_coeff_int,
164 const real artificial_diss_coeff_ext,
165 const std::array<real, nstate> &soln_int,
166 const std::array<real, nstate> &soln_ext,
167 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_int,
168 const std::array<dealii::Tensor<1,dim,real>, nstate> &soln_grad_ext,
169 const dealii::Tensor<1,dim,real> &normal_int,
171 const bool on_boundary =
false)
const override;
virtual 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 =0
Solution flux at the interface.
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.
virtual ~NumericalFluxDissipative()=default
Abstract class must have a virtual destructor and an implementation.
Base class of numerical flux associated with dissipation.
Files for the baseline physics.
const std::shared_ptr< ArtificialDissipationBase< dim, nstate > > artificial_dissip
Link to artificial dissipation.
SymmetricInternalPenalty(std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input, std::shared_ptr< ArtificialDissipationBase< dim, nstate >> artificial_dissipation_input)
Constructor.
BassiRebay2(std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input, std::shared_ptr< ArtificialDissipationBase< dim, nstate >> artificial_dissipation_input)
Constructor.
virtual 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 =0
Auxiliary flux at the interface.
CentralViscousNumericalFlux(std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input, std::shared_ptr< ArtificialDissipationBase< dim, nstate >> artificial_dissipation_input)
Constructor.
NumericalFluxDissipative(std::shared_ptr< Physics::PhysicsBase< dim, nstate, real >> physics_input, std::shared_ptr< ArtificialDissipationBase< dim, nstate >> artificial_dissipation_input)
Constructor.
Class to add artificial dissipation with an option to add one of the 3 dissipation types: 1...
const std::shared_ptr< Physics::PhysicsBase< dim, nstate, real > > pde_physics
Associated physics.
Symmetric interior penalty method.