1 #include "time_refinement_study_reference.h" 2 #include "flow_solver/flow_solver_factory.h" 3 #include "flow_solver/flow_solver_cases/periodic_1D_unsteady.h" 9 template <
int dim,
int nstate>
12 const dealii::ParameterHandler ¶meter_handler_input)
15 , reference_solution(this->calculate_reference_solution(parameters_input->flow_solver_param.final_time))
18 template <
int dim,
int nstate>
23 const double dt = final_time/number_of_timesteps;
34 this->
pcout <<
"Using timestep size dt = " << dt <<
" for reference solution." << std::endl;
39 template <
int dim,
int nstate>
41 const double final_time)
const 46 static_cast<void>(flow_solver_reference->run());
48 this->
pcout <<
" Actual final time: " << flow_solver_reference->ode_solver->current_time << std::endl;
50 return flow_solver_reference->dg->solution;
53 template <
int dim,
int nstate>
57 const double final_time_actual
62 if (abs(final_time_target-final_time_actual)<1E-13){
64 this->
pcout <<
"Comparing to reference solution at target final_time = " << final_time_target <<
" ..." << std::endl;
67 dealii::LinearAlgebra::distributed::Vector<double> cellwise_difference(this->reference_solution);
68 cellwise_difference.add(-1.0, dg->solution);
69 const double L2_error = cellwise_difference.l2_norm();
75 this->
pcout <<
" -------------------------------------------------------" << std::endl;
76 this->
pcout <<
" Calculating reference solution at actual final_time = " << final_time_actual <<
" ..."<<std::endl;
77 this->
pcout <<
" -------------------------------------------------------" << std::endl;
78 const dealii::LinearAlgebra::distributed::Vector<double> reference_solution_actual = calculate_reference_solution(final_time_actual);
80 dealii::LinearAlgebra::distributed::Vector<double> cellwise_difference(reference_solution_actual);
81 cellwise_difference.add(-1.0, dg->solution);
82 const double L2_error = cellwise_difference.l2_norm();
88 template <
int dim,
int nstate>
92 std::shared_ptr<dealii::ConvergenceTable> convergence_table,
94 const double expected_order)
const 99 std::unique_ptr<FlowSolver::Periodic1DUnsteady<dim, nstate>> flow_solver_case = std::make_unique<FlowSolver::Periodic1DUnsteady<dim,nstate>>(this->
all_parameters);
101 this->
pcout <<
"ERROR: the TimeRefinementStudyReference class is currently only valid " << std::endl
102 <<
"for flow_case_type = periodic_1D_unsteady. Aborting..." << std::endl;
105 const double energy_initial = flow_solver_case->compute_energy(flow_solver->dg);
106 const double final_time_actual = flow_solver->ode_solver->current_time;
107 this->
pcout <<
" Actual final time: " << final_time_actual << std::endl;
109 const int n_timesteps= flow_solver->ode_solver->current_iteration;
116 this->
pcout <<
"Computed error is " << L2_error << std::endl;
119 convergence_table->add_value(
"refinement", refinement);
120 convergence_table->add_value(
"dt", dt );
121 convergence_table->set_precision(
"dt", 16);
122 convergence_table->set_scientific(
"dt",
true);
123 convergence_table->add_value(
"final_time", final_time_actual );
124 convergence_table->set_precision(
"final_time", 16);
125 convergence_table->add_value(
"n_timesteps", n_timesteps);
126 convergence_table->add_value(
"L2_error",L2_error);
127 convergence_table->set_precision(
"L2_error", 16);
128 convergence_table->evaluate_convergence_rates(
"L2_error",
"dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
130 const double energy_end = flow_solver_case->compute_energy(flow_solver->dg);
131 const double energy_change = energy_initial - energy_end;
132 convergence_table->add_value(
"energy_change", energy_change);
133 convergence_table->set_precision(
"energy_change", 16);
134 convergence_table->evaluate_convergence_rates(
"energy_change",
"dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
140 const double gamma_aggregate_m1 = final_time_actual / (n_timesteps * dt)-1;
141 convergence_table->add_value(
"gamma_aggregate_m1", gamma_aggregate_m1);
142 convergence_table->set_precision(
"gamma_aggregate_m1", 16);
143 convergence_table->set_scientific(
"gamma_aggregate_m1",
true);
144 convergence_table->evaluate_convergence_rates(
"gamma_aggregate_m1",
"dt", dealii::ConvergenceTable::reduction_rate_log2, 1);
148 const double order_tolerance = 0.1;
149 if (refinement > 0) {
150 double L2_error_conv_rate = -log(L2_error_old/L2_error)/log(this->
refine_ratio);
151 this->
pcout <<
"Order at " << refinement <<
" is " << L2_error_conv_rate << std::endl;
152 if (abs(L2_error_conv_rate - expected_order) > order_tolerance){
154 this->
pcout <<
"Expected convergence order was not reached at refinement " << refinement <<std::endl;
158 return std::make_tuple(L2_error, testfail);
162 template <
int dim,
int nstate>
double final_time
Final solution time.
RefinementType
Type of refinement to run.
Advection time refinement study.
FlowCaseType flow_case_type
Selected FlowCaseType from the input file.
TimeRefinementStudyParam time_refinement_study_param
Contains the parameters for time refinement study.
FlowSolverParam flow_solver_param
Contains the parameters for simulation cases (flow solver test)
Selects which flow case to simulate.
Files for the baseline physics.
int run_refinement_study_and_write_result(const Parameters::AllParameters *parameters_in, const double expected_order, const bool check_only_last_refinement=false, const bool append_to_file=false) const
Run the refinements.
double calculate_L2_error_at_final_time_wrt_reference(std::shared_ptr< DGBase< dim, double >> dg, const Parameters::AllParameters parameters, double final_time_actual) const
Calculate L2 error at the final time in the passed parameters.
ODESolverEnum
Types of ODE solver.
int run_test() const override
Run test.
int number_of_timesteps_for_reference_solution
For time refinement study with reference solution, number of steps for reference solution.
Parameters::AllParameters reinit_params_for_reference_solution(int number_of_timesteps, double final_time) const
Reinitialize parameters and set initial_timestep according to reference solution and passed final tim...
Main parameter class that contains the various other sub-parameter classes.
const double refine_ratio
Ratio to refine by.
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.
std::tuple< double, int > process_and_write_conv_tables(std::shared_ptr< FlowSolver::FlowSolver< dim, nstate >> flow_solver, const Parameters::AllParameters params, double L2_error_old, std::shared_ptr< dealii::ConvergenceTable > convergence_table, int refinement, const double expected_order) const override
Calculate the L2 error and return local testfail for the converged flowsolver.
bool end_exactly_at_final_time
Flag to adjust the last timestep such that the simulation ends exactly at final_time.
TimeRefinementStudyReference(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler ¶meter_handler_input)
Constructor.
ODESolverEnum ode_solver_type
ODE solver type.
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
dealii::ConditionalOStream pcout
ConditionalOStream.
Time refinement study which compares to a reference solution.
DGBase is independent of the number of state variables.