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
Functions
SmoothnessEstimator::Legendre Namespace Reference

Functions

template<int dim, int spacedim, typename VectorType >
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)
 
template<int dim, int spacedim, typename VectorType >
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)
 
template<int dim, int spacedim>
FESeries::Legendre< dim, spacedim > default_fe_series (const hp::FECollection< dim, spacedim > &fe_collection, const unsigned int component=numbers::invalid_unsigned_int)
 

Detailed Description

Smoothness estimation strategy based on the decay of Legendre expansion coefficients.

In one dimension, the finite element solution on cell \(K\) with polynomial degree \(p\) can be written as

\begin{eqnarray*} u_h(x) &=& \sum_j u_j \varphi_j (x) \\ u_{h, k}(x) &=& \sum_{k=0}^{p} a_k \widetilde P_k (x), \quad a_k = \sum_j {\cal L}_{k,j} u_j \end{eqnarray*}

where \(u_j\) are degrees of freedom and \(\varphi_j\) are the corresponding shape functions. \(\{\widetilde P_k(x)\}\) are Legendre polynomials on cell \(K\). \(a_k\) and \({\cal L}_{k,j}\) are coefficients and transformation matrices from the Legendre expansion of each shape function. For practical reasons, we will perform the calculation of these matrices and coefficients only on the reference cell \(\hat K\). We only have to calculate the transformation matrices once this way. However, results are only applicable if the mapping from the reference cell to the actual cell is affine. We use the class FESeries::Legendre to determine all coefficients \(a_k\).

A function is analytic, i.e., representable by a power series, if and only if their Legendre expansion coefficients decay as (see [64])

\[ |a_k| \sim c \, \exp(-\sigma k) \]

We determine their decay rate \(\sigma\) by performing the linear regression fit of

\[ \ln |a_k| \sim C - \sigma k \]

for \(k=0,\ldots,p\), with \(p\) the polynomial degree of the finite element. The rate of the decay \(\sigma\) can be used to estimate the smoothness. For example, one strategy to implement hp-refinement criteria is to perform p-refinement if \(\sigma>1\) (see [127]).

Function Documentation

◆ coefficient_decay()

template<int dim, int spacedim, typename VectorType >
void SmoothnessEstimator::Legendre::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 
)

In this variant of the estimation strategy for higher dimensions, we will consider all mode vectors \(\bf k\) describing Legendre polynomials \(\widetilde P_{\bf k}\) and perform one least-squares fit over all coefficients at once. If there are multiple coefficients corresponding to the same absolute value of modes \(\|{\bf k}\|_1\), we take the maximum among those. Thus, the least-squares fit is performed on

\begin{eqnarray*} \widetilde P_{\bf k}({\bf x}) &=& \widetilde P_{k_1} (x_1) \ldots \widetilde P_{k_d} (x_d) \\ \ln \left( \max\limits_{\|{\bf k}\|_1} |a_{\bf k}| \right) &\sim& C - \sigma \|{\bf k}\|_1 \end{eqnarray*}

for \({\bf k}=(k_1,\ldots,k_d)\) and \(k_i=0,\ldots,p\), with \(p\) the polynomial degree of the finite element.

For a finite element approximation solution, this function writes the decay rate for every cell into the output vector smoothness_indicators.

Parameters
[in]fe_legendreFESeries::Legendre object to calculate coefficients. This object needs to be initialized to have at least \(p+1\) coefficients in each direction for every finite element in the collection, where \(p\) is its polynomial degree.
[in]dof_handlerA DoFHandler.
[in]solutionA solution vector.
[out]smoothness_indicatorsA vector for smoothness indicators.
[in]regression_strategyDetermines which norm will be used on the subset of coefficients \(\mathbf{k}\) with the same absolute value \(\|{\bf k}\|_1\). Default is VectorTools::Linfty_norm for a maximum approximation.
[in]smallest_abs_coefficientThe smallest absolute value of the coefficient to be used in linear regression. Note that Legendre coefficients of some functions may have a repeating pattern of zero coefficients (i.e. for functions that are locally symmetric or antisymmetric about the midpoint of the element in any coordinate direction). Thus this parameters allows to ignore small (in absolute value) coefficients within the linear regression fit. In case there are less than two nonzero coefficients, the returned value for this cell will be \(\sigma=\infty\).
[in]only_flagged_cellsSmoothness indicators are usually used to decide whether to perform h- or p-adaptation. So in most cases, we only need to calculate those indicators on cells flagged for refinement or coarsening. This parameter controls whether this particular subset or all cells will be considered. By default, all cells will be considered. When only flagged cells are supposed to be considered, smoothness indicators will only be set on those vector entries of flagged cells; the others will be set to a signaling NaN.

