Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
smoothness_estimator.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2018 - 2023 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_smoothness_estimator_h
17#define dealii_smoothness_estimator_h
18
19
20#include <deal.II/base/config.h>
21
23
25
26#include <functional>
27#include <vector>
28
29
31
32
33// forward declarations
34#ifndef DOXYGEN
35template <typename Number>
36class Vector;
37
38template <int dim, int spacedim>
40class DoFHandler;
41
42namespace FESeries
43{
44 template <int dim, int spacedim>
45 class Fourier;
46 template <int dim, int spacedim>
47 class Legendre;
48} // namespace FESeries
49
50namespace hp
51{
52 template <int dim, int spacedim>
53 class FECollection;
54} // namespace hp
55#endif
56
57
70{
107 namespace Legendre
108 {
160 template <int dim, int spacedim, typename VectorType>
161 void
163 const DoFHandler<dim, spacedim> & dof_handler,
164 const VectorType & solution,
165 Vector<float> & smoothness_indicators,
166 const VectorTools::NormType regression_strategy =
168 const double smallest_abs_coefficient = 1e-10,
169 const bool only_flagged_cells = false);
170
216 template <int dim, int spacedim, typename VectorType>
217 void
220 const DoFHandler<dim, spacedim> & dof_handler,
221 const VectorType & solution,
222 Vector<float> & smoothness_indicators,
223 const ComponentMask &coefficients_predicate = ComponentMask(),
224 const double smallest_abs_coefficient = 1e-10,
225 const bool only_flagged_cells = false);
226
247 template <int dim, int spacedim>
250 const hp::FECollection<dim, spacedim> &fe_collection,
251 const unsigned int component = numbers::invalid_unsigned_int);
252 } // namespace Legendre
253
254
255
360 namespace Fourier
361 {
405 template <int dim, int spacedim, typename VectorType>
406 void
408 const DoFHandler<dim, spacedim> & dof_handler,
409 const VectorType & solution,
410 Vector<float> & smoothness_indicators,
411 const VectorTools::NormType regression_strategy =
413 const double smallest_abs_coefficient = 1e-10,
414 const bool only_flagged_cells = false);
415
454 template <int dim, int spacedim, typename VectorType>
455 void
458 const DoFHandler<dim, spacedim> & dof_handler,
459 const VectorType & solution,
460 Vector<float> & smoothness_indicators,
461 const ComponentMask &coefficients_predicate = ComponentMask(),
462 const double smallest_abs_coefficient = 1e-10,
463 const bool only_flagged_cells = false);
464
485 template <int dim, int spacedim>
488 const hp::FECollection<dim, spacedim> &fe_collection,
489 const unsigned int component = numbers::invalid_unsigned_int);
490 } // namespace Fourier
491} // namespace SmoothnessEstimator
492
493
495
496#endif
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:160
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
void coefficient_decay_per_direction(FESeries::Fourier< dim, spacedim > &fe_fourier, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &smoothness_indicators, const ComponentMask &coefficients_predicate=ComponentMask(), const double smallest_abs_coefficient=1e-10, const bool only_flagged_cells=false)
FESeries::Fourier< dim, spacedim > default_fe_series(const hp::FECollection< dim, spacedim > &fe_collection, const unsigned int component=numbers::invalid_unsigned_int)
void coefficient_decay(FESeries::Fourier< dim, spacedim > &fe_fourier, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &smoothness_indicators, const VectorTools::NormType regression_strategy=VectorTools::Linfty_norm, const double smallest_abs_coefficient=1e-10, const bool only_flagged_cells=false)
void coefficient_decay_per_direction(FESeries::Legendre< dim, spacedim > &fe_legendre, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &smoothness_indicators, const ComponentMask &coefficients_predicate=ComponentMask(), const double smallest_abs_coefficient=1e-10, const bool only_flagged_cells=false)
FESeries::Legendre< dim, spacedim > default_fe_series(const hp::FECollection< dim, spacedim > &fe_collection, const unsigned int component=numbers::invalid_unsigned_int)
void coefficient_decay(FESeries::Legendre< dim, spacedim > &fe_legendre, const DoFHandler< dim, spacedim > &dof_handler, const VectorType &solution, Vector< float > &smoothness_indicators, const VectorTools::NormType regression_strategy=VectorTools::Linfty_norm, const double smallest_abs_coefficient=1e-10, const bool only_flagged_cells=false)
Definition hp.h:118
static const unsigned int invalid_unsigned_int
Definition types.h:213