[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
viscous_numerical_flux.hpp
1 #ifndef __VISCOUS_NUMERICAL_FLUX__
2 #define __VISCOUS_NUMERICAL_FLUX__
3 
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"
8 
9 namespace PHiLiP {
10 namespace NumericalFlux {
11 
13 template<int dim, int nstate, typename real>
15 {
16 public:
18 NumericalFluxDissipative(std::shared_ptr<Physics::PhysicsBase<dim, nstate, real>> physics_input, std::shared_ptr<ArtificialDissipationBase<dim, nstate>> artificial_dissipation_input)
19 : pde_physics(physics_input), artificial_dissip(artificial_dissipation_input)
20 {};
21 
23 virtual ~NumericalFluxDissipative() = default;
24 
26 virtual std::array<real, nstate> evaluate_solution_flux (
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;
30 
32 virtual std::array<real, nstate> evaluate_auxiliary_flux (
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,
42  const real &penalty,
43  const bool on_boundary = false) const = 0;
44 
45 protected:
46 const std::shared_ptr < Physics::PhysicsBase<dim, nstate, real> > pde_physics;
47 const std::shared_ptr < ArtificialDissipationBase<dim, nstate> > artificial_dissip;
48 };
49 
51 template<int dim, int nstate, typename real>
52 class CentralViscousNumericalFlux: public NumericalFluxDissipative<dim, nstate, real>
53 {
56 public:
58 CentralViscousNumericalFlux(std::shared_ptr<Physics::PhysicsBase<dim, nstate, real>> physics_input, std::shared_ptr < ArtificialDissipationBase<dim, nstate>> artificial_dissipation_input)
59 : NumericalFluxDissipative<dim,nstate,real>(physics_input,artificial_dissipation_input)
60 {};
61 
63 
65 std::array<real, nstate> evaluate_solution_flux (
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;
69 
71 
77 std::array<real, nstate> evaluate_auxiliary_flux (
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,
87  const real &penalty,
88  const bool on_boundary = false) const override;
89 
90 };
91 
93 template<int dim, int nstate, typename real>
94 class SymmetricInternalPenalty: public NumericalFluxDissipative<dim, nstate, real>
95 {
98 public:
100 SymmetricInternalPenalty(std::shared_ptr<Physics::PhysicsBase<dim, nstate, real>> physics_input, std::shared_ptr < ArtificialDissipationBase<dim, nstate>> artificial_dissipation_input)
101 : NumericalFluxDissipative<dim,nstate,real>(physics_input,artificial_dissipation_input)
102 {};
103 
105 
107 std::array<real, nstate> evaluate_solution_flux (
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;
111 
113 
119 std::array<real, nstate> evaluate_auxiliary_flux (
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,
129  const real &penalty,
130  const bool on_boundary = false) const override;
131 
132 };
133 
134 template<int dim, int nstate, typename real>
135 class BassiRebay2: public NumericalFluxDissipative<dim, nstate, real>
136 {
139 public:
141 BassiRebay2(std::shared_ptr<Physics::PhysicsBase<dim, nstate, real>> physics_input, std::shared_ptr<ArtificialDissipationBase<dim, nstate>> artificial_dissipation_input)
142 : NumericalFluxDissipative<dim,nstate,real>(physics_input,artificial_dissipation_input)
143 {};
144 
146 
148 std::array<real, nstate> evaluate_solution_flux (
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;
152 
154 
160 std::array<real, nstate> evaluate_auxiliary_flux (
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,
170  const real &penalty,
171  const bool on_boundary = false) const override;
172 
173 };
174 
175 } // NumericalFlux namespace
176 } // PHiLiP namespace
177 
178 #endif
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.
Definition: physics.h:34
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.
Definition: ADTypes.hpp:10
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.