1 #include "positivity_preserving_limiter.h" 2 #include "tvb_limiter.h" 3 #include <eigen/unsupported/Eigen/Polynomials> 4 #include <eigen/Eigen/Dense> 13 template <
int dim,
int nstate,
typename real>
17 , flow_solver_param(parameters_input->flow_solver_param)
18 , dx((flow_solver_param.grid_right_bound-flow_solver_param.grid_left_bound)/flow_solver_param.number_of_grid_elements_x)
19 , dy((flow_solver_param.grid_top_bound-flow_solver_param.grid_bottom_bound)/flow_solver_param.number_of_grid_elements_y)
20 , dz((flow_solver_param.grid_z_upper_bound-flow_solver_param.grid_z_lower_bound)/flow_solver_param.number_of_grid_elements_z)
24 PDE_enum pde_type = parameters_input->
pde_type;
26 std::shared_ptr< ManufacturedSolutionFunction<dim, real> > manufactured_solution_function
29 if (pde_type == PDE_enum::euler &&
nstate == dim + 2) {
30 euler_physics = std::make_shared < Physics::Euler<dim, nstate, real> >(
33 parameters_input->euler_param.gamma_gas,
34 parameters_input->euler_param.mach_inf,
35 parameters_input->euler_param.angle_of_attack,
36 parameters_input->euler_param.side_slip_angle,
37 manufactured_solution_function,
38 parameters_input->two_point_num_flux_type);
41 std::cout <<
"Error: Positivity-Preserving Limiter can only be applied for pde_type==euler" << std::endl;
48 tvbLimiter = std::make_shared < TVBLimiter<dim, nstate, real> >(parameters_input);
51 std::cout <<
"Error: Cannot create TVB limiter for dim > 1" << std::endl;
57 std::cout <<
"Error: number_of_grid_elements must be passed for all directions to use PPL Limiter." << std::endl;
62 std::cout <<
"Error: number_of_grid_elements must be passed for all directions to use PPL Limiter." << std::endl;
67 template <
int dim,
int nstate,
typename real>
69 const std::vector< real >& p_lim,
70 const std::array<real, nstate>& soln_cell_avg,
71 const std::array<std::vector<real>,
nstate>& soln_at_q,
72 const unsigned int n_quad_pts,
76 std::vector<real> theta2(n_quad_pts, 1);
77 Eigen::PolynomialSolver<double, 2> solver;
79 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
80 if (p_lim[iquad] >= eps) {
84 real s_coeff1 = (soln_at_q[
nstate - 1][iquad] - soln_cell_avg[
nstate - 1]) * (soln_at_q[0][iquad] - soln_cell_avg[0]) - 0.5 * pow(soln_at_q[1][iquad] - soln_cell_avg[1], 2);
86 real s_coeff2 = soln_cell_avg[0] * (soln_at_q[
nstate - 1][iquad] - soln_cell_avg[
nstate - 1]) + soln_cell_avg[
nstate - 1] * (soln_at_q[0][iquad] - soln_cell_avg[0])
87 - soln_cell_avg[1] * (soln_at_q[1][iquad] - soln_cell_avg[1]) - (eps / gamma) * (soln_at_q[0][iquad] - soln_cell_avg[0]);
89 real s_coeff3 = (soln_cell_avg[
nstate - 1] * soln_cell_avg[0]) - 0.5 * pow(soln_cell_avg[1], 2) - (eps / gamma) * soln_cell_avg[0];
92 s_coeff1 -= 0.5 * pow(soln_at_q[2][iquad] - soln_cell_avg[2], 2);
93 s_coeff2 -= soln_cell_avg[2] * (soln_at_q[2][iquad] - soln_cell_avg[2]);
94 s_coeff3 -= 0.5 * pow(soln_cell_avg[2], 2);
98 s_coeff1 -= 0.5 * pow(soln_at_q[3][iquad] - soln_cell_avg[3], 2);
99 s_coeff2 -= soln_cell_avg[3] * (soln_at_q[3][iquad] - soln_cell_avg[3]);
100 s_coeff3 -= 0.5 * pow(soln_cell_avg[3], 2);
103 Eigen::Vector3d coeff(s_coeff3, s_coeff2, s_coeff1);
104 solver.compute(coeff);
105 const Eigen::PolynomialSolver<double, 2>::RootType &r = solver.smallestRoot();
106 theta2[iquad] = r.real();
113 template <
int dim,
int nstate,
typename real>
115 const std::array<std::vector<real>,
nstate>& soln_at_q,
116 const unsigned int n_quad_pts,
119 std::vector<real> t2(n_quad_pts, 1);
121 std::array<real, nstate> soln_at_iquad;
124 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
125 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
126 soln_at_iquad[istate] = soln_at_q[istate][iquad];
130 if (nstate == dim + 2)
136 t2[iquad] = p_avg / (p_avg - p_lim);
138 if (t2[iquad] != 1) {
139 if (t2[iquad] >= 0 && t2[iquad] <= 1 && t2[iquad] < theta2) {
148 template <
int dim,
int nstate,
typename real>
150 const double density_avg,
151 const double density_min,
152 const double lower_bound,
156 real eps = std::min({ lower_bound, density_avg, p_avg });
157 if (eps < 0) eps = lower_bound;
160 if (density_avg - density_min > 1e-13)
161 theta = (density_avg - density_min == 0) ? 1.0 : std::min((density_avg - eps) / (density_avg - density_min), 1.0);
163 if (theta < 0 || theta > 1)
169 template <
int dim,
int nstate,
typename real>
171 dealii::LinearAlgebra::distributed::Vector<double>& solution,
172 const std::array<std::vector<real>,
nstate>& soln_coeff,
173 const unsigned int n_shape_fns,
174 const std::vector<dealii::types::global_dof_index>& current_dofs_indices)
177 for (
int istate = 0; istate <
nstate; istate++) {
178 for (
unsigned int ishape = 0; ishape < n_shape_fns; ++ishape) {
179 const unsigned int idof = istate * n_shape_fns + ishape;
180 solution[current_dofs_indices[idof]] = soln_coeff[istate][ishape];
183 if (istate == 0 && solution[current_dofs_indices[idof]] < 0) {
184 std::cout <<
"Error: Density is a negative value - Aborting... " << std::endl << solution[current_dofs_indices[idof]] << std::endl << std::flush;
189 if (istate == (nstate - 1) && solution[current_dofs_indices[idof]] < 0) {
190 std::cout <<
"Error: Total Energy is a negative value - Aborting... " << std::endl << solution[current_dofs_indices[idof]] << std::endl << std::flush;
196 if (isnan(solution[current_dofs_indices[idof]])) {
197 std::cout <<
"Error: Solution is NaN - Aborting... " << std::endl << solution[current_dofs_indices[idof]] << std::endl << std::flush;
204 template <
int dim,
int nstate,
typename real>
206 const std::array<std::array<std::vector<real>,
nstate>, dim>& soln_at_q,
207 const unsigned int n_quad_pts,
208 const std::vector<real>& quad_weights_GLL,
209 const std::vector<real>& quad_weights_GL,
212 std::array<real, nstate> soln_cell_avg;
217 }
else if (dim > 1) {
218 std::array<std::array<real, nstate>,dim> soln_cell_avg_dim;
220 for(
unsigned int idim = 0; idim < dim; ++idim) {
221 for(
unsigned int istate = 0; istate <
nstate; ++istate) {
222 soln_cell_avg_dim[idim][istate] = 0;
226 if constexpr (dim == 2) {
228 for(
unsigned int istate = 0; istate <
nstate; ++istate) {
229 unsigned int quad_pt = 0;
230 for(
unsigned int iquad=0; iquad<quad_weights_GLL.size(); ++iquad) {
231 for(
unsigned int jquad=0; jquad<quad_weights_GL.size(); ++jquad) {
232 soln_cell_avg_dim[0][istate] += quad_weights_GLL[iquad]*quad_weights_GL[jquad]*soln_at_q[0][istate][quad_pt];
239 for(
unsigned int istate = 0; istate <
nstate; ++istate) {
240 unsigned int quad_pt = 0;
241 for(
unsigned int iquad=0; iquad<quad_weights_GL.size(); ++iquad) {
242 for(
unsigned int jquad=0; jquad<quad_weights_GLL.size(); ++jquad) {
243 soln_cell_avg_dim[1][istate] += quad_weights_GL[iquad]*quad_weights_GLL[jquad]*soln_at_q[1][istate][quad_pt];
250 if constexpr (dim == 3) {
252 for(
unsigned int istate = 0; istate <
nstate; ++istate) {
253 unsigned int quad_pt = 0;
254 for(
unsigned int iquad=0; iquad<quad_weights_GLL.size(); ++iquad) {
255 for(
unsigned int jquad=0; jquad<quad_weights_GL.size(); ++jquad) {
256 for(
unsigned int kquad=0; kquad<quad_weights_GL.size(); ++kquad)
257 soln_cell_avg_dim[0][istate] += quad_weights_GLL[iquad]*quad_weights_GL[jquad]*quad_weights_GL[kquad]*soln_at_q[0][istate][quad_pt];
264 for(
unsigned int istate = 0; istate <
nstate; ++istate) {
265 unsigned int quad_pt = 0;
266 for(
unsigned int iquad=0; iquad<quad_weights_GL.size(); ++iquad) {
267 for(
unsigned int jquad=0; jquad<quad_weights_GLL.size(); ++jquad) {
268 for(
unsigned int kquad=0; kquad<quad_weights_GL.size(); ++kquad)
269 soln_cell_avg_dim[1][istate] += quad_weights_GL[iquad]*quad_weights_GLL[jquad]*quad_weights_GL[kquad]*soln_at_q[1][istate][quad_pt];
276 for(
unsigned int istate = 0; istate <
nstate; ++istate) {
277 unsigned int quad_pt = 0;
278 for(
unsigned int iquad=0; iquad<quad_weights_GL.size(); ++iquad) {
279 for(
unsigned int jquad=0; jquad<quad_weights_GL.size(); ++jquad) {
280 for(
unsigned int kquad=0; kquad<quad_weights_GLL.size(); ++kquad)
281 soln_cell_avg_dim[2][istate] += quad_weights_GL[iquad]*quad_weights_GL[jquad]*quad_weights_GLL[kquad]*soln_at_q[2][istate][quad_pt];
288 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
289 soln_cell_avg[istate] = 0;
293 const real lambda_1 = dt/this->
dx;
const real lambda_2 = dt/this->
dy; real lambda_3 = 0.0;
294 if constexpr(dim == 3)
295 lambda_3 = dt/this->
dz;
297 real max_local_wave_speed_1 = 0.0;
298 real max_local_wave_speed_2 = 0.0;
299 real max_local_wave_speed_3 = 0.0;
300 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
301 std::array<real,nstate> local_soln_at_q_1;
302 std::array<real,nstate> local_soln_at_q_2;
303 std::array<real,nstate> local_soln_at_q_3;
304 for(
unsigned int istate = 0; istate <
nstate; ++istate){
305 local_soln_at_q_1[istate] = soln_at_q[0][istate][iquad];
306 local_soln_at_q_2[istate] = soln_at_q[1][istate][iquad];
308 local_soln_at_q_3[istate] = soln_at_q[2][istate][iquad];
310 local_soln_at_q_3[istate] = 0.0;
313 const real local_wave_speed_1 = this->
euler_physics->max_convective_eigenvalue(local_soln_at_q_1);
314 const real local_wave_speed_2 = this->
euler_physics->max_convective_eigenvalue(local_soln_at_q_2);
316 real local_wave_speed_3 = 0.0;
318 local_wave_speed_3 = this->
euler_physics->max_convective_eigenvalue(local_soln_at_q_3);
320 if(local_wave_speed_1 > max_local_wave_speed_1) max_local_wave_speed_1 = local_wave_speed_1;
321 if(local_wave_speed_2 > max_local_wave_speed_2) max_local_wave_speed_2 = local_wave_speed_2;
322 if(dim == 3 && local_wave_speed_3 > max_local_wave_speed_3) max_local_wave_speed_3 = local_wave_speed_3;
326 real mu = max_local_wave_speed_1*lambda_1 + max_local_wave_speed_2*lambda_2 + max_local_wave_speed_3*lambda_3;
327 real avg_weight_1 = (max_local_wave_speed_1*lambda_1)/mu;
328 real avg_weight_2 = (max_local_wave_speed_2*lambda_2)/mu;
329 real avg_weight_3 = (max_local_wave_speed_3*lambda_3)/mu;
331 for (
unsigned int istate = 0; istate <
nstate; istate++) {
332 soln_cell_avg[istate] = avg_weight_1*soln_cell_avg_dim[0][istate] + avg_weight_2*soln_cell_avg_dim[1][istate];
334 soln_cell_avg[istate] += avg_weight_3*soln_cell_avg_dim[2][istate];
336 if (isnan(soln_cell_avg[istate])) {
337 std::cout <<
"Error: Solution Cell Avg is NaN - Aborting... " << std::endl << std::flush;
342 return soln_cell_avg;
345 template <
int dim,
int nstate,
typename real>
347 dealii::LinearAlgebra::distributed::Vector<double>& solution,
348 const dealii::DoFHandler<dim>& dof_handler,
349 const dealii::hp::FECollection<dim>& fe_collection,
350 const dealii::hp::QCollection<dim>& volume_quadrature_collection,
351 const unsigned int grid_degree,
352 const unsigned int max_degree,
353 const dealii::hp::FECollection<1> oneD_fe_collection_1state,
354 const dealii::hp::QCollection<1> oneD_quadrature_collection,
360 this->
tvbLimiter->limit(solution, dof_handler, fe_collection, volume_quadrature_collection, grid_degree, max_degree, oneD_fe_collection_1state, oneD_quadrature_collection, dt);
363 const unsigned int init_grid_degree = grid_degree;
366 dealii::QGauss<1> oneD_quad_GL(max_degree + 1);
367 dealii::QGaussLobatto<1> oneD_quad_GLL(max_degree + 1);
375 for (
auto soln_cell : dof_handler.active_cell_iterators()) {
376 if (!soln_cell->is_locally_owned())
continue;
378 std::vector<dealii::types::global_dof_index> current_dofs_indices;
380 const int i_fele = soln_cell->active_fe_index();
381 const dealii::FESystem<dim, dim>& current_fe_ref = fe_collection[i_fele];
382 const int poly_degree = current_fe_ref.tensor_degree();
384 const unsigned int n_dofs_curr_cell = current_fe_ref.n_dofs_per_cell();
387 current_dofs_indices.resize(n_dofs_curr_cell);
388 soln_cell->get_dof_indices(current_dofs_indices);
391 std::array<std::vector<real>,
nstate> soln_coeff;
393 const unsigned int n_shape_fns = n_dofs_curr_cell /
nstate;
394 real local_min_density = 1e6;
396 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
397 soln_coeff[istate].resize(n_shape_fns);
400 bool nan_check =
false;
402 for (
unsigned int idof = 0; idof < n_dofs_curr_cell; ++idof) {
403 const unsigned int istate = fe_collection[poly_degree].system_to_component_index(idof).first;
404 const unsigned int ishape = fe_collection[poly_degree].system_to_component_index(idof).second;
405 soln_coeff[istate][ishape] = solution[current_dofs_indices[idof]];
407 if (isnan(soln_coeff[istate][ishape])) {
412 const unsigned int n_quad_pts = n_shape_fns;
415 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
416 std::cout <<
"Error: Density passed to limiter is NaN - Aborting... " << std::endl;
418 for (
unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
419 std::cout << soln_coeff[istate][iquad] <<
" ";
422 std::cout << std::endl;
428 std::array<std::array<std::vector<real>, nstate>, dim> soln_at_q;
429 std::array<std::vector<real>, nstate> soln_at_q_dim;
431 for(
unsigned int idim = 0; idim < dim; idim++) {
432 for (
int istate = 0; istate <
nstate; istate++) {
433 soln_at_q_dim[istate].resize(n_quad_pts);
450 soln_at_q[idim] = soln_at_q_dim;
453 for (
unsigned int idim = 0; idim < dim; ++idim) {
454 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
455 if (soln_coeff[0][iquad] < local_min_density)
456 local_min_density = soln_coeff[0][iquad];
457 if (soln_at_q[idim][0][iquad] < local_min_density)
458 local_min_density = soln_at_q[idim][0][iquad];
462 std::vector< real > GLL_weights = oneD_quad_GLL.get_weights();
463 std::vector< real > GL_weights = oneD_quad_GL.get_weights();
464 std::array<real, nstate> soln_cell_avg;
466 soln_cell_avg =
get_soln_cell_avg_PPL(soln_at_q, n_quad_pts, oneD_quad_GLL.get_weights(), oneD_quad_GL.get_weights(), dt);
471 if (nstate == dim + 2) {
480 for (
unsigned int ishape = 0; ishape < n_shape_fns; ++ishape) {
481 soln_coeff[0][ishape] = theta*(soln_coeff[0][ishape] - soln_cell_avg[0]) + soln_cell_avg[0];
505 if (limiter_type == limiter_enum::positivity_preservingWang2012 && nstate == dim + 2) {
506 std::array<real, dim> theta2_quad;
507 for(
unsigned int idim = 0; idim < dim; ++idim) {
511 for(
unsigned int idim = 0; idim < dim; ++idim) {
512 if(theta2_quad[idim] < theta2)
513 theta2 = theta2_quad[idim];
517 if(theta2_soln < theta2)
518 theta2 = theta2_soln;
521 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
522 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
523 soln_coeff[istate][iquad] = theta2 * (soln_coeff[istate][iquad] - soln_cell_avg[istate])
524 + soln_cell_avg[istate];
529 if (limiter_type == limiter_enum::positivity_preservingZhang2010 && nstate == dim + 2) {
531 std::array<std::vector< real >, dim> p_lim_quad;
532 std::array<real, nstate> soln_at_iquad;
534 for(
unsigned int idim = 0; idim < dim; ++idim) {
535 p_lim_quad[idim].resize(n_quad_pts);
537 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
538 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
539 soln_at_iquad[istate] = soln_at_q[idim][istate][iquad];
541 p_lim_quad[idim][iquad] =
euler_physics->compute_pressure(soln_at_iquad);
545 std::array<std::vector< real >, dim> theta2_quad;
547 for(
unsigned int idim = 0; idim < dim; ++idim) {
548 theta2_quad[idim].resize(n_quad_pts);
553 std::vector< real > p_lim;
554 p_lim.resize(n_quad_pts);
555 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
556 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
557 soln_at_iquad[istate] = soln_coeff[istate][iquad];
559 p_lim[iquad] =
euler_physics->compute_pressure(soln_at_iquad);
564 for (
unsigned int istate = 0; istate <
nstate; ++istate) {
565 for (
unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
566 real min_theta2_quad = 1e6;
567 for(
unsigned int idim = 0; idim < dim; ++idim) {
568 if(theta2_quad[idim][iquad] < min_theta2_quad)
569 min_theta2_quad = theta2_quad[idim][iquad];
572 theta2 = std::min({ min_theta2_quad, theta2_soln[iquad] });
573 soln_coeff[istate][iquad] = theta2 * (soln_coeff[istate][iquad] - soln_cell_avg[istate])
574 + soln_cell_avg[istate];
580 std::cout <<
"Error: Theta2 is NaN - Aborting... " << std::endl << theta2 << std::endl << std::flush;
PartialDifferentialEquation pde_type
Store the PDE type to be solved.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
LimiterType
Limiter type to be applied on the solution.
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) override
LimiterParam limiter_param
Contains parameters for limiter.
real get_theta2_Wang2012(const std::array< std::vector< real >, nstate > &soln_at_q, const unsigned int n_quad_pts, const double p_avg)
unsigned int number_of_grid_elements_z
Number of subdivisions in z direction for a rectangle grid.
LimiterType bound_preserving_limiter
Variable to store specified limiter type.
void matrix_vector_mult(const std::vector< real > &input_vect, std::vector< real > &output_vect, const dealii::FullMatrix< double > &basis_x, const dealii::FullMatrix< double > &basis_y, const dealii::FullMatrix< double > &basis_z, const bool adding=false, const double factor=1.0) override
Computes a matrix-vector product using sum-factorization. Pass the one-dimensional basis...
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.
double min_density
Epsilon value for Positivity-Preserving Limiter.
PositivityPreservingLimiter(const Parameters::AllParameters *const parameters_input)
Constructor.
static std::shared_ptr< ManufacturedSolutionFunction< dim, real > > create_ManufacturedSolution(Parameters::AllParameters const *const param, int nstate)
Construct Manufactured solution object from global parameter file.
PartialDifferentialEquation
Possible Partial Differential Equations to solve.
const int nstate
Number of states.
Files for the baseline physics.
EulerParam euler_param
Contains parameters for the Euler equations non-dimensionalization.
Main parameter class that contains the various other sub-parameter classes.
std::array< real, nstate > get_soln_cell_avg(const std::array< std::vector< real >, nstate > &soln_at_q, const unsigned int n_quad_pts, const std::vector< real > &quad_weights)
Function to obtain the solution cell average.
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)
double ref_length
Reference length.
unsigned int number_of_grid_elements_x
Number of subdivisions in x direction for a rectangle grid.
real get_density_scaling_value(const double density_avg, const double density_min, const double pos_eps, const double p_avg)
real dz
Value required to compute solution cell average in 2D/3D, calculated using zmax and zmin parameters...
const Parameters::FlowSolverParam flow_solver_param
Flow solver parameters.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
Class for implementation of two forms of the Positivity-Preserving limiter derived from BoundPreservi...
unsigned int number_of_grid_elements_y
Number of subdivisions in y direction for a rectangle grid.
bool use_tvb_limiter
Flag for applying TVB Limiter.
std::array< real, nstate > get_soln_cell_avg_PPL(const std::array< std::array< std::vector< real >, nstate >, dim > &soln_at_q, const unsigned int n_quad_pts, const std::vector< real > &quad_weights_GLL, const std::vector< real > &quad_weights_GL, double &dt)
std::shared_ptr< Physics::Euler< dim, nstate, double > > euler_physics
Euler physics pointer. Used to compute pressure.
real dy
Value required to compute solution cell average in 2D/3D, calculated using ymax and ymin parameters...
const Parameters::AllParameters *const all_parameters
Pointer to parameters object.
real dx
Value required to compute solution cell average in 2D/3D, calculated using xmax and xmin parameters...