[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
positivity_preserving_limiter.h
1 #ifndef __POSITIVITY_PRESERVING_LIMITER__
2 #define __POSITIVITY_PRESERVING_LIMITER__
3 
4 #include "bound_preserving_limiter.h"
5 
6 namespace PHiLiP {
8 /**********************************
9 * Zhang, Xiangxiong, and Chi-Wang Shu.
10 * "On positivity-preserving high order discontinuous Galerkin schemes
11 * for compressible Euler equations on rectangular meshes."
12 * Journal of Computational Physics 229.23 (2010): 8918-8934.
13 **********************************/
14 /**********************************
15 * Wang, Cheng, et al.
16 * "Robust high order discontinuous Galerkin schemes for two-dimensional gaseous detonations."
17 * Journal of Computational Physics 231.2 (2012): 653-665.
18 **********************************/
19 template<int dim, int nstate, typename real>
20 class PositivityPreservingLimiter : public BoundPreservingLimiterState <dim, nstate, real>
21 {
22 public:
25  const Parameters::AllParameters* const parameters_input);
26 
28  ~PositivityPreservingLimiter() = default;
29 
31  std::shared_ptr<BoundPreservingLimiterState<dim, nstate, real>> tvbLimiter;
32 
34  std::shared_ptr < Physics::Euler<dim, nstate, double > > euler_physics;
35 
38 
42  void limit(
43  dealii::LinearAlgebra::distributed::Vector<double>& solution,
44  const dealii::DoFHandler<dim>& dof_handler,
45  const dealii::hp::FECollection<dim>& fe_collection,
46  const dealii::hp::QCollection<dim>& volume_quadrature_collection,
47  const unsigned int grid_degree,
48  const unsigned int max_degree,
49  const dealii::hp::FECollection<1> oneD_fe_collection_1state,
50  const dealii::hp::QCollection<1> oneD_quadrature_collection);
51 protected:
53  std::vector<real> get_theta2_Zhang2010(
54  const std::vector< real >& p_lim,
55  const std::array<real, nstate>& soln_cell_avg,
56  const std::array<std::vector<real>, nstate>& soln_at_q,
57  const unsigned int n_quad_pts,
58  const double eps,
59  const double gamma);
60 
63  const std::array<std::vector<real>, nstate>& soln_at_q,
64  const unsigned int n_quad_pts,
65  const double p_avg);
66 
70  const double density_avg,
71  const double density_min,
72  const double pos_eps,
73  const double p_avg);
74 
78  dealii::LinearAlgebra::distributed::Vector<double>& solution,
79  const std::array<std::vector<real>, nstate>& soln_coeff,
80  const unsigned int n_shape_fns,
81  const std::vector<dealii::types::global_dof_index>& current_dofs_indices);
82 
83 }; // End of PositivityPreservingLimiter Class
84 } // PHiLiP namespace
85 
86 #endif
87 
real get_theta2_Wang2012(const std::array< std::vector< real >, nstate > &soln_at_q, const unsigned int n_quad_pts, const double p_avg)
Obtain the theta value used to scale all the states using 3.7 in Wang, Shu 2012.
void write_limited_solution(dealii::LinearAlgebra::distributed::Vector< double > &solution, const std::array< std::vector< real >, nstate > &soln_coeff, const unsigned int n_shape_fns, const std::vector< dealii::types::global_dof_index > &current_dofs_indices)
std::shared_ptr< BoundPreservingLimiterState< dim, nstate, real > > tvbLimiter
Pointer to TVB limiter class (TVB limiter can be applied in conjunction with this limiter) ...
Base Class for bound preserving limiters templated on state.
PositivityPreservingLimiter(const Parameters::AllParameters *const parameters_input)
Constructor.
const int nstate
Number of states.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Main parameter class that contains the various other sub-parameter classes.
std::vector< real > get_theta2_Zhang2010(const std::vector< real > &p_lim, const std::array< real, nstate > &soln_cell_avg, const std::array< std::vector< real >, nstate > &soln_at_q, const unsigned int n_quad_pts, const double eps, const double gamma)
Obtain the theta value used to scale all the states using 3.16-3.18 in Zhang, Shu 2010...
~PositivityPreservingLimiter()=default
Destructor.
real get_density_scaling_value(const double density_avg, const double density_min, const double pos_eps, const double p_avg)
Class for implementation of two forms of the Positivity-Preserving limiter derived from BoundPreservi...
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)
std::shared_ptr< Physics::Euler< dim, nstate, double > > euler_physics
Euler physics pointer. Used to compute pressure.