1 #ifndef PHILIP_DG_BASE_STATE_HPP 2 #define PHILIP_DG_BASE_STATE_HPP 4 #include <deal.II/distributed/tria.h> 7 #include "parameters/all_parameters.h" 13 #if PHILIP_DIM==1 // dealii::parallel::distributed::Triangulation<dim> does not work for 1D 14 template <
int dim,
int nstate,
typename real,
typename MeshType = dealii::Triangulation<dim>>
16 template <
int dim,
int nstate,
typename real,
typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
30 const unsigned int degree,
31 const unsigned int max_degree_input,
32 const unsigned int grid_degree_input,
33 const std::shared_ptr<Triangulation> triangulation_input);
49 std::shared_ptr < Physics::ModelBase<dim, nstate, FadType > >
pde_model_fad;
51 std::unique_ptr < NumericalFlux::NumericalFluxConvective<dim, nstate, FadType > >
conv_num_flux_fad;
53 std::unique_ptr < NumericalFlux::NumericalFluxDissipative<dim, nstate, FadType > >
diss_num_flux_fad;
58 std::shared_ptr < Physics::ModelBase<dim, nstate, RadType > >
pde_model_rad;
60 std::unique_ptr < NumericalFlux::NumericalFluxConvective<dim, nstate, RadType > >
conv_num_flux_rad;
62 std::unique_ptr < NumericalFlux::NumericalFluxDissipative<dim, nstate, RadType > >
diss_num_flux_rad;
110 real
evaluate_CFL (std::vector< std::array<real,nstate> > soln_at_q,
const real artificial_dissipation,
const real cell_diameter,
const unsigned int cell_degree);
120 #endif // PHILIP_DG_BASE_STATE_HPP DGBaseState(const Parameters::AllParameters *const parameters_input, const unsigned int degree, const unsigned int max_degree_input, const unsigned int grid_degree_input, const std::shared_ptr< Triangulation > triangulation_input)
< Input parameters.
std::unique_ptr< NumericalFlux::NumericalFluxDissipative< dim, nstate, real > > diss_num_flux_double
Dissipative numerical flux with real type.
void set_use_auxiliary_eq()
Set use_auxiliary_eq flag.
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.
void update_model_variables()
Update the necessary variables declared in src/physics/model.h.
std::unique_ptr< NumericalFlux::NumericalFluxConvective< dim, nstate, FadFadType > > conv_num_flux_fad_fad
Convective numerical flux with FadFadType.
std::unique_ptr< NumericalFlux::NumericalFluxConvective< dim, nstate, real > > conv_num_flux_double
Convective numerical flux with real type.
std::unique_ptr< NumericalFlux::NumericalFluxConvective< dim, nstate, FadType > > conv_num_flux_fad
Convective numerical flux with FadType.
Files for the baseline physics.
std::shared_ptr< Physics::ModelBase< dim, nstate, FadFadType > > pde_model_fad_fad
Contains the model terms of the PDEType == PhysicsModel with FadFadType.
std::shared_ptr< Physics::ModelBase< dim, nstate, real > > pde_model_double
Contains the model terms of the PDEType == PhysicsModel with real type.
std::unique_ptr< NumericalFlux::NumericalFluxConvective< dim, nstate, RadFadType > > conv_num_flux_rad_fad
Convective numerical flux with RadFadDtype.
Main parameter class that contains the various other sub-parameter classes.
void reset_numerical_fluxes()
Reinitializes the numerical fluxes based on the current physics.
std::unique_ptr< NumericalFlux::NumericalFluxDissipative< dim, nstate, RadFadType > > diss_num_flux_rad_fad
Dissipative numerical flux with RadFadDtype.
std::shared_ptr< Physics::PhysicsBase< dim, nstate, FadFadType > > pde_physics_fad_fad
Contains the physics of the PDE with FadFadType.
std::shared_ptr< Physics::ModelBase< dim, nstate, RadFadType > > pde_model_rad_fad
Contains the model terms of the PDEType == PhysicsModel with RadFadType.
std::shared_ptr< Physics::PhysicsBase< dim, nstate, real > > pde_physics_double
Contains the physics of the PDE with real type.
void allocate_model_variables()
Allocate the necessary variables declared in src/physics/model.h.
std::shared_ptr< Physics::PhysicsBase< dim, nstate, RadType > > pde_physics_rad
Contains the physics of the PDE with RadType.
std::shared_ptr< ArtificialDissipationBase< dim, nstate > > artificial_dissip
Link to Artificial dissipation class (with three dissipation types, depending on the input)...
real evaluate_CFL(std::vector< std::array< real, nstate > > soln_at_q, const real artificial_dissipation, const real cell_diameter, const unsigned int cell_degree)
Evaluate the time it takes for the maximum wavespeed to cross the cell domain.
std::shared_ptr< Physics::PhysicsBase< dim, nstate, RadFadType > > pde_physics_rad_fad
Contains the physics of the PDE with RadFadDtype.
std::unique_ptr< NumericalFlux::NumericalFluxDissipative< dim, nstate, FadType > > diss_num_flux_fad
Dissipative numerical flux with FadType.
std::shared_ptr< Physics::ModelBase< dim, nstate, FadType > > pde_model_fad
Contains the model terms of the PDEType == PhysicsModel with FadType.
std::unique_ptr< NumericalFlux::NumericalFluxDissipative< dim, nstate, RadType > > diss_num_flux_rad
Dissipative numerical flux with RadType.
std::shared_ptr< Physics::ModelBase< dim, nstate, RadType > > pde_model_rad
Contains the model terms of the PDEType == PhysicsModel with RadType.
std::unique_ptr< NumericalFlux::NumericalFluxDissipative< dim, nstate, FadFadType > > diss_num_flux_fad_fad
Dissipative numerical flux with FadFadType.
std::shared_ptr< Physics::PhysicsBase< dim, nstate, FadType > > pde_physics_fad
Contains the physics of the PDE with FadType.
Abstract class templated on the number of state variables.
std::unique_ptr< NumericalFlux::NumericalFluxConvective< dim, nstate, RadType > > conv_num_flux_rad
Convective numerical flux with RadType.
DGBase is independent of the number of state variables.
void set_physics(std::shared_ptr< Physics::PhysicsBase< dim, nstate, real > > pde_physics_double_input, std::shared_ptr< Physics::PhysicsBase< dim, nstate, FadType > > pde_physics_fad_input, std::shared_ptr< Physics::PhysicsBase< dim, nstate, RadType > > pde_physics_rad_input, std::shared_ptr< Physics::PhysicsBase< dim, nstate, FadFadType > > pde_physics_fad_fad_input, std::shared_ptr< Physics::PhysicsBase< dim, nstate, RadFadType > > pde_physics_rad_fad_input)