[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
grid_refinement_study.h
1 #ifndef __GRID_REFINEMENT_STUDY_H__
2 #define __GRID_REFINEMENT_STUDY_H__
3 
4 #include "dg/dg_base.hpp"
5 #include "grid_refinement/gnu_out.h"
6 #include "parameters/all_parameters.h"
7 #include "physics/model.h"
8 #include "physics/physics.h"
9 #include "tests.h"
10 
11 namespace PHiLiP {
12 
13 namespace Tests {
14 
16 #if PHILIP_DIM==1 // dealii::parallel::distributed::Triangulation<dim> does not work for 1D
17 template <int dim, int nstate, typename MeshType = dealii::Triangulation<dim>>
18 #else
19 template <int dim, int nstate, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
20 #endif
22 {
23 public:
25  GridRefinementStudy() = delete;
27 
29  explicit GridRefinementStudy(
30  const Parameters::AllParameters *const parameters_input);
31 
32  int run_test() const;
33 
35  void get_grid(
36  const std::shared_ptr<MeshType>& grid,
37  const Parameters::GridRefinementStudyParam& grs_param) const;
38 
41  const std::shared_ptr<Physics::PhysicsBase<dim,nstate,double>>& physics_double,
42  const std::shared_ptr<Physics::PhysicsBase<dim,nstate,Sacado::Fad::DFad<double>>>& physics_adtype,
43  const Parameters::AllParameters& param,
44  const Parameters::GridRefinementStudyParam& grs_param) const;
45 };
46 
47 // constructs the mesh for the test
48 template <typename MeshType>
50 {
51 public:
52  static std::shared_ptr<MeshType>
53  create_MeshType(const MPI_Comm mpi_communicator);
54 };
55 
56 // function to perform the formatted output to gnuplot (of the solution error)
57 void output_gnufig_solution(
59 
60 // function to perform the formatted output to gnuplot (of the funcitonal error)
61 void output_gnufig_functional(
63 
64 } // Tests namespace
65 
66 } // PHiLiP namespace
67 
68 #endif // __GRID_REFINEMENT_STUDY_H__
const MPI_Comm mpi_communicator
MPI communicator.
Definition: tests.h:39
double approximate_exact_functional(const std::shared_ptr< Physics::PhysicsBase< dim, nstate, double >> &physics_double, const std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< double >>> &physics_adtype, const Parameters::AllParameters &param, const Parameters::GridRefinementStudyParam &grs_param) const
Approximates the exact functional using a uniformly refined grid.
Gnuplot utility class.
Definition: gnu_out.h:16
Files for the baseline physics.
Definition: ADTypes.hpp:10
Main parameter class that contains the various other sub-parameter classes.
Performs grid convergence for various polynomial degrees.
GridRefinementStudy()=delete
Constructor. Deleted the default constructor since it should not be used.
Parameters related to collection of grid refinement runs.
int run_test() const
Basically the main and only function of this class.
Base class of all the tests.
Definition: tests.h:17
void get_grid(const std::shared_ptr< MeshType > &grid, const Parameters::GridRefinementStudyParam &grs_param) const
gets the grid from the enum and reads file if neccesary