[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
mesh_adaptation.h
1 #ifndef __MESHADAPTATION_H__
2 #define __MESHADAPTATION_H__
3 
4 #include <deal.II/distributed/grid_refinement.h>
5 #include <deal.II/distributed/shared_tria.h>
6 #include <deal.II/distributed/tria.h>
7 #include <deal.II/grid/grid_refinement.h>
8 #include <deal.II/grid/tria.h>
9 
10 #include "dg/dg_base.hpp"
11 #include "mesh_error_estimate.h"
12 #include "mesh_error_factory.h"
13 #include "parameters/all_parameters.h"
14 
15 namespace PHiLiP {
16 
17 
18 #if PHILIP_DIM==1
19 template <int dim, typename real, typename MeshType = dealii::Triangulation<dim>>
20 #else
21 template <int dim, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
22 #endif
23 
28 {
29 public:
30 
32  MeshAdaptation(std::shared_ptr< DGBase<dim, real, MeshType> > dg_input, const Parameters::MeshAdaptationParam *const mesh_adaptation_param_input);
33 
35  ~MeshAdaptation() = default;
36 
38  std::unique_ptr<MeshErrorEstimateBase<dim, real, MeshType>> mesh_error;
39 
41  std::shared_ptr<DGBase<dim,real,MeshType>> dg;
42 
44  void adapt_mesh();
45 
48 
51 
52 protected:
53 
56 
59 
61  dealii::Vector<real> cellwise_errors;
62 
64  dealii::ConditionalOStream pcout;
65 
66 };
67 
68 } // namespace PHiLiP
69 
70 #endif
void adapt_mesh()
Function to adapt the mesh based on input parameters.
const Parameters::MeshAdaptationParam *const mesh_adaptation_param
Holds parameters of mesh adaptation.
MeshAdaptation(std::shared_ptr< DGBase< dim, real, MeshType > > dg_input, const Parameters::MeshAdaptationParam *const mesh_adaptation_param_input)
Constructor to initialize the class with a pointer to DG.
Files for the baseline physics.
Definition: ADTypes.hpp:10
~MeshAdaptation()=default
Destructor.
dealii::Vector< real > cellwise_errors
Stores errors in each cell.
std::unique_ptr< MeshErrorEstimateBase< dim, real, MeshType > > mesh_error
Pointer to the error estimator class.
std::shared_ptr< DGBase< dim, real, MeshType > > dg
Pointer to DGBase.
dealii::ConditionalOStream pcout
Parallel std::cout.
void smoothness_sensor_based_hp_refinement()
Decide whether to perform h or p refinement based on a smoothness indicator.
int current_mesh_adaptation_cycle
Stores the current adaptation cycle.
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
void fixed_fraction_isotropic_refinement_and_coarsening()
Performs fixed fraction refinement based on refinement and coarsening fractions.