For more theoretical details see [127] [101] [64].

Definition at line 103 of file smoothness_estimator.cc.

◆ coefficient_decay_per_direction()

template<int dim, int spacedim, typename VectorType >
void SmoothnessEstimator::Legendre::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 
)

In this variant of the estimation strategy for higher dimensions, we only consider modes in each coordinate direction, i.e., only mode vectors \(\bf k\) with one nonzero entry. We perform the least-squares fit in each coordinate direction separately and take the lowest decay rate \(\sigma\) among them.

For a finite element approximation solution, this function writes the decay rate for every cell into the output vector smoothness_indicators.

Parameters
[in]fe_legendreFESeries::Legendre object to calculate coefficients. This object needs to be initialized to have at least \(p+1\) coefficients in each direction, where \(p\) is the maximum polynomial degree to be used.
[in]dof_handlerA DoFHandler
[in]solutionA solution vector
[out]smoothness_indicatorsA vector for smoothness indicators
[in]coefficients_predicateA predicate to select Legendre coefficients \(a_j\), \(j=0,\ldots,p\) for linear regression in each coordinate direction. The user is responsible for updating the vector of flags provided to this function. Note that its size is \(p+1\), where \(p\) is the polynomial degree of the FE basis on a given element. The default implementation will use all Legendre coefficients in each coordinate direction, i.e. set all elements of the vector to true.
[in]smallest_abs_coefficientThe smallest absolute value of the coefficient to be used in linear regression in each coordinate direction. Note that Legendre coefficients of some functions may have a repeating pattern of zero coefficients (i.e. for functions that are locally symmetric or antisymmetric about the midpoint of the element in any coordinate direction). Thus this parameters allows to ignore small (in absolute value) coefficients within the linear regression fit. In case there are less than two nonzero coefficients for a coordinate direction, this direction will be skipped. If all coefficients are zero, the returned value for this cell will be \(\sigma=\infty\).
[in]only_flagged_cellsSmoothness indicators are usually used to decide whether to perform h- or p-adaptation. So in most cases, we only need to calculate those indicators on cells flagged for refinement or coarsening. This parameter controls whether this particular subset or all cells will be considered. By default, all cells will be considered. When only flagged cells are supposed to be considered, smoothness indicators will only be set on those vector entries of flagged cells; the others will be set to NaN.

For more theoretical details and the application within the deal.II library see [59].

Definition at line 183 of file smoothness_estimator.cc.

◆ default_fe_series()

template<int dim, int spacedim>
FESeries::Legendre< dim, spacedim > SmoothnessEstimator::Legendre::default_fe_series ( const hp::FECollection< dim, spacedim > &  fe_collection,
const unsigned int  component = numbers::invalid_unsigned_int 
)

Returns a FESeries::Legendre object for Legendre series expansions with the default configuration for smoothness estimation purposes.

For each finite element of the provided fe_collection, we use as many modes as its polynomial degree plus two. This includes the first Legendre polynomial which is just a constant. Further for each element, we use a Gaussian quadrature designed to yield exact results for the highest order Legendre polynomial used.

As the Legendre expansion can only be performed on scalar fields, this class does not operate on vector-valued finite elements and will therefore throw an assertion. However, each component of a finite element field can be treated as a scalar field, respectively, on which Legendre expansions are again possible. For this purpose, the optional parameter component defines which component of each FiniteElement will be used. The default value of component only applies to scalar FEs, in which case it indicates that the sole component is to be decomposed. For vector-valued FEs, a non-default value must be explicitly provided.

Definition at line 290 of file smoothness_estimator.cc.