![]() |
Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/fe_series.h>
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 () |
| Subscriptor & | operator= (const Subscriptor &) |
| Subscriptor & | operator= (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 ::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
. 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 91 of file fe_series.h.
| 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.
| 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 .
|
private |
hp::FECollection for which transformation matrices will be calculated.
Definition at line 122 of file fe_series.h.
|
private |
hp::QCollection used in calculation of transformation matrices.
Definition at line 127 of file fe_series.h.
|
private |
Angular frequencies
.
Definition at line 132 of file fe_series.h.
|
private |
Transformation matrices for each FiniteElement.
Definition at line 137 of file fe_series.h.
|
private |
Auxiliary vector to store unrolled coefficients.
Definition at line 142 of file fe_series.h.
1.8.16