[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
physics_model.h
1 #ifndef __PHYSICS_MODEL__
2 #define __PHYSICS_MODEL__
3 
5 #include "physics.h"
6 #include "navier_stokes.h"
7 #include "model.h"
8 
9 namespace PHiLiP {
10 namespace Physics {
11 
13 template <int dim, int nstate, typename real, int nstate_baseline_physics>
14 class PhysicsModel : public PhysicsBase <dim, nstate, real>
15 {
16 public:
19  const Parameters::AllParameters *const parameters_input,
21  std::shared_ptr< ModelBase<dim,nstate,real> > model_input,
23  const bool has_nonzero_diffusion,
24  const bool has_nonzero_physical_source);
25 
27  const int n_model_equations;
28 
30  std::shared_ptr< PhysicsBase<dim,nstate_baseline_physics,real> > physics_baseline;
31 
33  std::shared_ptr< ModelBase<dim,nstate,real> > model;
34 
36  std::array<dealii::Tensor<1,dim,real>,nstate> convective_flux (
37  const std::array<real,nstate> &conservative_soln) const;
38 
40  std::array<dealii::Tensor<1,dim,real>,nstate> dissipative_flux (
41  const std::array<real,nstate> &conservative_soln,
42  const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
43  const dealii::types::global_dof_index cell_index) const;
44 
46  std::array<real,nstate> physical_source_term (
47  const dealii::Point<dim,real> &pos,
48  const std::array<real,nstate> &conservative_soln,
49  const std::array<dealii::Tensor<1,dim,real>,nstate> &solution_gradient,
50  const dealii::types::global_dof_index cell_index) const;
51 
53  std::array<real,nstate> source_term (
54  const dealii::Point<dim,real> &pos,
55  const std::array<real,nstate> &conservative_soln,
56  const real current_time,
57  const dealii::types::global_dof_index cell_index) const;
58 
59  //===========================================================================================
60  // All other functions required by PhysicsBase:
61  //===========================================================================================
63  std::array<dealii::Tensor<1,dim,real>,nstate> convective_numerical_split_flux (
64  const std::array<real,nstate> &conservative_soln1,
65  const std::array<real,nstate> &conservative_soln2) const;
66 
68  std::array<real,nstate> compute_entropy_variables (
69  const std::array<real,nstate> &conservative_soln) const;
70 
73  const std::array<real,nstate> &entropy_var) const;
74 
78  std::array<real,nstate> convective_eigenvalues (
79  const std::array<real,nstate> &/*solution*/,
80  const dealii::Tensor<1,dim,real> &/*normal*/) const;
81 
83  real max_convective_eigenvalue (const std::array<real,nstate> &soln) const;
84 
87  const std::array<real,nstate> &soln,
88  const dealii::Tensor<1,dim,real> &normal) const override;
89 
91  real max_viscous_eigenvalue (const std::array<real,nstate> &soln) const;
92 
95  const int /*boundary_type*/,
96  const dealii::Point<dim, real> &/*pos*/,
97  const dealii::Tensor<1,dim,real> &/*normal*/,
98  const std::array<real,nstate> &/*soln_int*/,
99  const std::array<dealii::Tensor<1,dim,real>,nstate> &/*soln_grad_int*/,
100  std::array<real,nstate> &/*soln_bc*/,
101  std::array<dealii::Tensor<1,dim,real>,nstate> &/*soln_grad_bc*/) const;
102 
104 
106  dealii::Vector<double> post_compute_derived_quantities_vector (
107  const dealii::Vector<double> &uh,
108  const std::vector<dealii::Tensor<1,dim> > &duh,
109  const std::vector<dealii::Tensor<2,dim> > &dduh,
110  const dealii::Tensor<1,dim> &normals,
111  const dealii::Point<dim> &evaluation_points) const;
112 
114 
116  std::vector<std::string> post_get_names () const;
117 
119 
121  std::vector<dealii::DataComponentInterpretation::DataComponentInterpretation> post_get_data_component_interpretation () const;
122 
124 
126  dealii::UpdateFlags post_get_needed_update_flags () const;
127 
128 protected:
129  const MPI_Comm mpi_communicator;
130  const int mpi_rank;
131  const int n_mpi;
132 
135  dealii::ConditionalOStream pcout;
136 };
137 
138 } // Physics namespace
139 } // PHiLiP namespace
140 
141 #endif
const int n_model_equations
Number of model equations (i.e. those additional to the baseline physics)
Definition: physics_model.h:27
dealii::Vector< double > post_compute_derived_quantities_vector(const dealii::Vector< double > &uh, const std::vector< dealii::Tensor< 1, dim > > &duh, const std::vector< dealii::Tensor< 2, dim > > &dduh, const dealii::Tensor< 1, dim > &normals, const dealii::Point< dim > &evaluation_points) const
Returns current vector solution to be used by PhysicsPostprocessor to output current solution...
std::array< dealii::Tensor< 1, dim, real >, nstate > dissipative_flux(const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const
Dissipative (i.e. viscous) flux: .
const bool has_nonzero_diffusion
Flag to signal that diffusion term is non-zero.
Definition: physics.h:59
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.
Definition: physics.h:34
Manufactured solution used for grid studies to check convergence orders.
const int mpi_rank
MPI rank.
Physics Model equations. Derived from PhysicsBase, holds a baseline physics and model terms and equat...
Definition: physics_model.h:14
void boundary_face_values(const int, const dealii::Point< dim, real > &, const dealii::Tensor< 1, dim, real > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &, std::array< real, nstate > &, std::array< dealii::Tensor< 1, dim, real >, nstate > &) const
Evaluates boundary values and gradients on the other side of the face.
PartialDifferentialEquation
Possible Partial Differential Equations to solve.
Files for the baseline physics.
Definition: ADTypes.hpp:10
std::array< real, nstate > physical_source_term(const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const std::array< dealii::Tensor< 1, dim, real >, nstate > &solution_gradient, const dealii::types::global_dof_index cell_index) const
Physical source term.
std::vector< dealii::DataComponentInterpretation::DataComponentInterpretation > post_get_data_component_interpretation() const
Returns DataComponentInterpretation of the solution to be used by PhysicsPostprocessor to output curr...
Physics model additional terms and equations to the baseline physics.
Definition: model.h:18
Main parameter class that contains the various other sub-parameter classes.
std::array< dealii::Tensor< 1, dim, real >, nstate > convective_flux(const std::array< real, nstate > &conservative_soln) const
Convective flux: .
real max_viscous_eigenvalue(const std::array< real, nstate > &soln) const
Maximum viscous eigenvalue.
std::vector< std::string > post_get_names() const
Returns names of the solution to be used by PhysicsPostprocessor to output current solution...
std::array< real, nstate > source_term(const dealii::Point< dim, real > &pos, const std::array< real, nstate > &conservative_soln, const real current_time, const dealii::types::global_dof_index cell_index) const
Source term that does not require differentiation.
std::array< real, nstate > convective_eigenvalues(const std::array< real, nstate > &, const dealii::Tensor< 1, dim, real > &) const
real max_convective_eigenvalue(const std::array< real, nstate > &soln) const
Maximum convective eigenvalue.
std::shared_ptr< PhysicsBase< dim, nstate_baseline_physics, real > > physics_baseline
Baseline physics object with nstate==nstate_baseline_physics.
Definition: physics_model.h:30
std::array< real, nstate > compute_entropy_variables(const std::array< real, nstate > &conservative_soln) const
Computes the entropy variables.
std::shared_ptr< ModelBase< dim, nstate, real > > model
Model object.
Definition: physics_model.h:33
dealii::ConditionalOStream pcout
ConditionalOStream.
const bool has_nonzero_physical_source
Flag to signal that physical source term is non-zero.
Definition: physics.h:62
dealii::UpdateFlags post_get_needed_update_flags() const
Returns required update flags of the solution to be used by PhysicsPostprocessor to output current so...
const MPI_Comm mpi_communicator
MPI communicator.
std::array< dealii::Tensor< 1, dim, real >, nstate > convective_numerical_split_flux(const std::array< real, nstate > &conservative_soln1, const std::array< real, nstate > &conservative_soln2) const
Convective Numerical Split Flux for split form.
real max_convective_normal_eigenvalue(const std::array< real, nstate > &soln, const dealii::Tensor< 1, dim, real > &normal) const override
Maximum convective normal eigenvalue (used in Lax-Friedrichs)
std::array< real, nstate > compute_conservative_variables_from_entropy_variables(const std::array< real, nstate > &entropy_var) const
Computes the conservative variables from the entropy variables.
std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function
Manufactured solution function.
Definition: physics.h:71
PhysicsModel(const Parameters::AllParameters *const parameters_input, Parameters::AllParameters::PartialDifferentialEquation baseline_physics_type, std::shared_ptr< ModelBase< dim, nstate, real > > model_input, std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function, const bool has_nonzero_diffusion, const bool has_nonzero_physical_source)
Constructor.