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" 12 template <
int dim,
int nstate>
15 const dealii::ParameterHandler ¶meter_handler_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)
24 template <
int dim,
int nstate>
40 template <
int dim,
int nstate>
44 std::shared_ptr<ExactSolutionFunction<dim,nstate,double>> exact_solution_function;
46 this->
pcout <<
"End time: " << final_time << std::endl;
48 int overintegrate = 10;
51 dealii::Vector<double> difference_per_cell(dg->solution.size());
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;
58 pcout <<
"Norm not defined. Aborting... " << std::endl;
62 dealii::VectorTools::integrate_difference(dg->dof_handler,
64 *exact_solution_function,
66 dealii::QGauss<dim>(poly_degree + overintegrate),
68 Lp_error = dealii::VectorTools::compute_global_error(*dg->triangulation,
74 template <
int dim,
int nstate>
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;
89 dealii::ConvergenceTable convergence_table;
90 double L2_error_old = 0;
91 double L2_error_conv_rate=0;
95 pcout <<
"\n\n---------------------------------------------" << std::endl;
96 pcout <<
"Refinement number " << refinement <<
" of " << n_time_calculations - 1 << std::endl;
97 pcout <<
"---------------------------------------------" << std::endl;
101 static_cast<void>(flow_solver->run());
103 pcout <<
"Finished flowsolver " << std::endl;
105 const double final_time_actual = flow_solver->ode_solver->current_time;
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;
117 const int n_timesteps= flow_solver->ode_solver->current_iteration;
118 pcout <<
" at dt = " << dt << std::endl;
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){
128 pcout <<
"Expected L2 error for RK3(2)5F[3S*+] using an atol = rtol = 1e-4 was not reached " << refinement <<std::endl;
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);
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);
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){
164 pcout <<
"Expected convergence order was not reached at refinement " << refinement <<std::endl;
168 if (refinement < n_time_calculations-1 &&
pcout.is_active()) convergence_table.write_text(
pcout.get_stream());
170 L2_error_old = L2_error;
175 if (
pcout.is_active()) convergence_table.write_text(
pcout.get_stream());
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();
double final_time
Final solution time.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
Files for the baseline physics.
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.
ODESolverEnum
Types of ODE solver.
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 ¶meter_handler_input)
Factory to return the correct flow solver given input file.
const Parameters::AllParameters *const all_parameters
Pointer to all parameters.
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 ¶meter_handler_input)
Constructor.
ODESolverEnum ode_solver_type
ODE solver type.
dealii::ConditionalOStream pcout
ConditionalOStream.
DGBase is independent of the number of state variables.
unsigned int poly_degree
Initial solution polynomial degree.
const double refine_ratio
Ratio to refine by.
Base class of all the tests.
int run_test() const override
Run test.