[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
grid_study.h
1 #ifndef __GRID_STUDY_H__
2 #define __GRID_STUDY_H__
3 
4 #include <deal.II/base/convergence_table.h>
5 
6 #include "dg/dg_base.hpp"
7 #include "parameters/all_parameters.h"
8 #include "physics/physics.h"
9 #include "tests.h"
10 
11 namespace PHiLiP {
12 namespace Tests {
13 
15 template <int dim, int nstate>
16 class GridStudy: public TestsBase
17 {
18 public:
20  GridStudy () = delete;
22 
24  explicit GridStudy(const Parameters::AllParameters *const parameters_input);
25 
27 
37  int run_test () const;
38 
39 protected:
41  void print_mesh_info(const dealii::Triangulation<dim> &triangulation,
42  const std::string &filename) const;
44 
46  static dealii::Point<dim> warp (const dealii::Point<dim> &p);
47 
49 
54 
57 
59  std::string get_convergence_tables_baseline_filename(const Parameters::AllParameters *const parameters_input) const;
60 
63  const std::string error_filename_baseline,
64  const dealii::ConvergenceTable convergence_table,
65  const unsigned int poly_degree) const;
66 };
67 
68 
69 // /// Manufactured grid convergence
70 // /** Currently the main function as all my test cases simply
71 // * check for optimal convergence of the solution
72 // */
73 // template<int dim>
74 // int manufactured_grid_convergence (Parameters::AllParameters &parameters);
75 
76 } // Tests namespace
77 } // PHiLiP namespace
78 #endif
std::string get_convergence_tables_baseline_filename(const Parameters::AllParameters *const parameters_input) const
Returns the baseline filename for outputted convergence tables.
Definition: grid_study.cpp:118
double integrate_solution_over_domain(DGBase< dim, double > &dg) const
L2-Integral of the solution over the entire domain.
Definition: grid_study.cpp:65
Performs grid convergence for various polynomial degrees.
Definition: grid_study.h:16
int run_test() const
Manufactured grid convergence.
Definition: grid_study.cpp:153
Files for the baseline physics.
Definition: ADTypes.hpp:10
Main parameter class that contains the various other sub-parameter classes.
void write_convergence_table_to_output_file(const std::string error_filename_baseline, const dealii::ConvergenceTable convergence_table, const unsigned int poly_degree) const
Writes the convergence table output file.
Definition: grid_study.cpp:138
GridStudy()=delete
Constructor. Deleted the default constructor since it should not be used.
void print_mesh_info(const dealii::Triangulation< dim > &triangulation, const std::string &filename) const
Prints our mesh info and generates eps file if 2D grid.
Definition: grid_study.cpp:651
void initialize_perturbed_solution(DGBase< dim, double > &dg, const Physics::PhysicsBase< dim, nstate, double > &physics) const
Initialize the solution with the exact solution.
Definition: grid_study.cpp:49
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
Base class of all the tests.
Definition: tests.h:17
static dealii::Point< dim > warp(const dealii::Point< dim > &p)
Warps mesh into sinusoidal.
Definition: grid_study.cpp:640