[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_fixed_fraction.h
1 #ifndef __GRID_REFINEMENT_FIXED_FRACTION_H__
2 #define __GRID_REFINEMENT_FIXED_FRACTION_H__
3 
4 #include "grid_refinement/grid_refinement.h"
5 
6 namespace PHiLiP {
7 
8 namespace GridRefinement {
9 
11 
22 #if PHILIP_DIM==1
23 template <int dim, int nstate, typename real, typename MeshType = dealii::Triangulation<dim>>
24 #else
25 template <int dim, int nstate, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
26 #endif
27 class GridRefinement_FixedFraction : public GridRefinementBase<dim,nstate,real,MeshType>
28 {
29 public:
31 
33 
45  void refine_grid() override;
46 
47 protected:
48 
49  // specified refinement functions for different cases
50 
52 
60  void refine_grid_h();
61 
63 
67  void refine_grid_p();
68 
70 
76  void refine_grid_hp();
77 
79 
83  void refine_boundary_h();
84 
86 
93  void smoothness_indicator();
94 
96 
101  void anisotropic_h();
102 
104 
113 
115 
123 
125 
130  void error_indicator();
131 
132  // error distribution for each indicator type
133 
135 
145  void error_indicator_error();
146 
148 
163 
165 
177 
179 
198 
200 
202  std::vector< std::pair<dealii::Vector<real>, std::string> > output_results_vtk_method() override;
203 
204 protected:
205  dealii::Vector<real> indicator;
206  dealii::Vector<real> smoothness;
207 };
208 
209 } // namespace GridRefinement
210 
211 } // namespace PHiLiP
212 
213 #endif // __GRID_REFINEMENT_H__
void refine_grid() override
Perform call to the grid refinement object of choice.
void error_indicator_error()
Compute error indicator based on Lq norm relative to exact manufactured solution. ...
void error_indicator()
Performs call to proper error indcator function based on type of error_indicator parameter.
void anisotropic_h_jump_based()
Sets anisotropic refinement flags for cells based on discontinuity between neighbors.
Files for the baseline physics.
Definition: ADTypes.hpp:10
void anisotropic_h_reconstruction_based()
Sets anisotropic refinement flags for cells based on directional derivatives reconstructed along the ...
void refine_grid_hp()
(NOT IMPLEMENTED) Based on error and smoothness indicator, perform fixed-fraction flagging and decisi...
void error_indicator_adjoint()
Compute error indicator for goal-oriented approach using a dual-weighted residual.
void refine_grid_p()
Based on error indicator, performs fixed-fraction flagging for polynomial order enrichment.
void smoothness_indicator()
(NOT IMPLEMENTED) Computes smoothness indicator for -refinement decisions
void anisotropic_h()
Performs anisotropic adaptation of the mesh.
std::vector< std::pair< dealii::Vector< real >, std::string > > output_results_vtk_method() override
Output refinement method dependent results.
void error_indicator_hessian()
Compute error indicator based on reconstructed directional derivatives.
void error_indicator_residual()
(NOT IMPLEMENTED) Compute error indicator based on fine grid residual distribution ...
void refine_boundary_h()
Flags the domain boundaries for refinement.
void refine_grid_h()
Based on error indicator, performs fixed-fraction flagging for element splitting. ...