[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_tests.cpp
1 #include "positivity_preserving_tests.h"
2 #include "mesh/grids/positivity_preserving_tests_grid.h"
3 #include "mesh/grids/straight_periodic_cube.hpp"
4 #include <deal.II/grid/grid_generator.h>
5 #include "physics/physics_factory.h"
6 
7 namespace PHiLiP {
8 namespace FlowSolver {
9 
10 template <int dim, int nstate>
12  : CubeFlow_UniformGrid<dim, nstate>(parameters_input)
13  , unsteady_data_table_filename_with_extension(this->all_param.flow_solver_param.unsteady_data_table_filename+".txt")
14 {
15  this->euler_physics = std::dynamic_pointer_cast<Physics::Euler<dim,dim+2,double>>(
17 }
18 
19 template <int dim, int nstate>
20 std::shared_ptr<Triangulation> PositivityPreservingTests<dim,nstate>::generate_grid() const
21 {
22  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation> (
23  #if PHILIP_DIM!=1
24  this->mpi_communicator
25  #endif
26  );
27 
28  if(dim >= 1) {
29  if(this->all_param.flow_solver_param.grid_right_bound == this->all_param.flow_solver_param.grid_left_bound) {
30  std::cout << "Error: xmax and xmin need to be provided as parameters - Aborting... " << std::endl << std::flush;
31  std::abort();
32  }
33  }
34 
35  if(dim >= 2) {
36  if(this->all_param.flow_solver_param.grid_top_bound == this->all_param.flow_solver_param.grid_bottom_bound) {
37  std::cout << "Error: ymax and ymin need to be provided as parameters - Aborting... " << std::endl << std::flush;
38  std::abort();
39  }
40  }
41 
42  if(dim == 3) {
43  if(this->all_param.flow_solver_param.grid_z_upper_bound == this->all_param.flow_solver_param.grid_z_lower_bound) {
44  std::cout << "Error: zmax and zmin need to be provided as parameters - Aborting... " << std::endl << std::flush;
45  std::abort();
46  }
47  }
48  using flow_case_enum = Parameters::FlowSolverParam::FlowCaseType;
49  flow_case_enum flow_case_type = this->all_param.flow_solver_param.flow_case_type;
50 
51  if (dim==1 && (flow_case_type == flow_case_enum::sod_shock_tube
52  || flow_case_type == flow_case_enum::leblanc_shock_tube
53  || flow_case_type == flow_case_enum::shu_osher_problem)) {
54  Grids::shock_tube_1D_grid<dim>(*grid, &this->all_param.flow_solver_param);
55  }
56  else if (dim==2 && flow_case_type == flow_case_enum::shock_diffraction) {
57  Grids::shock_diffraction_grid<dim>(*grid, &this->all_param.flow_solver_param);
58  }
59  else if (dim==2 && flow_case_type == flow_case_enum::astrophysical_jet) {
60  Grids::astrophysical_jet_grid<dim>(*grid, &this->all_param.flow_solver_param);
61  }
62  else if (dim==2 && flow_case_type == flow_case_enum::double_mach_reflection) {
63  Grids::double_mach_reflection_grid<dim>(*grid, &this->all_param.flow_solver_param);
64  }
65  else if (dim==2 && flow_case_type == flow_case_enum::strong_vortex_shock_wave) {
66  Grids::svsw_grid<dim>(*grid, &this->all_param.flow_solver_param);
67  }
68  return grid;
69 }
70 
71 template <int dim, int nstate>
73 {
74  this->pcout << "- - Courant-Friedrichs-Lewy number: " << this->all_param.flow_solver_param.courant_friedrichs_lewy_number << std::endl;
75 }
76 
77 template<int dim, int nstate>
79 {
80  //create 1D solution polynomial basis functions and corresponding projection operator
81  //to interpolate the solution to the quadrature nodes, and to project it back to the
82  //modal coefficients.
83  const unsigned int init_grid_degree = dg.max_grid_degree;
84  const unsigned int poly_degree = this->all_param.flow_solver_param.poly_degree;
85  //Constructor for the operators
86  OPERATOR::basis_functions<dim, 2 * dim, double> soln_basis(1, poly_degree, init_grid_degree);
87  OPERATOR::vol_projection_operator<dim, 2 * dim, double> soln_basis_projection_oper(1, dg.max_degree, init_grid_degree);
88 
89 
90  // Build the oneD operator to perform interpolation/projection
91  soln_basis.build_1D_volume_operator(dg.oneD_fe_collection_1state[poly_degree], dg.oneD_quadrature_collection[poly_degree]);
92  soln_basis_projection_oper.build_1D_volume_operator(dg.oneD_fe_collection_1state[poly_degree], dg.oneD_quadrature_collection[poly_degree]);
93 
94  for (auto soln_cell = dg.dof_handler.begin_active(); soln_cell != dg.dof_handler.end(); ++soln_cell) {
95  if (!soln_cell->is_locally_owned()) continue;
96 
97 
98  std::vector<dealii::types::global_dof_index> current_dofs_indices;
99  // Current reference element related to this physical cell
100  const int i_fele = soln_cell->active_fe_index();
101  const dealii::FESystem<dim, dim>& current_fe_ref = dg.fe_collection[i_fele];
102  const int poly_degree = current_fe_ref.tensor_degree();
103 
104  const unsigned int n_dofs_curr_cell = current_fe_ref.n_dofs_per_cell();
105 
106  // Obtain the mapping from local dof indices to global dof indices
107  current_dofs_indices.resize(n_dofs_curr_cell);
108  soln_cell->get_dof_indices(current_dofs_indices);
109 
110  // Extract the local solution dofs in the cell from the global solution dofs
111  std::array<std::vector<double>, nstate> soln_coeff;
112  const unsigned int n_shape_fns = n_dofs_curr_cell / nstate;
113 
114  for (unsigned int istate = 0; istate < nstate; ++istate) {
115  soln_coeff[istate].resize(n_shape_fns);
116  }
117 
118  // Allocate solution dofs and set local max and min
119  for (unsigned int idof = 0; idof < n_dofs_curr_cell; ++idof) {
120  const unsigned int istate = dg.fe_collection[poly_degree].system_to_component_index(idof).first;
121  const unsigned int ishape = dg.fe_collection[poly_degree].system_to_component_index(idof).second;
122  soln_coeff[istate][ishape] = dg.solution[current_dofs_indices[idof]];
123  }
124 
125  const unsigned int n_quad_pts = dg.volume_quadrature_collection[poly_degree].size();
126 
127  std::array<std::vector<double>, nstate> soln_at_q;
128 
129  // Interpolate solution dofs to quadrature pts.
130  for (int istate = 0; istate < nstate; istate++) {
131  soln_at_q[istate].resize(n_quad_pts);
132  soln_basis.matrix_vector_mult_1D(soln_coeff[istate], soln_at_q[istate], soln_basis.oneD_vol_operator);
133  }
134 
135  for (unsigned int iquad = 0; iquad < n_quad_pts; ++iquad) {
136  // Verify that positivity of density is preserved
137  if (soln_at_q[0][iquad] < 0) {
138  std::cout << "Flow Solver Error: Density is negative - Aborting... " << std::endl << std::flush;
139  std::abort();
140  }
141 
142  if (soln_at_q[nstate - 1][iquad] < 0) {
143  std::cout << "Flow Solver Error: Total Energy is negative - Aborting... " << std::endl << std::flush;
144  std::abort();
145  }
146 
147  if ((isnan(soln_at_q[0][iquad]))) {
148  std::cout << "Flow Solver Error: Density is NaN - Aborting... " << std::endl << std::flush;
149  std::abort();
150  }
151  }
152  }
153 }
154 
155 template<int dim, int nstate>
157 {
158  // Check that poly_degree is uniform everywhere
159  if (dg.get_max_fe_degree() != dg.get_min_fe_degree()) {
160  // Note: This function may have issues with nonuniform p. Should test in debug mode if developing in the future.
161  this->pcout << "ERROR: compute_integrated_quantities() is untested for nonuniform p. Aborting..." << std::endl;
162  std::abort();
163  }
164 
165  double integrated_quantity = 0.0;
166 
167  const unsigned int grid_degree = dg.high_order_grid->fe_system.tensor_degree();
168  const unsigned int poly_degree = dg.max_degree;
169 
170  // Set the quadrature of size dim and 1D for sum-factorization.
171  dealii::Quadrature<1> quad_1D = dg.oneD_quadrature_collection[poly_degree];
172  std::vector<double> quad_weights = dg.volume_quadrature_collection[poly_degree].get_weights();
173  unsigned int n_quad_pts = dg.volume_quadrature_collection[poly_degree].size();
174 
175  // Construct the basis functions and mapping shape functions.
176  OPERATOR::basis_functions<dim,2*dim,double> soln_basis(1, poly_degree, grid_degree);
177  OPERATOR::mapping_shape_functions<dim,2*dim,double> mapping_basis(1, poly_degree, grid_degree);
178  // Build basis function volume operator from 1D finite element for 1 state.
179  soln_basis.build_1D_volume_operator(dg.oneD_fe_collection_1state[poly_degree], quad_1D);
180  // Build mapping shape functions operators using the oneD high_ordeR_grid finite element
181  mapping_basis.build_1D_shape_functions_at_grid_nodes(dg.high_order_grid->oneD_fe_system, dg.high_order_grid->oneD_grid_nodes);
182  mapping_basis.build_1D_shape_functions_at_flux_nodes(dg.high_order_grid->oneD_fe_system, quad_1D, dg.oneD_face_quadrature);
183  // If in the future we need the physical quadrature node location, turn these flags to true and the constructor will
184  // automatically compute it for you. Currently set to false as to not compute extra unused terms.
185  const bool store_vol_flux_nodes = false;//currently doesn't need the volume physical nodal position
186  const bool store_surf_flux_nodes = false;//currently doesn't need the surface physical nodal position
187 
188  const unsigned int n_dofs = dg.fe_collection[poly_degree].n_dofs_per_cell();
189  const unsigned int n_shape_fns = n_dofs / nstate;
190  std::vector<dealii::types::global_dof_index> dofs_indices (n_dofs);
191  auto metric_cell = dg.high_order_grid->dof_handler_grid.begin_active();
192  // Changed for loop to update metric_cell.
193  for (auto cell = dg.dof_handler.begin_active(); cell!= dg.dof_handler.end(); ++cell, ++metric_cell) {
194  if (!cell->is_locally_owned()) continue;
195  cell->get_dof_indices (dofs_indices);
196 
197  // We first need to extract the mapping support points (grid nodes) from high_order_grid.
198  const dealii::FESystem<dim> &fe_metric = dg.high_order_grid->fe_system;
199  const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
200  const unsigned int n_grid_nodes = n_metric_dofs / dim;
201  std::vector<dealii::types::global_dof_index> metric_dof_indices(n_metric_dofs);
202  metric_cell->get_dof_indices (metric_dof_indices);
203  std::array<std::vector<double>,dim> mapping_support_points;
204  for(int idim=0; idim<dim; idim++){
205  mapping_support_points[idim].resize(n_grid_nodes);
206  }
207  // Get the mapping support points (physical grid nodes) from high_order_grid.
208  // Store it in such a way we can use sum-factorization on it with the mapping basis functions.
209  const std::vector<unsigned int > &index_renumbering = dealii::FETools::hierarchic_to_lexicographic_numbering<dim>(grid_degree);
210  for (unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
211  const double val = (dg.high_order_grid->volume_nodes[metric_dof_indices[idof]]);
212  const unsigned int istate = fe_metric.system_to_component_index(idof).first;
213  const unsigned int ishape = fe_metric.system_to_component_index(idof).second;
214  const unsigned int igrid_node = index_renumbering[ishape];
215  mapping_support_points[istate][igrid_node] = val;
216  }
217  // Construct the metric operators.
218  OPERATOR::metric_operators<double, dim, 2*dim> metric_oper(nstate, poly_degree, grid_degree, store_vol_flux_nodes, store_surf_flux_nodes);
219  // Build the metric terms to compute the gradient and volume node positions.
220  // This functions will compute the determinant of the metric Jacobian and metric cofactor matrix.
221  // If flags store_vol_flux_nodes and store_surf_flux_nodes set as true it will also compute the physical quadrature positions.
222  metric_oper.build_volume_metric_operators(
223  n_quad_pts, n_grid_nodes,
224  mapping_support_points,
225  mapping_basis,
226  dg.all_parameters->use_invariant_curl_form);
227 
228  // Fetch the modal soln coefficients
229  // We immediately separate them by state as to be able to use sum-factorization
230  // in the interpolation operator. If we left it by n_dofs_cell, then the matrix-vector
231  // mult would sum the states at the quadrature point.
232  // That is why the basis functions are based off the 1state oneD fe_collection.
233  std::array<std::vector<double>,nstate> soln_coeff;
234  for (unsigned int idof = 0; idof < n_dofs; ++idof) {
235  const unsigned int istate = dg.fe_collection[poly_degree].system_to_component_index(idof).first;
236  const unsigned int ishape = dg.fe_collection[poly_degree].system_to_component_index(idof).second;
237  if(ishape == 0){
238  soln_coeff[istate].resize(n_shape_fns);
239  }
240 
241  soln_coeff[istate][ishape] = dg.solution(dofs_indices[idof]);
242  }
243  // Interpolate each state to the quadrature points using sum-factorization
244  // with the basis functions in each reference direction.
245  std::array<std::vector<double>,nstate> soln_at_q_vect;
246  for(int istate=0; istate<nstate; istate++){
247  soln_at_q_vect[istate].resize(n_quad_pts);
248  // Interpolate soln coeff to volume cubature nodes.
249  soln_basis.matrix_vector_mult_1D(soln_coeff[istate], soln_at_q_vect[istate],
250  soln_basis.oneD_vol_operator);
251  }
252 
253  // Loop over quadrature nodes, compute quantities to be integrated, and integrate them.
254  for (unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
255 
256  std::array<double,nstate> soln_at_q;
257  // Extract solution in a way that the physics ca n use them.
258  for(int istate=0; istate<nstate; istate++){
259  soln_at_q[istate] = soln_at_q_vect[istate][iquad];
260  }
261 
262  //#####################################################################
263  // Compute integrated quantities here
264  //#####################################################################
265  const double quadrature_entropy = this->euler_physics->compute_numerical_entropy_function(soln_at_q);
266  //Using std::cout because of cell->is_locally_owned check
267  if (isnan(quadrature_entropy)) std::cout << "WARNING: NaN entropy detected at a node!" << std::endl;
268  integrated_quantity += quadrature_entropy * quad_weights[iquad] * metric_oper.det_Jac_vol[iquad];
269  //#####################################################################
270  }
271  }
272 
273  //MPI
274  integrated_quantity = dealii::Utilities::MPI::sum(integrated_quantity, this->mpi_communicator);
275 
276  return integrated_quantity;
277 }
278 
279 template <int dim, int nstate>
281  const std::shared_ptr<ODE::ODESolverBase<dim, double>> ode_solver,
282  const std::shared_ptr <DGBase<dim, double>> dg,
283  const std::shared_ptr <dealii::TableHandler> unsteady_data_table)
284 {
285  //unpack current iteration and current time from ode solver
286  const unsigned int current_iteration = ode_solver->current_iteration;
287  const double current_time = ode_solver->current_time;
288 
289  // All discrete proofs use solution nodes, therefore it is best to report
290  // entropy on the solution nodes rather than by overintegrating.
291  const double current_numerical_entropy = this->compute_integrated_entropy(*dg); // no overintegration
292  if (current_iteration==0) this->previous_numerical_entropy = current_numerical_entropy;
293  const double entropy = current_numerical_entropy - previous_numerical_entropy + ode_solver->FR_entropy_contribution_RRK_solver;
294  this->previous_numerical_entropy = current_numerical_entropy;
295 
296  if (std::isnan(entropy)){
297  this->pcout << "Entropy is nan. Aborting flow simulation..." << std::endl << std::flush;
298  std::abort();
299  }
300  if (current_iteration == 0) initial_entropy = current_numerical_entropy;
301 
302  if(nstate == dim + 2)
303  this->check_positivity_density(*dg);
304 
305  if (this->mpi_rank == 0) {
306 
307  unsteady_data_table->add_value("iteration", current_iteration);
308  // Add values to data table
309  this->add_value_to_data_table(current_time, "time", unsteady_data_table);
310  this->add_value_to_data_table(entropy,"entropy",unsteady_data_table);
311  unsteady_data_table->set_scientific("entropy", false);
312  this->add_value_to_data_table(current_numerical_entropy,"current_numerical_entropy",unsteady_data_table);
313  unsteady_data_table->set_scientific("current_numerical_entropy", false);
314  this->add_value_to_data_table(entropy/initial_entropy,"U/Uo",unsteady_data_table);
315  unsteady_data_table->set_scientific("U/Uo", false);
316 
317 
318  // Write to file
319  std::ofstream unsteady_data_table_file(this->unsteady_data_table_filename_with_extension);
320  unsteady_data_table->write_text(unsteady_data_table_file);
321  }
322 
323  if (current_iteration % this->all_param.ode_solver_param.print_iteration_modulo == 0) {
324  // Print to console
325  this->pcout << " Iter: " << current_iteration
326  << " Time: " << std::setprecision(16) << current_time
327  << " Current Numerical Entropy: " << current_numerical_entropy
328  << " Entropy: " << entropy
329  << " (U-Uo)/Uo: " << entropy/initial_entropy;
330 
331  this->pcout << std::endl;
332  }
333 
334  // Update local maximum wave speed before calculating next time step
336 }
337 
339 
340 } // FlowSolver namespace
341 } // PHiLiP namespace
const Parameters::AllParameters all_param
All parameters.
const MPI_Comm mpi_communicator
MPI communicator.
void check_positivity_density(DGBase< dim, double > &dg)
Check positivity of density and total energy + verify that density is not NaN.
FlowCaseType
Selects the flow case to be simulated.
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
Definition: operators.cpp:1082
void compute_unsteady_data_and_write_to_table(const std::shared_ptr< ODE::ODESolverBase< dim, double >> ode_solver, const std::shared_ptr< DGBase< dim, double >> dg, const std::shared_ptr< dealii::TableHandler > unsteady_data_table) override
Compute the desired unsteady data and write it to a table.
FlowCaseType flow_case_type
Selected FlowCaseType from the input file.
double courant_friedrichs_lewy_number
Courant-Friedrich-Lewy (CFL) number for constant time step.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
double initial_entropy
Storing entropy at first step.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: dg_base.hpp:91
void build_1D_volume_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
Definition: operators.cpp:1732
dealii::QGauss< 0 > oneD_face_quadrature
1D surface quadrature is always one single point for all poly degrees.
Definition: dg_base.hpp:634
Files for the baseline physics.
Definition: ADTypes.hpp:10
Base class ODE solver.
unsigned int print_iteration_modulo
If ode_output==verbose, print every print_iteration_modulo iterations.
double compute_integrated_entropy(DGBase< dim, double > &dg) const
Updates the maximum local wave speed.
void build_1D_shape_functions_at_grid_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Constructs the volume operator and gradient operator.
Definition: operators.cpp:2155
unsigned int poly_degree
Polynomial order (P) of the basis functions for DG.
Main parameter class that contains the various other sub-parameter classes.
const dealii::hp::FECollection< 1 > oneD_fe_collection_1state
1D Finite Element Collection for p-finite-element to represent the solution for a single state...
Definition: dg_base.hpp:627
double grid_z_upper_bound
Maximum z bound of domain for a rectangle grid.
dealii::LinearAlgebra::distributed::Vector< double > solution
Current modal coefficients of the solution.
Definition: dg_base.hpp:409
std::vector< real > det_Jac_vol
The determinant of the metric Jacobian at volume cubature nodes.
Definition: operators.h:1171
void add_value_to_data_table(const double value, const std::string value_string, const std::shared_ptr< dealii::TableHandler > data_table) const
Add a value to a given data table with scientific format.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
PositivityPreservingTests(const Parameters::AllParameters *const parameters_input)
Constructor.
void update_maximum_local_wave_speed(DGBase< dim, double > &dg)
Updates the maximum local wave speed.
const unsigned int max_degree
Maximum degree used for p-refi1nement.
Definition: dg_base.hpp:104
Base metric operators class that stores functions used in both the volume and on surface.
Definition: operators.h:1086
double grid_top_bound
Maximum y bound of domain for a rectangle grid.
unsigned int get_max_fe_degree()
Gets the maximum value of currently active FE degree.
Definition: dg_base.cpp:305
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
Definition: operators.h:1026
Euler equations. Derived from PhysicsBase.
Definition: euler.h:78
std::shared_ptr< HighOrderGrid< dim, real, MeshType > > high_order_grid
High order grid that will provide the MappingFEField.
Definition: dg_base.hpp:655
const unsigned int max_grid_degree
Maximum grid degree used for hp-refi1nement.
Definition: dg_base.hpp:109
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
Definition: operators.h:355
static std::shared_ptr< PhysicsBase< dim, nstate, real > > create_Physics(const Parameters::AllParameters *const parameters_input, std::shared_ptr< ModelBase< dim, nstate, real > > model_input=nullptr)
Factory to return the correct physics given input file.
unsigned int get_min_fe_degree()
Gets the minimum value of currently active FE degree.
Definition: dg_base.cpp:317
dealii::DoFHandler< dim > dof_handler
Finite Element Collection to represent the high-order grid.
Definition: dg_base.hpp:652
const std::string unsteady_data_table_filename_with_extension
Filename (with extension) for the unsteady data table.
void build_volume_metric_operators(const unsigned int n_quad_pts, const unsigned int n_metric_dofs, const std::array< std::vector< real >, dim > &mapping_support_points, mapping_shape_functions< dim, n_faces, real > &mapping_basis, const bool use_invariant_curl_form=false)
Builds the volume metric operators.
Definition: operators.cpp:2294
double grid_right_bound
Right bound of domain for hyper_cube mesh based cases.
dealii::hp::QCollection< 1 > oneD_quadrature_collection
1D quadrature to generate Lagrange polynomials for the sake of flux interpolation.
Definition: dg_base.hpp:632
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
dealii::ConditionalOStream pcout
ConditionalOStream.
dealii::hp::QCollection< dim > volume_quadrature_collection
Finite Element Collection to represent the high-order grid.
Definition: dg_base.hpp:608
void build_1D_shape_functions_at_flux_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature, const dealii::Quadrature< 0 > &face_quadrature)
Constructs the volume, gradient, surface, and surface gradient operator.
Definition: operators.cpp:2164
void matrix_vector_mult_1D(const std::vector< real > &input_vect, std::vector< real > &output_vect, const dealii::FullMatrix< double > &basis_x, const bool adding=false, const double factor=1.0)
Apply the matrix vector operation using the 1D operator in each direction.
Definition: operators.cpp:308
void display_additional_flow_case_specific_parameters() const override
Display additional more specific flow case parameters.
const dealii::hp::FECollection< dim > fe_collection
Finite Element Collection for p-finite-element to represent the solution.
Definition: dg_base.hpp:597
Projection operator corresponding to basis functions onto M-norm (L2).
Definition: operators.h:694
std::shared_ptr< Triangulation > generate_grid() const override
Function to generate the grid.