1 #ifndef __SIZE_FIELD_H__ 2 #define __SIZE_FIELD_H__ 4 #include <deal.II/grid/tria.h> 6 #include <deal.II/fe/fe.h> 7 #include <deal.II/fe/mapping.h> 9 #include "grid_refinement/field.h" 13 namespace GridRefinement {
63 template <
int dim,
typename real>
104 const real & complexity,
105 const dealii::Vector<real> & B,
106 const dealii::DoFHandler<dim> & dof_handler,
108 const real & poly_degree
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,
153 const dealii::Vector<real> & p_field
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,
237 dealii::Vector<real> & p_field
246 const real complexity,
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,
256 const real & poly_degree
273 const real complexity,
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,
283 const dealii::Vector<real> & p_field
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,
305 const real & poly_degree
326 const dealii::Vector<real> & eta,
327 const dealii::DoFHandler<dim> & dof_handler,
329 const real & poly_degree
350 const dealii::Vector<real> & B,
351 const dealii::DoFHandler<dim> & dof_handler,
353 const dealii::Vector<real> & p_field
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,
390 const dealii::Vector<real> & p_field
413 const dealii::Vector<real>& eta,
419 const dealii::DoFHandler<dim>& dof_handler,
420 const dealii::Vector<real>& I_c,
470 const std::function<real(real)> func,
473 real rel_tolerance = 1e-6,
474 real abs_tolerance = 1.0
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 ...
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...
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.
Files for the baseline physics.
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.
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.
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, .
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...
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.
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...
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.
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 ...
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)