1 #ifndef __POSITIVITY_PRESERVING_LIMITER__ 2 #define __POSITIVITY_PRESERVING_LIMITER__ 4 #include "bound_preserving_limiter.h" 19 template<
int dim,
int nstate,
typename real>
31 std::shared_ptr<BoundPreservingLimiterState<dim, nstate, real>>
tvbLimiter;
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);
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,
63 const std::array<std::vector<real>,
nstate>& soln_at_q,
64 const unsigned int n_quad_pts,
70 const double density_avg,
71 const double density_min,
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);
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 > ¤t_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.
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.