Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/fe_series.h>
Public Member Functions | |
Fourier (const unsigned int size_in_each_direction, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim > &q_collection) | |
void | calculate (const ::Vector< double > &local_dof_values, const unsigned int cell_active_fe_index, Table< dim, std::complex< double > > &fourier_coefficients) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
void | ensure_existence (const unsigned int fe_index) |
Private Attributes | |
SmartPointer< const hp::FECollection< dim > > | fe_collection |
SmartPointer< const hp::QCollection< dim > > | q_collection |
Table< dim, Tensor< 1, dim > > | k_vectors |
std::vector< FullMatrix< std::complex< double > > > | fourier_transform_matrices |
std::vector< std::complex< double > > | unrolled_coefficients |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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 \) .
Definition at line 87 of file fe_series.h.
FESeries::Fourier< dim >::Fourier | ( | const unsigned int | size_in_each_direction, |
const hp::FECollection< dim > & | 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 74 of file fe_series.cc.
void FESeries::Fourier< dim >::calculate | ( | const ::Vector< double > & | local_dof_values, |
const unsigned int | cell_active_fe_index, | ||
Table< dim, std::complex< double > > & | fourier_coefficients | ||
) |
Calculate fourier_coefficients
of the cell vector field given by local_dof_values
corresponding to FiniteElement with cell_active_fe_index
.
Definition at line 87 of file fe_series.cc.
|
private |
Ensure that the transformation matrix for FiniteElement index fe_index
is calculated. If not, calculate it.
|
private |
hp::FECollection for which transformation matrices will be calculated.
Definition at line 114 of file fe_series.h.
|
private |
hp::QCollection used in calculation of transformation matrices.
Definition at line 119 of file fe_series.h.
|
private |
Angular frequencies \( 2 \pi {\bf k} \) .
Definition at line 130 of file fe_series.h.
|
private |
Transformation matrices for each FiniteElement.
Definition at line 135 of file fe_series.h.
|
private |
Auxiliary vector to store unrolled coefficients.
Definition at line 140 of file fe_series.h.