1 #ifndef __PHYSICS_MODEL__ 2 #define __PHYSICS_MODEL__ 6 #include "navier_stokes.h" 13 template <
int dim,
int nstate,
typename real,
int nstate_baseline_physics>
33 std::shared_ptr< ModelBase<dim,nstate,real> >
model;
37 const std::array<real,nstate> &conservative_soln)
const;
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;
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;
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;
64 const std::array<real,nstate> &conservative_soln1,
65 const std::array<real,nstate> &conservative_soln2)
const;
69 const std::array<real,nstate> &conservative_soln)
const;
73 const std::array<real,nstate> &entropy_var)
const;
79 const std::array<real,nstate> &,
80 const dealii::Tensor<1,dim,real> &)
const;
87 const std::array<real,nstate> &soln,
88 const dealii::Tensor<1,dim,real> &normal)
const override;
96 const dealii::Point<dim, real> &,
97 const dealii::Tensor<1,dim,real> &,
98 const std::array<real,nstate> &,
99 const std::array<dealii::Tensor<1,dim,real>,nstate> &,
100 std::array<real,nstate> &,
101 std::array<dealii::Tensor<1,dim,real>,nstate> &)
const;
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;
const int n_model_equations
Number of model equations (i.e. those additional to the baseline physics)
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.
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.
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...
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.
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.
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.
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.
dealii::ConditionalOStream pcout
ConditionalOStream.
const bool has_nonzero_physical_source
Flag to signal that physical source term is non-zero.
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.
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.