Reference documentation for deal.II version 9.4.1
\(\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
Public Types | Public Member Functions | Private Attributes | List of all members
FESeries::Fourier< dim, spacedim > Class Template Reference

#include <deal.II/fe/fe_series.h>

Inheritance diagram for FESeries::Fourier< dim, spacedim >:
[legend]

Public Types

using CoefficientType = typename std::complex< double >
 

Public Member Functions

 Fourier (const std::vector< unsigned int > &n_coefficients_per_direction, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim > &q_collection, const unsigned int component=numbers::invalid_unsigned_int)
 
template<typename Number >
void calculate (const ::Vector< Number > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, CoefficientType > &fourier_coefficients)
 
unsigned int get_n_coefficients_per_direction (const unsigned int index) const
 
void precalculate_all_transformation_matrices ()
 
template<class Archive >
void save_transformation_matrices (Archive &ar, const unsigned int version)
 
template<class Archive >
void load_transformation_matrices (Archive &ar, const unsigned int version)
 
bool operator== (const Fourier< dim, spacedim > &fourier) const
 
template<typename Number >
void calculate (const Vector< Number > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, CoefficientType > &fourier_coefficients)
 

Private Attributes

const std::vector< unsigned intn_coefficients_per_direction
 
SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
 
const hp::QCollection< dim > q_collection
 
Table< dim, Tensor< 1, dim > > k_vectors
 
std::vector< FullMatrix< CoefficientType > > fourier_transform_matrices
 
std::vector< CoefficientTypeunrolled_coefficients
 
const unsigned int component
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
static std::mutex mutex
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
void check_no_subscribers () const noexcept
 

Detailed Description

template<int dim, int spacedim = dim>
class FESeries::Fourier< dim, spacedim >

A class to calculate expansion of a scalar FE (or a single component of vector-valued FE) field into Fourier series on a reference element. The exponential form of the Fourier series is based on completeness and Hermitian orthogonality of the set of exponential functions \( \phi_{\bf k}({\bf x}) = \exp(2 \pi i\, {\bf k} \cdot {\bf x})\). For example in 1D the L2-orthogonality condition reads

\[ \int_0^1 \phi_k(x) \phi_l^\ast(x) dx=\delta_{kl}. \]

Note that \( \phi_{\bf k} = \phi_{-\bf k}^\ast \).

The arbitrary scalar FE field on the reference element can be expanded in the complete orthogonal exponential basis as

\[ u({\bf x}) = \sum_{\bf k} c_{\bf k} \phi_{\bf k}({\bf x}). \]

From the orthogonality property of the basis, it follows that

\[ c_{\bf k} = \int_{[0,1]^d} u({\bf x}) \phi_{\bf k}^\ast ({\bf x}) d{\bf x}\,. \]

It is this complex-valued expansion coefficients, that are calculated by this class. Note that \( u({\bf x}) = \sum_i u_i N_i({\bf x})\), where \( N_i({\bf x}) \) are real-valued FiniteElement shape functions. Consequently \( c_{\bf k} \equiv c_{-\bf k}^\ast \) and we only need to compute \( c_{\bf k} \) for positive indices \( \bf k \) .

Definition at line 88 of file fe_series.h.

Member Typedef Documentation

◆ CoefficientType

template<int dim, int spacedim = dim>
using FESeries::Fourier< dim, spacedim >::CoefficientType = typename std::complex<double>

Definition at line 91 of file fe_series.h.

Constructor & Destructor Documentation

◆ Fourier()

template<int dim, int spacedim>
FESeries::Fourier< dim, spacedim >::Fourier ( const std::vector< unsigned int > &  n_coefficients_per_direction,
const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim > &  q_collection,
const unsigned int  component = numbers::invalid_unsigned_int 
)

Constructor that initializes all required data structures.

The n_coefficients_per_direction defines the number of coefficients in each direction, fe_collection is the hp::FECollection for which expansion will be used and q_collection is the hp::QCollection used to integrate the expansion for each FiniteElement in fe_collection.

As the Fourier 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 Fourier 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 192 of file fe_series_fourier.cc.

Member Function Documentation

◆ calculate() [1/2]

template<int dim, int spacedim = dim>
template<typename Number >
void FESeries::Fourier< dim, spacedim >::calculate ( const ::Vector< Number > &  local_dof_values,
const unsigned int  cell_active_fe_index,
Table< dim, CoefficientType > &  fourier_coefficients 
)

Calculate fourier_coefficients of the cell vector field given by local_dof_values corresponding to FiniteElement with cell_active_fe_index .

◆ get_n_coefficients_per_direction()

template<int dim, int spacedim>
unsigned int FESeries::Fourier< dim, spacedim >::get_n_coefficients_per_direction ( const unsigned int  index) const

Return the number of coefficients in each coordinate direction for the finite element associated with index in the provided hp::FECollection.

Definition at line 267 of file fe_series_fourier.cc.

◆ precalculate_all_transformation_matrices()

template<int dim, int spacedim>
void FESeries::Fourier< dim, spacedim >::precalculate_all_transformation_matrices

Calculate all transformation matrices to transfer the finite element solution to the series expansion representation.

These matrices will be generated on demand by calling calculate() and stored for recurring purposes. Usually, this operation consumes a lot of workload. With this function, all matrices will be calculated in advance. This way, we can separate their costly generation from the actual application.

Definition at line 246 of file fe_series_fourier.cc.

◆ save_transformation_matrices()

template<int dim, int spacedim = dim>
template<class Archive >
void FESeries::Fourier< dim, spacedim >::save_transformation_matrices ( Archive &  ar,
const unsigned int  version 
)

Write all transformation matrices of this object to a stream for the purpose of serialization.

Since any of its transformation matrices has to be generated only once for a given scenario, it is common practice to determine them in advance calling precalculate_all_transformation_matrices() and keep them via serialization.

◆ load_transformation_matrices()

template<int dim, int spacedim = dim>
template<class Archive >
void FESeries::Fourier< dim, spacedim >::load_transformation_matrices ( Archive &  ar,
const unsigned int  version 
)

Read all transformation matrices from a stream and recover them for this object.

◆ operator==()

template<int dim, int spacedim>
bool FESeries::Fourier< dim, spacedim >::operator== ( const Fourier< dim, spacedim > &  fourier) const
inline

Test for equality of two series expansion objects.

Definition at line 230 of file fe_series_fourier.cc.

◆ calculate() [2/2]

template<int dim, int spacedim = dim>
template<typename Number >
void FESeries::Fourier< dim, spacedim >::calculate ( const Vector< Number > &  local_dof_values,
const unsigned int  cell_active_fe_index,
Table< dim, CoefficientType > &  fourier_coefficients 
)

Definition at line 278 of file fe_series_fourier.cc.

Member Data Documentation

◆ n_coefficients_per_direction

template<int dim, int spacedim = dim>
const std::vector<unsigned int> FESeries::Fourier< dim, spacedim >::n_coefficients_per_direction
private

Number of coefficients in each direction for each finite element in the registered hp::FECollection.

Definition at line 179 of file fe_series.h.

◆ fe_collection

template<int dim, int spacedim = dim>
SmartPointer<const hp::FECollection<dim, spacedim> > FESeries::Fourier< dim, spacedim >::fe_collection
private

hp::FECollection for which transformation matrices will be calculated.

Definition at line 184 of file fe_series.h.

◆ q_collection

template<int dim, int spacedim = dim>
const hp::QCollection<dim> FESeries::Fourier< dim, spacedim >::q_collection
private

hp::QCollection used in calculation of transformation matrices.

Definition at line 189 of file fe_series.h.

◆ k_vectors

template<int dim, int spacedim = dim>
Table<dim, Tensor<1, dim> > FESeries::Fourier< dim, spacedim >::k_vectors
private

Angular frequencies \( 2 \pi {\bf k} \) .

Definition at line 194 of file fe_series.h.

◆ fourier_transform_matrices

template<int dim, int spacedim = dim>
std::vector<FullMatrix<CoefficientType> > FESeries::Fourier< dim, spacedim >::fourier_transform_matrices
private

Transformation matrices for each FiniteElement.

Definition at line 199 of file fe_series.h.

◆ unrolled_coefficients

template<int dim, int spacedim = dim>
std::vector<CoefficientType> FESeries::Fourier< dim, spacedim >::unrolled_coefficients
private

Auxiliary vector to store unrolled coefficients.

Definition at line 204 of file fe_series.h.

◆ component

template<int dim, int spacedim = dim>
const unsigned int FESeries::Fourier< dim, spacedim >::component
private

Which component of FiniteElement should be used to calculate the expansion.

Definition at line 210 of file fe_series.h.


The documentation for this class was generated from the following files: