[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 
33  if (nb_c_value > 0){
35  }
36 
37  return parameters;
38 }
39 
40 template <int dim, int nstate>
42  double &Lp_error_pressure,
43  std::shared_ptr<DGBase<dim,double>> dg,
44  const Parameters::AllParameters parameters,
45  double final_time,
46  int norm_p) const
47 {
48  //generate exact solution at final time
49  std::shared_ptr<ExactSolutionFunction<dim,nstate,double>> exact_solution_function;
50  exact_solution_function = ExactSolutionFactory<dim,nstate,double>::create_ExactSolutionFunction(parameters.flow_solver_param, final_time);
51  int overintegrate = 10;
52 
53  // For Euler, compare only density or pressure
54  // deal.ii compute_global_error() does not interface simply
55  // Therefore, overintegrated error calculation is coded here
56 
57  // Need to do MPI sum manaully, so declaring a new local error
58  double Lp_error_density_local = 0;
59  double Lp_error_pressure_local = 0;
60 
61  std::shared_ptr< Physics::Euler<dim,dim+2,double> > euler_physics = std::dynamic_pointer_cast<Physics::Euler<dim,dim+2,double>>(
63  dealii::QGauss<dim> quad_extra(dg->max_degree+1+overintegrate);
64  dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->fe_collection[dg->max_degree], quad_extra,
65  dealii::update_values | dealii::update_gradients | dealii::update_JxW_values | dealii::update_quadrature_points);
66 
67  const unsigned int n_quad_pts = fe_values_extra.n_quadrature_points;
68  std::array<double,nstate> soln_at_q;
69 
70  std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
71  for (auto cell : dg->dof_handler.active_cell_iterators()) {
72  if (!cell->is_locally_owned()) continue;
73  fe_values_extra.reinit (cell);
74  cell->get_dof_indices (dofs_indices);
75 
76  for (unsigned int iquad=0; iquad<n_quad_pts; ++iquad) {
77 
78  std::fill(soln_at_q.begin(), soln_at_q.end(), 0.0);
79 
80  for (unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof) {
81  const unsigned int istate = fe_values_extra.get_fe().system_to_component_index(idof).first;
82  soln_at_q[istate] += dg->solution[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, istate);
83  }
84 
85  const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
86 
87  //Compute Lp error for istate=0, which is density
88  std::array<double,nstate> exact_soln_at_q;
89  for (unsigned int istate = 0; istate < nstate; ++istate) {
90  exact_soln_at_q[istate] = exact_solution_function->value(qpoint, istate);
91  }
92 
93  const double exact_pressure = euler_physics->compute_pressure(exact_soln_at_q);
94  const double pressure = euler_physics->compute_pressure(soln_at_q);
95 
96  if (norm_p > 0){
97  Lp_error_density_local += pow(abs(soln_at_q[0] - exact_soln_at_q[0]), norm_p) * fe_values_extra.JxW(iquad);
98  Lp_error_pressure_local += pow(abs(pressure-exact_pressure), norm_p) * fe_values_extra.JxW(iquad);
99  } else{
100  //L-infinity norm
101  Lp_error_density_local = std::max(abs(soln_at_q[0]-exact_soln_at_q[0]), Lp_error_density_local);
102  Lp_error_pressure_local = std::max(abs(pressure-exact_pressure), Lp_error_pressure_local);
103  }
104  }
105  }
106 
107  //MPI sum
108  if (norm_p > 0) {
109  Lp_error_density= dealii::Utilities::MPI::sum(Lp_error_density_local, this->mpi_communicator);
110  Lp_error_density= pow(Lp_error_density, 1.0/((double)norm_p));
111  Lp_error_pressure = dealii::Utilities::MPI::sum(Lp_error_pressure_local, this->mpi_communicator);
112  Lp_error_pressure= pow(Lp_error_pressure, 1.0/((double)norm_p));
113  } else {
114  // L-infinity norm
115  Lp_error_density = dealii::Utilities::MPI::max(Lp_error_density_local, this->mpi_communicator);
116  Lp_error_pressure = dealii::Utilities::MPI::max(Lp_error_pressure_local, this->mpi_communicator);
117  }
118 
119 }
120 
121 template <int dim, int nstate>
123 {
124 
125  const double final_time = this->all_parameters->flow_solver_param.final_time;
126  const double initial_time_step = this->all_parameters->ode_solver_param.initial_time_step;
127  const int n_steps = floor(final_time/initial_time_step);
128  if (n_steps * initial_time_step != final_time){
129  pcout << "WARNING: final_time is not evenly divisible by initial_time_step!" << std::endl
130  << "Remainder is " << fmod(final_time, initial_time_step)
131  << ". Consider modifying parameters." << std::endl;
132  }
133 
134  int testfail = 0;
135  dealii::ConvergenceTable convergence_table_density;
136  dealii::ConvergenceTable convergence_table_pressure;
137 
138  std::ofstream conv_tab_file;
139  const std::string fname = "convergence_tables.txt";
140  conv_tab_file.open(fname);
141 // cESFR range test -------------------------------------
142  const unsigned int nb_c_value = this->all_parameters->flow_solver_param.number_ESFR_parameter_values;
144  const double c_max = this->all_parameters->flow_solver_param.ESFR_parameter_values_end;
145  const double log_c_min = std::log10(c_min);
146  const double log_c_max = std::log10(c_max);
147  std::vector<double> c_array(nb_c_value+1);
148 
149 
150  // Create log space array of c_value
151  for (unsigned int ic = 0; ic < nb_c_value; ic++) {
152  double log_c = log_c_min + (log_c_max - log_c_min) / (nb_c_value - 1) * ic;
153  c_array[ic] = std::pow(10.0, log_c);
154  }
155  if (this->all_parameters->flux_reconstruction_type == PHiLiP::Parameters::AllParameters::Flux_Reconstruction::user_specified_value){
156  // Add cPlus value at the end
158  }
159  else{
160  // default value to ensure it is not empty, it won't be used later.
161  c_array[nb_c_value] = 0.0;
162  }
163 
164  // Loop over c_array to compute slope
165  for (unsigned int ic = 0; ic < nb_c_value+1; ic++) {
166  double c_value = c_array[ic];
167 // --------------------------------------------------------
168 
169  double L2_error_pressure_old = 0;
170  double L2_error_pressure_conv_rate=0;
171 
172 
173  for (int refinement = 0; refinement < n_calculations; ++refinement){
174 
175  pcout << "\n\n---------------------------------------------" << std::endl;
176  pcout << "Refinement number " << refinement + 1 << " of " << n_calculations << ", flux reconstruction parameter c = " << c_value << std::endl;
177  pcout << "---------------------------------------------" << std::endl;
178 
179  const Parameters::AllParameters params = reinit_params_and_refine(refinement, c_value, nb_c_value);
180  auto params_modified = params;
181 
182  std::unique_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver = FlowSolver::FlowSolverFactory<dim,nstate>::select_flow_case(&params_modified, parameter_handler);
183  std::unique_ptr<FlowSolver::PeriodicEntropyTests<dim, nstate>> flow_solver_case = std::make_unique<FlowSolver::PeriodicEntropyTests<dim,nstate>>(&params_modified);
184 
185  static_cast<void>(flow_solver->run());
186  pcout << "Finished flowsolver " << std::endl;
187 
188  const double final_time_actual = flow_solver->ode_solver->current_time;
189 
190  double L1_error_density=0;
191  double L1_error_pressure=0;
192  calculate_Lp_error_at_final_time_wrt_function(L1_error_density, L1_error_pressure, flow_solver->dg, params_modified,final_time_actual, 1);
193  double L2_error_density =0;
194  double L2_error_pressure=0;
195  calculate_Lp_error_at_final_time_wrt_function(L2_error_density,L2_error_pressure, flow_solver->dg, params_modified,final_time_actual, 2);
196  double Linfty_error_density = 0;
197  double Linfty_error_pressure = 0;
198  calculate_Lp_error_at_final_time_wrt_function(Linfty_error_density,Linfty_error_pressure, flow_solver->dg, params_modified,final_time_actual, -1);
199  pcout << "Computed density errors are: " << std::endl
200  << " L1: " << L1_error_density << std::endl
201  << " L2: " << L2_error_density << std::endl
202  << " Linfty: " << Linfty_error_density << std::endl;
203 
204  const double dt = flow_solver_case->get_constant_time_step(flow_solver->dg);
205  const int n_cells = pow(params_modified.flow_solver_param.number_of_grid_elements_per_dimension, PHILIP_DIM);
206  pcout << " at dt = " << dt << std::endl;
207 
208  // Convergence for density
209  if (this->all_parameters->flux_reconstruction_type == PHiLiP::Parameters::AllParameters::Flux_Reconstruction::user_specified_value){
210  convergence_table_density.add_value("cESFR", params_modified.FR_user_specified_correction_parameter_value );
211  convergence_table_density.set_precision("cESFR", 16);
212  convergence_table_density.set_scientific("cESFR", true);
213  }
214  convergence_table_density.add_value("refinement", refinement);
215  convergence_table_density.add_value("dt", dt );
216  convergence_table_density.set_precision("dt", 16);
217  convergence_table_density.set_scientific("dt", true);
218  convergence_table_density.add_value("n_cells",n_cells);
219  convergence_table_density.add_value("L1_error_density",L1_error_density);
220  convergence_table_density.set_precision("L1_error_density", 16);
221  convergence_table_density.evaluate_convergence_rates("L1_error_density", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
222  convergence_table_density.add_value("L2_error_density",L2_error_density);
223  convergence_table_density.set_precision("L2_error_density", 16);
224  convergence_table_density.evaluate_convergence_rates("L2_error_density", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
225  convergence_table_density.add_value("Linfty_error_density",Linfty_error_density);
226  convergence_table_density.set_precision("Linfty_error_density", 16);
227  convergence_table_density.evaluate_convergence_rates("Linfty_error_density", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
228 
229  if (params_modified.ode_solver_param.ode_solver_type == Parameters::ODESolverParam::ODESolverEnum::rrk_explicit_solver) {
230  const double gamma_agg = final_time_actual / (dt * flow_solver->ode_solver->current_iteration);
231 
232  convergence_table_density.add_value("gamma_agg",gamma_agg-1.0);
233  convergence_table_density.set_precision("gamma_agg", 16);
234  convergence_table_density.evaluate_convergence_rates("gamma_agg", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
235  }
236 
237  // Convergence for pressure
238  if (this->all_parameters->flux_reconstruction_type == PHiLiP::Parameters::AllParameters::Flux_Reconstruction::user_specified_value){
239  convergence_table_pressure.add_value("cESFR", params_modified.FR_user_specified_correction_parameter_value );
240  convergence_table_pressure.set_precision("cESFR", 16);
241  convergence_table_pressure.set_scientific("cESFR", true);
242  }
243  convergence_table_pressure.add_value("refinement", refinement);
244  convergence_table_pressure.add_value("dt", dt );
245  convergence_table_pressure.set_precision("dt", 16);
246  convergence_table_pressure.set_scientific("dt", true);
247  convergence_table_pressure.add_value("n_cells",n_cells);
248  convergence_table_pressure.add_value("L1_error_pressure",L1_error_pressure);
249  convergence_table_pressure.set_precision("L1_error_pressure", 16);
250  convergence_table_pressure.evaluate_convergence_rates("L1_error_pressure", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
251  convergence_table_pressure.add_value("L2_error_pressure",L2_error_pressure);
252  convergence_table_pressure.set_precision("L2_error_pressure", 16);
253  convergence_table_pressure.evaluate_convergence_rates("L2_error_pressure", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
254  convergence_table_pressure.add_value("Linfty_error_pressure",Linfty_error_pressure);
255  convergence_table_pressure.set_precision("Linfty_error_pressure", 16);
256  convergence_table_pressure.evaluate_convergence_rates("Linfty_error_pressure", "n_cells", dealii::ConvergenceTable::reduction_rate_log2, PHILIP_DIM);
257 
258  //Checking convergence order
259  double expected_order = params_modified.flow_solver_param.poly_degree + 1;
260  if(c_value > c_array[nb_c_value])
261  expected_order -= 1;
262  //set tolerance to make test pass for ctest. Note that the grids are very coarse (not in asymptotic range)
263  const double order_tolerance = 1.0;
264  if (refinement > 0) {
265  L2_error_pressure_conv_rate = -log(L2_error_pressure_old/L2_error_pressure)/log(refine_ratio);
266  pcout << "Order for L2 pressure at " << refinement << " is " << L2_error_pressure_conv_rate << std::endl;
267  if (abs(L2_error_pressure_conv_rate - expected_order) > order_tolerance){
268  testfail = 1;
269  pcout << "Expected convergence order for L2 pressure was not reached at refinement " << refinement <<std::endl;
270  }
271  if (refinement < n_calculations-1 && pcout.is_active()){
272  // Print current convergence results for solution monitoring
273  convergence_table_density.write_text(pcout.get_stream());
274  convergence_table_pressure.write_text(pcout.get_stream());
275  }
276  }
277  L2_error_pressure_old = L2_error_pressure;
278  }// refinement loop
279  //Printing and writing convergence tables
280  pcout << std::endl;
281  if (pcout.is_active()){
282  convergence_table_density.write_text(pcout.get_stream());
283  convergence_table_pressure.write_text(pcout.get_stream());
284  }
285 
286  convergence_table_density.write_text(conv_tab_file);
287  convergence_table_pressure.write_text(conv_tab_file);
288 
289  convergence_table_density.clear();
290  convergence_table_pressure.clear();
291  }// c ESFR range loop
292  conv_tab_file.close();
293  return testfail;
294 }
295 #if PHILIP_DIM!=1
297 #endif
298 } // Tests namespace
299 } // 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)
double ESFR_parameter_values_start
For user defined FR parameter tests, value of starting FR param.
const MPI_Comm mpi_communicator
MPI communicator.
Definition: tests.h:39
const int n_calculations
Number of times to solve for convergence summary.
int number_ESFR_parameter_values
For user defined FR parameter tests, number of values to be tested.
double ESFR_parameter_values_end
For user defined FR parameter tests, value of final FR param.
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.
Main parameter class that contains the various other sub-parameter classes.
Flux_Reconstruction flux_reconstruction_type
Store flux reconstruction type.
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, double cvalue, int nb_c_value) const
Reinitialize parameters while refining the timestep. Necessary because all_parameters is constant...
double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
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.
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.