Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
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 Member Functions

 Fourier (const unsigned int size_in_each_direction, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim > &q_collection)
 
template<typename Number >
void calculate (const ::Vector< Number > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, CoefficientType > &fourier_coefficients)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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)
 

Private Attributes

SmartPointer< const hp::FECollection< dim, spacedim > > fe_collection
 
SmartPointer< const hp::QCollection< dim > > q_collection
 
Table< dim, Tensor< 1, dim > > k_vectors
 
std::vector< FullMatrix< CoefficientType > > fourier_transform_matrices
 
std::vector< CoefficientType > unrolled_coefficients
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

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

A class to calculate expansion of a scalar 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 \) .

Author
Denis Davydov, 2016.

Definition at line 91 of file fe_series.h.

Constructor & Destructor Documentation

◆ Fourier()

template<int dim, int spacedim>
FESeries::Fourier< dim, spacedim >::Fourier ( const unsigned int  size_in_each_direction,
const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim > &  q_collection 
)

A non-default constructor. The size_in_each_direction defines the number of modes 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.

Definition at line 170 of file fe_series_fourier.cc.

Member Function Documentation

◆ calculate()

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 .

Member Data Documentation

◆ 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 122 of file fe_series.h.

◆ q_collection

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

hp::QCollection used in calculation of transformation matrices.

Definition at line 127 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 132 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 137 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 142 of file fe_series.h.


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