[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
bound_preserving_limiter.h
1 #ifndef __BOUND_PRESERVING_LIMITER__
2 #define __BOUND_PRESERVING_LIMITER__
3 
4 #include <deal.II/dofs/dof_handler.h>
5 #include <deal.II/base/parameter_handler.h>
6 
7 #include "parameters/all_parameters.h"
8 #include "operators/operators.h"
9 #include "physics/euler.h"
10 
11 namespace PHiLiP {
12 
14 
22 template<int dim, typename real>
24 {
25 public:
27  explicit BoundPreservingLimiter(
28  const int nstate_input,//number of states input
29  const Parameters::AllParameters* const parameters_input);//pointer to parameters
30 
32  ~BoundPreservingLimiter() = default;
33 
35  const int nstate;
36 
39 
41  virtual void limit(
42  dealii::LinearAlgebra::distributed::Vector<double>& solution,
43  const dealii::DoFHandler<dim>& dof_handler,
44  const dealii::hp::FECollection<dim>& fe_collection,
45  const dealii::hp::QCollection<dim>& volume_quadrature_collection,
46  const unsigned int grid_degree,
47  const unsigned int max_degree,
48  const dealii::hp::FECollection<1> oneD_fe_collection_1state,
49  const dealii::hp::QCollection<1> oneD_quadrature_collection,
50  double dt) = 0;
51 }; // End of BoundPreservingLimiter Class
52 
54 template<int dim, int nstate, typename real>
56 {
57 public:
60 
63  const Parameters::AllParameters* const parameters_input);
64 
66  ~BoundPreservingLimiterState() = default;
67 
69  virtual void limit(
70  dealii::LinearAlgebra::distributed::Vector<double>& solution,
71  const dealii::DoFHandler<dim>& dof_handler,
72  const dealii::hp::FECollection<dim>& fe_collection,
73  const dealii::hp::QCollection<dim>& volume_quadrature_collection,
74  const unsigned int grid_degree,
75  const unsigned int max_degree,
76  const dealii::hp::FECollection<1> oneD_fe_collection_1state,
77  const dealii::hp::QCollection<1> oneD_quadrature_collection,
78  double dt) = 0;
79 
81  std::array<real, nstate> get_soln_cell_avg(
82  const std::array<std::vector<real>, nstate>& soln_at_q,
83  const unsigned int n_quad_pts,
84  const std::vector<real>& quad_weights);
85 
86 }; // End of BoundPreservingLimiterState Class
87 } // PHiLiP namespace
88 
89 #endif
90 
~BoundPreservingLimiter()=default
Destructor.
Base Class for bound preserving limiters templated on state.
const int nstate
Number of states.
Files for the baseline physics.
Definition: ADTypes.hpp:10
BoundPreservingLimiter(const int nstate_input, const Parameters::AllParameters *const parameters_input)
Constructor.
Main parameter class that contains the various other sub-parameter classes.
Base Class for implementation of bound preserving limiters.
virtual void limit(dealii::LinearAlgebra::distributed::Vector< double > &solution, const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &volume_quadrature_collection, const unsigned int grid_degree, const unsigned int max_degree, const dealii::hp::FECollection< 1 > oneD_fe_collection_1state, const dealii::hp::QCollection< 1 > oneD_quadrature_collection, double dt)=0
Function to limit the solution.
const Parameters::AllParameters *const all_parameters
Pointer to parameters object.