[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
size_field.h
1 #ifndef __SIZE_FIELD_H__
2 #define __SIZE_FIELD_H__
3 
4 #include <deal.II/grid/tria.h>
5 
6 #include <deal.II/fe/fe.h>
7 #include <deal.II/fe/mapping.h>
8 
9 #include "grid_refinement/field.h"
10 
11 namespace PHiLiP {
12 
13 namespace GridRefinement {
14 
16 
63 template <int dim, typename real>
64 class SizeField
65 {
66 public:
68 
103  static void isotropic_uniform(
104  const real & complexity,
105  const dealii::Vector<real> & B,
106  const dealii::DoFHandler<dim> & dof_handler,
107  std::unique_ptr<Field<dim,real>> & h_field,
108  const real & poly_degree
109  );
110 
112 
144  static void isotropic_h(
145  const real complexity,
146  const dealii::Vector<real> & B,
147  const dealii::DoFHandler<dim> & dof_handler,
148  const dealii::hp::MappingCollection<dim> & mapping_collection,
149  const dealii::hp::FECollection<dim> & fe_collection,
150  const dealii::hp::QCollection<dim> & quadrature_collection,
151  const dealii::UpdateFlags & update_flags,
152  std::unique_ptr<Field<dim,real>> & h_field,
153  const dealii::Vector<real> & p_field
154  );
155 
156  // computes updated p-field with a constant h-field
157  // NOT IMPLEMENTED yet
158  /*
159  static void isotropic_p(
160  const dealii::Vector<real> & Bm, ///< Continuous error model for \f$p\f$ directional derivatives
161  const dealii::Vector<real> & B, ///< Continuous error model for \f$p+1\f$ directional derivatives
162  const dealii::Vector<real> & Bp, ///< Continuous error model for \f$p+2\f$ directional derivatives
163  const dealii::DoFHandler<dim> & dof_handler, ///< DoFHandler describing the mesh
164  const dealii::hp::MappingCollection<dim> & mapping_collection, ///< Element mapping collection
165  const dealii::hp::FECollection<dim> & fe_collection, ///< Finite element collection
166  const dealii::hp::QCollection<dim> & quadrature_collection, ///< Quadrature rules collection
167  const dealii::UpdateFlags & update_flags, ///< Update flags for the volume finite elements
168  const std::unique_ptr<Field<dim,real>> & h_field, ///< (Input) Current size-field
169  dealii::Vector<real> & p_field); ///< (Output) Target polynomial field
170  */
171 
173 
226  static void isotropic_hp(
227  const real complexity,
228  const dealii::Vector<real> & Bm,
229  const dealii::Vector<real> & B,
230  const dealii::Vector<real> & Bp,
231  const dealii::DoFHandler<dim> & dof_handler,
232  const dealii::hp::MappingCollection<dim> & mapping_collection,
233  const dealii::hp::FECollection<dim> & fe_collection,
234  const dealii::hp::QCollection<dim> & quadrature_collection,
235  const dealii::UpdateFlags & update_flags,
236  std::unique_ptr<Field<dim,real>> & h_field,
237  dealii::Vector<real> & p_field
238  );
239 
241 
245  static void adjoint_uniform_balan(
246  const real complexity,
247  const real r_max,
248  const real c_max,
249  const dealii::Vector<real> & eta,
250  const dealii::DoFHandler<dim> & dof_handler,
251  const dealii::hp::MappingCollection<dim> & mapping_collection,
252  const dealii::hp::FECollection<dim> & fe_collection,
253  const dealii::hp::QCollection<dim> & quadrature_collection,
254  const dealii::UpdateFlags & update_flags,
255  std::unique_ptr<Field<dim,real>>& h_field,
256  const real & poly_degree
257  );
258 
260 
272  static void adjoint_h_balan(
273  const real complexity,
274  const real r_max,
275  const real c_max,
276  const dealii::Vector<real> & eta,
277  const dealii::DoFHandler<dim> & dof_handler,
278  const dealii::hp::MappingCollection<dim> & mapping_collection,
279  const dealii::hp::FECollection<dim> & fe_collection,
280  const dealii::hp::QCollection<dim> & quadrature_collection,
281  const dealii::UpdateFlags & update_flags,
282  std::unique_ptr<Field<dim,real>>& h_field,
283  const dealii::Vector<real> & p_field
284  );
285 
287 
296  static void adjoint_h_equal(
297  const real complexity,
298  const dealii::Vector<real> & eta,
299  const dealii::DoFHandler<dim> & dof_handler,
300  const dealii::hp::MappingCollection<dim> & mapping_collection,
301  const dealii::hp::FECollection<dim> & fe_collection,
302  const dealii::hp::QCollection<dim> & quadrature_collection,
303  const dealii::UpdateFlags & update_flags,
304  std::unique_ptr<Field<dim,real>>& h_field,
305  const real & poly_degree
306  );
307 
308 protected:
309 
311 
324  static void update_h_dwr(
325  const real tau,
326  const dealii::Vector<real> & eta,
327  const dealii::DoFHandler<dim> & dof_handler,
328  std::unique_ptr<Field<dim,real>>& h_field,
329  const real & poly_degree
330  );
331 
333 
348  static void update_h_optimal(
349  const real lambda,
350  const dealii::Vector<real> & B,
351  const dealii::DoFHandler<dim> & dof_handler,
352  std::unique_ptr<Field<dim,real>> & h_field,
353  const dealii::Vector<real> & p_field
354  );
355 
357 
383  static real evaluate_complexity(
384  const dealii::DoFHandler<dim> & dof_handler,
385  const dealii::hp::MappingCollection<dim> & mapping_collection,
386  const dealii::hp::FECollection<dim> & fe_collection,
387  const dealii::hp::QCollection<dim> & quadrature_collection,
388  const dealii::UpdateFlags & update_flags,
389  const std::unique_ptr<Field<dim,real>> & h_field,
390  const dealii::Vector<real> & p_field
391  );
392 
394 
412  static void update_alpha_vector_balan(
413  const dealii::Vector<real>& eta,
414  const real r_max,
415  const real c_max,
416  const real eta_min,
417  const real eta_max,
418  const real eta_ref,
419  const dealii::DoFHandler<dim>& dof_handler,
420  const dealii::Vector<real>& I_c,
421  std::unique_ptr<Field<dim,real>>& h_field
422  );
423 
425 
453  static real update_alpha_k_balan(
454  const real eta_k,
455  const real r_max,
456  const real c_max,
457  const real eta_min,
458  const real eta_max,
459  const real eta_ref
460  );
461 
463 
469  static real bisection(
470  const std::function<real(real)> func,
471  real lower_bound,
472  real upper_bound,
473  real rel_tolerance = 1e-6,
474  real abs_tolerance = 1.0
475  );
476 
477 };
478 
479 } // namespace GridRefinement
480 
481 } //namespace PHiLiP
482 
483 #endif // __SIZE_FIELD_H__
static void adjoint_uniform_balan(const real complexity, const real r_max, const real c_max, const dealii::Vector< real > &eta, const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::MappingCollection< dim > &mapping_collection, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &quadrature_collection, const dealii::UpdateFlags &update_flags, std::unique_ptr< Field< dim, real >> &h_field, const real &poly_degree)
Performs adjoint based size field adaptation for uniform polynomial distribution based on the method ...
Definition: size_field.cpp:272
static void isotropic_hp(const real complexity, const dealii::Vector< real > &Bm, const dealii::Vector< real > &B, const dealii::Vector< real > &Bp, const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::MappingCollection< dim > &mapping_collection, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &quadrature_collection, const dealii::UpdateFlags &update_flags, std::unique_ptr< Field< dim, real >> &h_field, dealii::Vector< real > &p_field)
Computes the size field (element scale) and updated polynomial distrubition for high-order error func...
Definition: size_field.cpp:206
static void isotropic_h(const real complexity, const dealii::Vector< real > &B, const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::MappingCollection< dim > &mapping_collection, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &quadrature_collection, const dealii::UpdateFlags &update_flags, std::unique_ptr< Field< dim, real >> &h_field, const dealii::Vector< real > &p_field)
Computes the size field (element scale) for a potentially non-uniform - distribution.
Definition: size_field.cpp:70
Files for the baseline physics.
Definition: ADTypes.hpp:10
static real bisection(const std::function< real(real)> func, real lower_bound, real upper_bound, real rel_tolerance=1e-6, real abs_tolerance=1.0)
Bisect function based on starting bounds.
Definition: size_field.cpp:619
static void update_h_optimal(const real lambda, const dealii::Vector< real > &B, const dealii::DoFHandler< dim > &dof_handler, std::unique_ptr< Field< dim, real >> &h_field, const dealii::Vector< real > &p_field)
Based on a chosen bisection paraemeter, updates the - field to an optimal distribution.
Definition: size_field.cpp:143
static void update_h_dwr(const real tau, const dealii::Vector< real > &eta, const dealii::DoFHandler< dim > &dof_handler, std::unique_ptr< Field< dim, real >> &h_field, const real &poly_degree)
Updates the size field based on the DWR disitribution and a reference value, .
Definition: size_field.cpp:507
static real evaluate_complexity(const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::MappingCollection< dim > &mapping_collection, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &quadrature_collection, const dealii::UpdateFlags &update_flags, const std::unique_ptr< Field< dim, real >> &h_field, const dealii::Vector< real > &p_field)
Evaluates the continuous complexity (approximate degrees of freedom) for the target mesh representat...
Definition: size_field.cpp:99
Field element class.
Definition: field.h:452
static void isotropic_uniform(const real &complexity, const dealii::Vector< real > &B, const dealii::DoFHandler< dim > &dof_handler, std::unique_ptr< Field< dim, real >> &h_field, const real &poly_degree)
Computes the size field (element scale) for a uniform - distribution.
Definition: size_field.cpp:38
static void adjoint_h_balan(const real complexity, const real r_max, const real c_max, const dealii::Vector< real > &eta, const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::MappingCollection< dim > &mapping_collection, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &quadrature_collection, const dealii::UpdateFlags &update_flags, std::unique_ptr< Field< dim, real >> &h_field, const dealii::Vector< real > &p_field)
Performs adjoint based size field adaptation for non-uniform polynomial distribution based on the met...
Definition: size_field.cpp:307
static void update_alpha_vector_balan(const dealii::Vector< real > &eta, const real r_max, const real c_max, const real eta_min, const real eta_max, const real eta_ref, const dealii::DoFHandler< dim > &dof_handler, const dealii::Vector< real > &I_c, std::unique_ptr< Field< dim, real >> &h_field)
Updates the size field for the mesh through the area measure based on reference value.
Definition: size_field.cpp:535
Size Field Static Class.
Definition: size_field.h:64
static void adjoint_h_equal(const real complexity, const dealii::Vector< real > &eta, const dealii::DoFHandler< dim > &dof_handler, const dealii::hp::MappingCollection< dim > &mapping_collection, const dealii::hp::FECollection< dim > &fe_collection, const dealii::hp::QCollection< dim > &quadrature_collection, const dealii::UpdateFlags &update_flags, std::unique_ptr< Field< dim, real >> &h_field, const real &poly_degree)
Performs adjoint based size field adaptation with a uniform - distribution based on equidistribution ...
Definition: size_field.cpp:419
static real update_alpha_k_balan(const real eta_k, const real r_max, const real c_max, const real eta_min, const real eta_max, const real eta_ref)
Determines local sizing factor (from adjoint estimates)
Definition: size_field.cpp:577