Reference documentation for deal.II version 9.2.0
|
#include <deal.II/fe/fe_series.h>
Public Types | |
using | CoefficientType = double |
Public Member Functions | |
Legendre (const std::vector< unsigned int > &n_coefficients_per_direction, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim > &q_collection) | |
Legendre (const unsigned int n_coefficients_per_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 > &legendre_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 Legendre< dim, spacedim > &legendre) const |
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 | |
const std::vector< unsigned int > | n_coefficients_per_direction |
SmartPointer< const hp::FECollection< dim, spacedim > > | fe_collection |
const hp::QCollection< dim > | q_collection |
std::vector< FullMatrix< CoefficientType > > | legendre_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 series of Legendre functions on a reference element.
Legendre functions are solutions to Legendre's differential equation
\[ \frac{d}{dx}\left([1-x^2] \frac{d}{dx} P_n(x)\right) + n[n+1] P_n(x) = 0 \]
and can be expressed using Rodrigues' formula
\[ P_n(x) = \frac{1}{2^n n!} \frac{d^n}{dx^n}[x^2-1]^n. \]
These polynomials are orthogonal with respect to the \( L^2 \) inner product on the interval \( [-1;1] \)
\[ \int_{-1}^1 P_m(x) P_n(x) = \frac{2}{2n + 1} \delta_{mn} \]
and are complete. A family of \( L^2 \)-orthogonal polynomials on \( [0;1] \) can be constructed via
\[ \widetilde P_m = \sqrt{2} P_m(2x-1). \]
An arbitrary scalar FE field on the reference element \( [0;1] \) can be expanded in the complete orthogonal basis as
\[ u(x) = \sum_{m} c_m \widetilde P_{m}(x). \]
From the orthogonality property of the basis, it follows that
\[ c_m = \frac{2m+1}{2} \int_0^1 u(x) \widetilde P_m(x) dx . \]
This class calculates coefficients \( c_{\bf k} \) using \( dim \)-dimensional Legendre polynomials constructed from \( \widetilde P_m(x) \) using tensor product rule.
Definition at line 259 of file fe_series.h.
using FESeries::Legendre< dim, spacedim >::CoefficientType = double |
Definition at line 262 of file fe_series.h.
FESeries::Legendre< dim, spacedim >::Legendre | ( | const std::vector< unsigned int > & | n_coefficients_per_direction, |
const hp::FECollection< dim, spacedim > & | fe_collection, | ||
const hp::QCollection< dim > & | q_collection | ||
) |
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
.
Definition at line 198 of file fe_series_legendre.cc.
FESeries::Legendre< dim, spacedim >::Legendre | ( | const unsigned int | n_coefficients_per_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 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
.
Definition at line 222 of file fe_series_legendre.cc.
void FESeries::Legendre< dim, spacedim >::calculate | ( | const ::Vector< Number > & | local_dof_values, |
const unsigned int | cell_active_fe_index, | ||
Table< dim, CoefficientType > & | legendre_coefficients | ||
) |
Calculate legendre_coefficients
of the cell vector field given by local_dof_values
corresponding to FiniteElement with cell_active_fe_index
.
Definition at line 281 of file fe_series_legendre.cc.
unsigned int FESeries::Legendre< 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 270 of file fe_series_legendre.cc.
void FESeries::Legendre< 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 251 of file fe_series_legendre.cc.
void FESeries::Legendre< 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.
void FESeries::Legendre< dim, spacedim >::load_transformation_matrices | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read all transformation matrices from a stream and recover them for this object.
|
inline |
Test for equality of two series expansion objects.
Definition at line 238 of file fe_series_legendre.cc.
|
private |
Number of coefficients in each direction for each finite element in the registered hp::FECollection.
Definition at line 352 of file fe_series.h.
|
private |
hp::FECollection for which transformation matrices will be calculated.
Definition at line 357 of file fe_series.h.
|
private |
hp::QCollection used in calculation of transformation matrices.
Definition at line 362 of file fe_series.h.
|
private |
Transformation matrices for each FiniteElement.
Definition at line 367 of file fe_series.h.
|
private |
Auxiliary vector to store unrolled coefficients.
Definition at line 372 of file fe_series.h.