[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_continuous.h
1 #ifndef __GRID_REFINEMENT_CONTINUOUS_H__
2 #define __GRID_REFINEMENT_CONTINUOUS_H__
3 
4 #include "grid_refinement/grid_refinement.h"
5 
6 namespace PHiLiP {
7 
8 namespace GridRefinement {
9 
11 
32 #if PHILIP_DIM==1
33 template <int dim, int nstate, typename real, typename MeshType = dealii::Triangulation<dim>>
34 #else
35 template <int dim, int nstate, typename real, typename MeshType = dealii::parallel::distributed::Triangulation<dim>>
36 #endif
37 class GridRefinement_Continuous : public GridRefinementBase<dim,nstate,real,MeshType>
38 {
39 public:
40  // overriding the other constructors to call delegated constructor for this class
41 
45  std::shared_ptr< PHiLiP::Adjoint<dim, nstate, real, MeshType> > adj_input,
46  std::shared_ptr< PHiLiP::Physics::PhysicsBase<dim,nstate,real> > physics_input);
47 
51  std::shared_ptr< PHiLiP::DGBase<dim, real, MeshType> > dg_input,
52  std::shared_ptr< PHiLiP::Physics::PhysicsBase<dim,nstate,real> > physics_input,
53  std::shared_ptr< PHiLiP::Functional<dim, nstate, real, MeshType> > functional_input);
54 
58  std::shared_ptr< PHiLiP::DGBase<dim, real, MeshType> > dg_input,
59  std::shared_ptr< PHiLiP::Physics::PhysicsBase<dim,nstate,real> > physics_input);
60 
64  // PHiLiP::Parameters::AllParameters const *const param_input,
65  std::shared_ptr< PHiLiP::DGBase<dim, real, MeshType> > dg_input);
66 
67 protected:
71  std::shared_ptr< PHiLiP::Adjoint<dim, nstate, real, MeshType> > adj_input,
72  std::shared_ptr< PHiLiP::Functional<dim, nstate, real, MeshType> > functional_input,
73  std::shared_ptr< PHiLiP::DGBase<dim, real, MeshType> > dg_input,
74  std::shared_ptr< PHiLiP::Physics::PhysicsBase<dim,nstate,real> > physics_input);
75 
77 
79 
90  void refine_grid() override;
91 
92 protected:
93 
94  // specified refinement functions for different cases
95 
97 
106  void refine_grid_h();
107 
109 
114  void refine_grid_p();
115 
117 
121  void refine_grid_hp();
122 
124 
128  std::vector< std::pair<dealii::Vector<real>, std::string> > output_results_vtk_method() override;
129 
130  // getting the size or tensor fields based on indicator
131 
133 
139  void field();
140 
141  // based on exact error function from manufactured solution
142 
144 
157  void field_h_error();
159  void field_p_error();
161  void field_hp_error();
162 
163  // based on hessian or feature-based error
164 
166 
185  void field_h_hessian();
187  void field_p_hessian();
189  void field_hp_hessian();
190 
191  // based on high-order solution residual
192 
194  void field_h_residual();
196  void field_p_residual();
198  void field_hp_residual();
199 
200  // based on adjoint solution and dual-weighted residual
201 
203 
221  void field_h_adjoint();
223  void field_p_adjoint();
225  void field_hp_adjoint();
226 
227  // performing output to appropriate mesh generator
228 
230 
242  void refine_grid_gmsh();
243 
245 
272  void refine_grid_msh();
273 
274  // scheduling of complexity growth
275 
277 
288  real current_complexity();
289 
291 
298  void target_complexity();
299 
302  std::vector<real> complexity_vector;
303 
305 
310  void get_current_field_h();
311 
313 
316  void get_current_field_p();
317 
318  std::unique_ptr<Field<dim,real>> h_field;
319  dealii::Vector<real> p_field;
320 };
321 
322 } // namespace GridRefinement
323 
324 } // namespace PHiLiP
325 
326 #endif // __GRID_REFINEMENT_CONTINUOUS_H__
void field_h_error()
Generate mesh distribution based on manufactured solution.
void target_complexity()
Updates the complexity target based on the current refinement iteration.
std::unique_ptr< Field< dim, real > > h_field
Continuous representation of the mesh size and anisotropy distribution.
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.
Definition: physics.h:34
Adjoint class.
Definition: adjoint.h:39
void refine_grid_msh()
Generates an output .msh file with matrix information about the target frame field.
Files for the baseline physics.
Definition: ADTypes.hpp:10
void refine_grid_p()
Updates the global polynomial distribution based on the target -field.
GridRefinement_Continuous(PHiLiP::Parameters::GridRefinementParam gr_param_input, std::shared_ptr< PHiLiP::Adjoint< dim, nstate, real, MeshType > > adj_input, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, real > > physics_input)
Constructor. Stores the adjoint object, physics and parameters.
void field_hp_error()
(NOT IMPLEMENTED) Generate mesh and polynomial distribution based on manufactured solution ...
void field_hp_residual()
(NOT IMPLEMENTED) Generate mesh and polynomial distribution based on fine-grid residual distribution ...
void field_hp_hessian()
(NOT IMPLEMENTED) Generate mesh and polynomial distribution based on directional derivatives ...
void field()
Updates the continuous mesh and polynomial distribution.
void get_current_field_p()
Evaluates the polynomial distribution for the current mesh.
void refine_grid_gmsh()
Generates a new mesh based on GMSH using .pos and .geo files for i/o.
real current_complexity()
Evaluates the current complexity of the mesh.
void field_p_hessian()
(NOT IMPLEMENTED) Generate polynomial distribution based on directinal derivatives ...
void field_hp_adjoint()
(NOT IMPLEMENTED) Generate mesh and polynomial distribution based on dual-weighted residual distribut...
void field_p_residual()
(NOT IMPLEMENTED) Generate polynomial distribution based on fine-grid residual distribution ...
void refine_grid() override
Perform call to the grid refinement object of choice.
void field_p_error()
(NOT IMPLEMENTED) Generate polynomial distribution based on manufactured solution ...
dealii::Vector< real > p_field
Continuous representation of the polynomial distribution.
Parameters related to individual grid refinement run.
real complexity_initial
Initial mesh complexity at construction.
void field_p_adjoint()
(NOT IMPLEMENTED) Generate polynomial distribution based on dual-weighted residual distribution ...
void field_h_adjoint()
Generate mesh distribution based on dual-weighted residual distribution.
std::vector< std::pair< dealii::Vector< real >, std::string > > output_results_vtk_method() override
Output refinement method dependent results.
void refine_grid_h()
Performs call to global remeshing method.
void field_h_hessian()
Generate mesh distribution based on directional derivatives.
Functional base class.
Definition: functional.h:43
void refine_grid_hp()
(NOT IMPLEMENTED) Performs call to global remeshing method with updated polynomial distribution ...
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
void get_current_field_h()
Evaluates the mesh size and ansitropy description for the current mesh.
std::vector< real > complexity_vector
Vector of complexity target goals for each iteration.
void field_h_residual()
(NOT IMPLEMENTED) Generate mesh distribution based on fine-grid residual distribution ...