[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
time_refinement_study.cpp
1 #include "time_refinement_study.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 "cmath"
7 //#include "ode_solver/runge_kutta_ode_solver.h"
8 
9 namespace PHiLiP {
10 namespace Tests {
11 
12 template <int dim, int nstate>
14  const PHiLiP::Parameters::AllParameters *const parameters_input,
15  const dealii::ParameterHandler &parameter_handler_input)
16  : TestsBase::TestsBase(parameters_input),
17  parameter_handler(parameter_handler_input),
18  n_time_calculations(parameters_input->time_refinement_study_param.number_of_times_to_solve),
19  refine_ratio(parameters_input->time_refinement_study_param.refinement_ratio)
20 {}
21 
22 
23 
24 template <int dim, int nstate>
26 {
28 
29  parameters.ode_solver_param.initial_time_step *= pow(refine_ratio,refinement);
30 
31  //For RRK, do not end at exact time because of how relaxation parameter convergence is calculated
32  using ODESolverEnum = Parameters::ODESolverParam::ODESolverEnum;
33  if (parameters.ode_solver_param.ode_solver_type == ODESolverEnum::rrk_explicit_solver){
35  }
36 
37  return parameters;
38 }
39 
40 template <int dim, int nstate>
41 double TimeRefinementStudy<dim,nstate>::calculate_Lp_error_at_final_time_wrt_function(std::shared_ptr<DGBase<dim,double>> dg, const Parameters::AllParameters parameters, double final_time, 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  this->pcout << "End time: " << final_time << std::endl;
47  double Lp_error=0;
48  int overintegrate = 10;
49 
50  int poly_degree = parameters.grid_refinement_study_param.poly_degree;
51  dealii::Vector<double> difference_per_cell(dg->solution.size());
52 
53  dealii::VectorTools::NormType norm_type;
54  if (norm_p == 1) norm_type = dealii::VectorTools::L1_norm;
55  else if (norm_p == 2) norm_type = dealii::VectorTools::L2_norm;
56  else if (norm_p<0) norm_type = dealii::VectorTools::Linfty_norm;
57  else {
58  pcout << "Norm not defined. Aborting... " << std::endl;
59  std::abort();
60  }
61 
62  dealii::VectorTools::integrate_difference(dg->dof_handler,
63  dg->solution,
64  *exact_solution_function,
65  difference_per_cell,
66  dealii::QGauss<dim>(poly_degree + overintegrate), //overintegrating
67  norm_type);
68  Lp_error = dealii::VectorTools::compute_global_error(*dg->triangulation,
69  difference_per_cell,
70  norm_type);
71  return Lp_error;
72 }
73 
74 template <int dim, int nstate>
76 {
77 
78  const double final_time = this->all_parameters->flow_solver_param.final_time;
79  const double initial_time_step = this->all_parameters->ode_solver_param.initial_time_step;
80  const int n_steps = floor(final_time/initial_time_step);
81  if (n_steps * initial_time_step != final_time){
82  pcout << "WARNING: final_time is not evenly divisible by initial_time_step!" << std::endl
83  << "Remainder is " << fmod(final_time, initial_time_step)
84  << ". Consider modifying parameters." << std::endl;
85  }
86 
87  int testfail = 0;
88 
89  dealii::ConvergenceTable convergence_table;
90  double L2_error_old = 0;
91  double L2_error_conv_rate=0;
92 
93  for (int refinement = 0; refinement < n_time_calculations; ++refinement){
94 
95  pcout << "\n\n---------------------------------------------" << std::endl;
96  pcout << "Refinement number " << refinement << " of " << n_time_calculations - 1 << std::endl;
97  pcout << "---------------------------------------------" << std::endl;
98 
100  std::unique_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver = FlowSolver::FlowSolverFactory<dim,nstate>::select_flow_case(&params, parameter_handler);
101  static_cast<void>(flow_solver->run());
102 
103  pcout << "Finished flowsolver " << std::endl;
104 
105  const double final_time_actual = flow_solver->ode_solver->current_time;
106 
107  //check Lp error
108  const double L1_error = calculate_Lp_error_at_final_time_wrt_function(flow_solver->dg, params,final_time_actual, 1);
109  const double L2_error = calculate_Lp_error_at_final_time_wrt_function(flow_solver->dg, params,final_time_actual, 2);
110  const double Linfty_error = calculate_Lp_error_at_final_time_wrt_function(flow_solver->dg, params,final_time_actual, -1);
111  pcout << "Computed errors are: " << std::endl
112  << " L1: " << L1_error << std::endl
113  << " L2: " << L2_error << std::endl
114  << " Linfty: " << Linfty_error << std::endl;
115 
116  const double dt = params.ode_solver_param.initial_time_step;
117  const int n_timesteps= flow_solver->ode_solver->current_iteration;
118  pcout << " at dt = " << dt << std::endl;
119 
120 
122  && this->all_parameters->ode_solver_param.atol == 1e-4 && this->all_parameters->ode_solver_param.rtol == 1e-4
123  && this->all_parameters->time_refinement_study_param.number_of_times_to_solve == 1){
124  double L2_error_expected = 2.14808703658e-5;
125  pcout << " Expected L2 error is: " << L2_error_expected << std::endl;
126  if (L2_error > L2_error_expected + 1e-9 || L2_error < L2_error_expected - 1e-9){
127  testfail = 1;
128  pcout << "Expected L2 error for RK3(2)5F[3S*+] using an atol = rtol = 1e-4 was not reached " << refinement <<std::endl;
129  }
130  }
131 
132  convergence_table.add_value("refinement", refinement);
133  convergence_table.add_value("dt", dt );
134  convergence_table.set_precision("dt", 16);
135  convergence_table.set_scientific("dt", true);
136  convergence_table.add_value("L1_error",L1_error);
137  convergence_table.set_precision("L1_error", 16);
138  convergence_table.evaluate_convergence_rates("L1_error", "dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
139  convergence_table.add_value("L2_error",L2_error);
140  convergence_table.set_precision("L2_error", 16);
141  convergence_table.evaluate_convergence_rates("L2_error", "dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
142  convergence_table.add_value("Linfty_error",Linfty_error);
143  convergence_table.set_precision("Linfty_error", 16);
144  convergence_table.evaluate_convergence_rates("Linfty_error", "dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
145 
146  using ODESolverEnum = Parameters::ODESolverParam::ODESolverEnum;
147  if(params.ode_solver_param.ode_solver_type == ODESolverEnum::rrk_explicit_solver){
148  //for burgers, this is the average gamma over the runtime
149  const double gamma_aggregate_m1 = final_time_actual / (n_timesteps * dt)-1;
150  convergence_table.add_value("gamma_aggregate_m1", gamma_aggregate_m1);
151  convergence_table.set_precision("gamma_aggregate_m1", 16);
152  convergence_table.set_scientific("gamma_aggregate_m1", true);
153  convergence_table.evaluate_convergence_rates("gamma_aggregate_m1", "dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
154  }
155 
156  //Checking convergence order
157  const double expected_order = params.ode_solver_param.rk_order;
158  const double order_tolerance = 0.1;
159  if (refinement > 0) {
160  L2_error_conv_rate = -log(L2_error_old/L2_error)/log(refine_ratio);
161  pcout << "L2 order at " << refinement << " is " << L2_error_conv_rate << std::endl;
162  if (abs(L2_error_conv_rate - expected_order) > order_tolerance){
163  testfail = 1;
164  pcout << "Expected convergence order was not reached at refinement " << refinement <<std::endl;
165  }
166 
167  // output current time refinement results to console
168  if (refinement < n_time_calculations-1 && pcout.is_active()) convergence_table.write_text(pcout.get_stream());
169  }
170  L2_error_old = L2_error;
171  }
172 
173  //Printing and writing convergence table
174  pcout << std::endl;
175  if (pcout.is_active()) convergence_table.write_text(pcout.get_stream());
176 
177  std::ofstream conv_tab_file;
178  const std::string fname = "temporal_convergence_table.txt";
179  conv_tab_file.open(fname);
180  convergence_table.write_text(conv_tab_file);
181  conv_tab_file.close();
182 
183  return testfail;
184 }
185 
186 #if PHILIP_DIM==1
188 #endif
189 } // Tests namespace
190 } // PHiLiP namespace
double final_time
Final solution time.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
Files for the baseline physics.
Definition: ADTypes.hpp:10
GridRefinementStudyParam grid_refinement_study_param
Contains the parameters for grid refinement study.
Parameters::AllParameters reinit_params_and_refine_timestep(int refinement) const
Reinitialize parameters while refining the timestep. Necessary because all_parameters is constant...
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.
int rk_order
Order of the RK method; assigned based on runge_kutta_method.
ODESolverParam ode_solver_param
Contains parameters for ODE solver.
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
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
Fourth-order, three register low-storage Runge-Kutta method.
const int n_time_calculations
Number of times to solve for convergence summary.
bool end_exactly_at_final_time
Flag to adjust the last timestep such that the simulation ends exactly at final_time.
Advection time refinement study.
RKMethodEnum runge_kutta_method
Runge-kutta method.
double calculate_Lp_error_at_final_time_wrt_function(std::shared_ptr< DGBase< dim, double >> dg, const Parameters::AllParameters parameters, double final_time, int norm_p) const
TimeRefinementStudy(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Constructor.
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
unsigned int poly_degree
Initial solution polynomial degree.
const double refine_ratio
Ratio to refine by.
Base class of all the tests.
Definition: tests.h:17
int run_test() const override
Run test.