[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
h_refinement_study_isentropic_vortex.cpp
1 #include "h_refinement_study_isentropic_vortex.h"
2 #include "flow_solver/flow_solver_factory.h"
3 #include "flow_solver/flow_solver_cases/periodic_1D_unsteady.h"
4 #include "flow_solver/flow_solver_cases/periodic_entropy_tests.h"
5 #include "physics/exact_solutions/exact_solution.h"
6 #include "physics/euler.h"
7 #include "cmath"
8 //#include "ode_solver/runge_kutta_ode_solver.h"
9 
10 namespace PHiLiP {
11 namespace Tests {
12 
13 template <int dim, int nstate>
15  const PHiLiP::Parameters::AllParameters *const parameters_input,
16  const dealii::ParameterHandler &parameter_handler_input)
17  : TestsBase::TestsBase(parameters_input),
18  parameter_handler(parameter_handler_input),
19  n_calculations(parameters_input->time_refinement_study_param.number_of_times_to_solve),
20  refine_ratio(parameters_input->time_refinement_study_param.refinement_ratio)
21 {}
22 
23 
24 
25 template <int dim, int nstate>
27 {
29 
30  parameters.flow_solver_param.number_of_grid_elements_per_dimension *= pow(2, refinement);
32  return parameters;
33 }
34 
35 template <int dim, int nstate>
37  double &Lp_error_pressure,
38  std::shared_ptr<DGBase<dim,double>> dg,
39  const Parameters::AllParameters parameters,
40  double final_time,
41  int norm_p) const
42 {
43  //generate exact solution at final time
44  std::shared_ptr<ExactSolutionFunction<dim,nstate,double>> exact_solution_function;
45  exact_solution_function = ExactSolutionFactory<dim,nstate,double>::create_ExactSolutionFunction(parameters.flow_solver_param, final_time);
46  int overintegrate = 10;
47 
48  // For Euler, compare only density or pressure
49  // deal.ii compute_global_error() does not interface simply
50  // Therefore, overintegrated error calculation is coded here
51 
52  // Need to do MPI sum manaully, so declaring a new local error
53  double Lp_error_density_local = 0;
54  double Lp_error_pressure_local = 0;
55 
56  std::shared_ptr< Physics::Euler<dim,dim+2,double> > euler_physics = std::dynamic_pointer_cast<Physics::Euler<dim,dim+2,double>>(
58  dealii::QGauss<dim> quad_extra(dg->max_degree+1+overintegrate);
59  dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[dg->max_degree], quad_extra,
60  dealii::update_values | dealii::update_gradients | dealii::update_JxW_values | dealii::update_quadrature_points);
61 
62  const unsigned int n_quad_pts = fe_values_extra.n_quadrature_points;
63  std::array<double,nstate> soln_at_q;
64 
65  std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
66  for (auto cell : dg->dof_handler.active_cell_iterators()) {
67  if (!cell->is_locally_owned()) continue;
68  fe_values_extra.reinit (cell);
69  cell->get_dof_indices (dofs_indices);
70 
71  for (unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
72 
73  std::fill(soln_at_q.begin(), soln_at_q.end(), 0.0);
74 
75  for (unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof) {
76  const unsigned int istate = fe_values_extra.get_fe().system_to_component_index(idof).first;
77  soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, istate);
78  }
79 
80  const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
81 
82  //Compute Lp error for istate=0, which is density
83  std::array<double,nstate> exact_soln_at_q;
84  for (unsigned int istate = 0; istate < nstate; ++istate) {
85  exact_soln_at_q[istate] = exact_solution_function->value(qpoint, istate);
86  }
87 
88  const double exact_pressure = euler_physics->compute_pressure(exact_soln_at_q);
89  const double pressure = euler_physics->compute_pressure(soln_at_q);
90 
91  if (norm_p > 0){
92  Lp_error_density_local += pow(abs(soln_at_q[0] - exact_soln_at_q[0]), norm_p) * fe_values_extra.JxW(iquad);
93  Lp_error_pressure_local += pow(abs(pressure-exact_pressure), norm_p) * fe_values_extra.JxW(iquad);
94  } else{
95  //L-infinity norm
96  Lp_error_density_local = std::max(abs(soln_at_q[0]-exact_soln_at_q[0]), Lp_error_density_local);
97  Lp_error_pressure_local = std::max(abs(pressure-exact_pressure), Lp_error_pressure_local);
98  }
99  }
100  }
101 
102  //MPI sum
103  if (norm_p > 0) {
104  Lp_error_density= dealii::Utilities::MPI::sum(Lp_error_density_local, this->mpi_communicator);
105  Lp_error_density= pow(Lp_error_density, 1.0/((double)norm_p));
106  Lp_error_pressure = dealii::Utilities::MPI::sum(Lp_error_pressure_local, this->mpi_communicator);
107  Lp_error_pressure= pow(Lp_error_pressure, 1.0/((double)norm_p));
108  } else {
109  // L-infinity norm
110  Lp_error_density = dealii::Utilities::MPI::max(Lp_error_density_local, this->mpi_communicator);
111  Lp_error_pressure = dealii::Utilities::MPI::max(Lp_error_pressure_local, this->mpi_communicator);
112  }
113 
114 }
115 
116 template <int dim, int nstate>
118 {
119 
120  const double final_time = this->all_parameters->flow_solver_param.final_time;
121  const double initial_time_step = this->all_parameters->ode_solver_param.initial_time_step;
122  const int n_steps = floor(final_time/initial_time_step);
123  if (n_steps * initial_time_step != final_time){
124  pcout << "WARNING: final_time is not evenly divisible by initial_time_step!" << std::endl
125  << "Remainder is " << fmod(final_time, initial_time_step)
126  << ". Consider modifying parameters." << std::endl;
127  }
128 
129  int testfail = 0;
130 
131  dealii::ConvergenceTable convergence_table_density;
132  dealii::ConvergenceTable convergence_table_pressure;
133  double L2_error_pressure_old = 0;
134  double L2_error_pressure_conv_rate=0;
135 
136 
137  for (int refinement = 0; refinement < n_calculations; ++refinement){
138 
139  pcout << "\n\n---------------------------------------------" << std::endl;
140  pcout << "Refinement number " << refinement << " of " << n_calculations - 1 << std::endl;
141  pcout << "---------------------------------------------" << std::endl;
142 
143  const Parameters::AllParameters params = reinit_params_and_refine(refinement);
144  std::unique_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver = FlowSolver::FlowSolverFactory<dim,nstate>::select_flow_case(&params, parameter_handler);
145  std::unique_ptr<FlowSolver::PeriodicEntropyTests<dim, nstate>> flow_solver_case = std::make_unique<FlowSolver::PeriodicEntropyTests<dim,nstate>>(&params);
146 
147  static_cast<void>(flow_solver->run());
148  pcout << "Finished flowsolver " << std::endl;
149 
150  const double final_time_actual = flow_solver->ode_solver->current_time;
151 
152  double L1_error_density=0;
153  double L1_error_pressure=0;
154  calculate_Lp_error_at_final_time_wrt_function(L1_error_density, L1_error_pressure, flow_solver->dg, params,final_time_actual, 1);
155  double L2_error_density =0;
156  double L2_error_pressure=0;
157  calculate_Lp_error_at_final_time_wrt_function(L2_error_density,L2_error_pressure, flow_solver->dg, params,final_time_actual, 2);
158  double Linfty_error_density = 0;
159  double Linfty_error_pressure = 0;
160  calculate_Lp_error_at_final_time_wrt_function(Linfty_error_density,Linfty_error_pressure, flow_solver->dg, params,final_time_actual, -1);
161  pcout << "Computed density errors are: " << std::endl
162  << " L1: " << L1_error_density << std::endl
163  << " L2: " << L2_error_density << std::endl
164  << " Linfty: " << Linfty_error_density << std::endl;
165 
166  const double dt = flow_solver_case->get_constant_time_step(flow_solver->dg);
167  const int n_cells = pow(params.flow_solver_param.number_of_grid_elements_per_dimension, PHILIP_DIM);
168  pcout << " at dt = " << dt << std::endl;
169 
170  // Convergence for density
171  convergence_table_density.add_value("refinement", refinement);
172  convergence_table_density.add_value("dt", dt );
173  convergence_table_density.set_precision("dt", 16);
174  convergence_table_density.set_scientific("dt", true);
175  convergence_table_density.add_value("n_cells",n_cells);
176  convergence_table_density.add_value("L1_error_density",L1_error_density);
177  convergence_table_density.set_precision("L1_error_density", 16);
178  convergence_table_density.evaluate_convergence_rates("L1_error_density", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
179  convergence_table_density.add_value("L2_error_density",L2_error_density);
180  convergence_table_density.set_precision("L2_error_density", 16);
181  convergence_table_density.evaluate_convergence_rates("L2_error_density", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
182  convergence_table_density.add_value("Linfty_error_density",Linfty_error_density);
183  convergence_table_density.set_precision("Linfty_error_density", 16);
184  convergence_table_density.evaluate_convergence_rates("Linfty_error_density", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
185 
186  if (params.ode_solver_param.ode_solver_type == Parameters::ODESolverParam::ODESolverEnum::rrk_explicit_solver) {
187  const double gamma_agg = final_time_actual / (dt * flow_solver->ode_solver->current_iteration);
188 
189  convergence_table_density.add_value("gamma_agg",gamma_agg-1.0);
190  convergence_table_density.set_precision("gamma_agg", 16);
191  convergence_table_density.evaluate_convergence_rates("gamma_agg", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
192  }
193 
194  // Convergence for pressure
195  convergence_table_pressure.add_value("refinement", refinement);
196  convergence_table_pressure.add_value("dt", dt );
197  convergence_table_pressure.set_precision("dt", 16);
198  convergence_table_pressure.set_scientific("dt", true);
199  convergence_table_pressure.add_value("n_cells",n_cells);
200  convergence_table_pressure.add_value("L1_error_pressure",L1_error_pressure);
201  convergence_table_pressure.set_precision("L1_error_pressure", 16);
202  convergence_table_pressure.evaluate_convergence_rates("L1_error_pressure", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
203  convergence_table_pressure.add_value("L2_error_pressure",L2_error_pressure);
204  convergence_table_pressure.set_precision("L2_error_pressure", 16);
205  convergence_table_pressure.evaluate_convergence_rates("L2_error_pressure", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
206  convergence_table_pressure.add_value("Linfty_error_pressure",Linfty_error_pressure);
207  convergence_table_pressure.set_precision("Linfty_error_pressure", 16);
208  convergence_table_pressure.evaluate_convergence_rates("Linfty_error_pressure", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
209 
210  //Checking convergence order
211  const double expected_order = params.flow_solver_param.poly_degree + 1;
212  //set tolerance to make test pass for ctest. Note that the grids are very coarse (not in asymptotic range)
213  const double order_tolerance = 1.0;
214  if (refinement > 0) {
215  L2_error_pressure_conv_rate = -log(L2_error_pressure_old/L2_error_pressure)/log(refine_ratio);
216  pcout << "Order for L2 pressure at " << refinement << " is " << L2_error_pressure_conv_rate << std::endl;
217  if (abs(L2_error_pressure_conv_rate - expected_order) > order_tolerance){
218  testfail = 1;
219  pcout << "Expected convergence order for L2 pressure was not reached at refinement " << refinement <<std::endl;
220  }
221  if (refinement < n_calculations-1 && pcout.is_active()){
222  // Print current convergence results for solution monitoring
223  convergence_table_density.write_text(pcout.get_stream());
224  convergence_table_pressure.write_text(pcout.get_stream());
225  }
226  }
227  L2_error_pressure_old = L2_error_pressure;
228  }
229 
230  //Printing and writing convergence tables
231  pcout << std::endl;
232  if (pcout.is_active()){
233  convergence_table_density.write_text(pcout.get_stream());
234  convergence_table_pressure.write_text(pcout.get_stream());
235  }
236  std::ofstream conv_tab_file;
237  const std::string fname = "convergence_tables.txt";
238  conv_tab_file.open(fname);
239  convergence_table_density.write_text(conv_tab_file);
240  convergence_table_pressure.write_text(conv_tab_file);
241  conv_tab_file.close();
242 
243  return testfail;
244 }
245 #if PHILIP_DIM!=1
247 #endif
248 } // Tests namespace
249 } // PHiLiP namespace
double final_time
Final solution time.
unsigned int number_of_grid_elements_per_dimension
Number of grid elements per dimension for hyper_cube mesh based cases.
h refinement test for the isentropic vortex advection test case.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
const MPI_Comm mpi_communicator
MPI communicator.
Definition: tests.h:39
const int n_calculations
Number of times to solve for convergence summary.
Files for the baseline physics.
Definition: ADTypes.hpp:10
static std::shared_ptr< ExactSolutionFunction< dim, nstate, real > > create_ExactSolutionFunction(const Parameters::FlowSolverParam &flow_solver_parameters, const double time_compare)
Construct ExactSolutionFunction object from global parameter file.
unsigned int poly_degree
Polynomial order (P) of the basis functions for DG.
Main parameter class that contains the various other sub-parameter classes.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
double initial_time_step
Time step used in ODE solver.
static std::unique_ptr< FlowSolver< dim, nstate > > select_flow_case(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Factory to return the correct flow solver given input file.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
Definition: tests.h:20
Parameters::AllParameters reinit_params_and_refine(int refinement) const
Reinitialize parameters while refining the timestep. Necessary because all_parameters is constant...
Euler equations. Derived from PhysicsBase.
Definition: euler.h:78
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.
ODESolverEnum ode_solver_type
ODE solver type.
dealii::ConditionalOStream pcout
ConditionalOStream.
Definition: tests.h:45
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
void calculate_Lp_error_at_final_time_wrt_function(double &Lp_error_density, double &Lp_error_pressure, std::shared_ptr< DGBase< dim, double >> dg, const Parameters::AllParameters parameters, double final_time, int norm_p) const
Base class of all the tests.
Definition: tests.h:17
HRefinementStudyIsentropicVortex(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.