[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
general_refinement_study.h
1 #ifndef __TIME_REFINEMENT_STUDY__
2 #define __TIME_REFINEMENT_STUDY__
3 
4 #include <deal.II/base/convergence_table.h>
5 
6 #include "flow_solver/flow_solver.h"
7 #include "dg/dg_base.hpp"
8 #include "tests.h"
9 
10 namespace PHiLiP {
11 namespace Tests {
12 
14 template <int dim, int nstate>
16 {
17 public:
19  enum RefinementType { timestep, cell_length }; // in the future, can also add refinement by poly_degree here
20 
23  const Parameters::AllParameters *const parameters_input,
24  const dealii::ParameterHandler &parameter_handler_input,
25  const RefinementType refinement_type_input
26  );
27 
29  const dealii::ParameterHandler &parameter_handler;
30 
32  int run_test () const override;
33 protected:
36 
38  const int n_calculations;
39 
41  const double refine_ratio;
42 
43 
45  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;
46 
51  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;
52 
54  virtual std::tuple<double,int> process_and_write_conv_tables(std::shared_ptr<FlowSolver::FlowSolver<dim,nstate>> flow_solver,
55  const Parameters::AllParameters params,
56  double L2_error_old,
57  std::shared_ptr<dealii::ConvergenceTable> convergence_table,
58  int refinement,
59  const double expected_order) const;
61  Parameters::AllParameters reinit_params_and_refine(const Parameters::AllParameters *parameters_in, int refinement, const RefinementType how) const;
62 
63 };
64 
65 } // End of Tests namespace
66 } // End of PHiLiP namespace
67 
68 #endif
GeneralRefinementStudy(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input, const RefinementType refinement_type_input)
Constructor.
Advection time refinement study.
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
Selects which flow case to simulate.
Definition: flow_solver.h:64
Parameters::AllParameters reinit_params_and_refine(const Parameters::AllParameters *parameters_in, int refinement, const RefinementType how) const
Reinitialize parameters while refining the timestep. Necessary because all_parameters is constant...
Files for the baseline physics.
Definition: ADTypes.hpp:10
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.
const RefinementType refinement_type
Type of refinement to run.
const int n_calculations
Number of times to solve for convergence summary.
Main parameter class that contains the various other sub-parameter classes.
const double refine_ratio
Ratio to refine by.
const dealii::ParameterHandler & parameter_handler
Parameter handler for storing the .prm file being ran.
virtual 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
Calculate the L2 error and return local testfail for the converged flowsolver.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
Base class of all the tests.
Definition: tests.h:17