|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_hp_refinement_h
17 #define dealii_hp_refinement_h
33 template <
typename Number>
38 template <
int dim,
int spacedim>
141 template <
typename Number>
143 std::function<
bool(
const Number &,
const Number &)>;
158 template <
int dim,
int spacedim>
175 template <
int dim,
int spacedim>
178 const std::vector<bool> & p_flags);
204 template <
int dim,
typename Number,
int spacedim>
208 const Vector<Number> & criteria,
209 const Number p_refine_threshold,
210 const Number p_coarsen_threshold,
212 &compare_refine = std::greater_equal<Number>(),
214 &compare_coarsen = std::less_equal<Number>());
246 template <
int dim,
typename Number,
int spacedim>
250 const Vector<Number> & criteria,
251 const double p_refine_fraction = 0.5,
252 const double p_coarsen_fraction = 0.5,
254 &compare_refine = std::greater_equal<Number>(),
256 &compare_coarsen = std::less_equal<Number>());
289 template <
int dim,
typename Number,
int spacedim>
293 const Vector<Number> & criteria,
294 const double p_refine_fraction = 0.5,
295 const double p_coarsen_fraction = 0.5,
297 &compare_refine = std::greater_equal<Number>(),
299 &compare_coarsen = std::less_equal<Number>());
326 template <
int dim,
typename Number,
int spacedim>
330 const Vector<Number> & sobolev_indices);
351 template <
int dim,
typename Number,
int spacedim>
355 const Vector<Number> & criteria,
356 const Vector<Number> & references,
561 template <
int dim,
typename Number,
int spacedim>
564 const Vector<Number> & error_indicators,
565 Vector<Number> & predicted_errors,
566 const double gamma_p = std::sqrt(0.4),
567 const double gamma_h = 2.,
568 const double gamma_n = 1.);
589 template <
int dim,
int spacedim>
636 template <
int dim,
int spacedim>
649 #endif // dealii_hp_refinement_h
void p_adaptivity_from_relative_threshold(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const double p_refine_fraction=0.5, const double p_coarsen_fraction=0.5, const ComparisonFunction< typename identity< Number >::type > &compare_refine=std::greater_equal< Number >(), const ComparisonFunction< typename identity< Number >::type > &compare_coarsen=std::less_equal< Number >())
void p_adaptivity_from_regularity(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &sobolev_indices)
void predict_error(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &error_indicators, Vector< Number > &predicted_errors, const double gamma_p=std::sqrt(0.4), const double gamma_h=2., const double gamma_n=1.)
void full_p_adaptivity(const hp::DoFHandler< dim, spacedim > &dof_handler)
std::function< bool(const Number &, const Number &)> ComparisonFunction
void p_adaptivity_from_flags(const hp::DoFHandler< dim, spacedim > &dof_handler, const std::vector< bool > &p_flags)
#define DEAL_II_NAMESPACE_OPEN
void choose_p_over_h(const hp::DoFHandler< dim, spacedim > &dof_handler)
void p_adaptivity_fixed_number(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const double p_refine_fraction=0.5, const double p_coarsen_fraction=0.5, const ComparisonFunction< typename identity< Number >::type > &compare_refine=std::greater_equal< Number >(), const ComparisonFunction< typename identity< Number >::type > &compare_coarsen=std::less_equal< Number >())
void p_adaptivity_from_absolute_threshold(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const Number p_refine_threshold, const Number p_coarsen_threshold, const ComparisonFunction< typename identity< Number >::type > &compare_refine=std::greater_equal< Number >(), const ComparisonFunction< typename identity< Number >::type > &compare_coarsen=std::less_equal< Number >())
void p_adaptivity_from_reference(const hp::DoFHandler< dim, spacedim > &dof_handler, const Vector< Number > &criteria, const Vector< Number > &references, const ComparisonFunction< typename identity< Number >::type > &compare_refine, const ComparisonFunction< typename identity< Number >::type > &compare_coarsen)
void force_p_over_h(const hp::DoFHandler< dim, spacedim > &dof_handler)
#define DEAL_II_NAMESPACE_CLOSE