![]() |
Reference documentation for deal.II version 8.5.1
|
#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) |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) 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::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static::ExceptionBase & | ExcNoSubscriber (char *arg1, char *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 . For example in 1D the L2-orthogonality condition reads
Note that .
The arbitrary scalar FE field on the reference element can be expanded in the complete orthogonal exponential basis as
From the orthogonality property of the basis, it follows that
It is this complex-valued expansion coefficients, that are calculated by this class. Note that , where
are real-valued FiniteElement shape functions. Consequently
and we only need to compute
for positive indices
.
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 .
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.