Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/fe_series.h>
Public Member Functions | |
Legendre (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 > &legendre_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 | |
const unsigned int | N |
SmartPointer< const hp::FECollection< dim, spacedim > > | fe_collection |
SmartPointer< 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 189 of file fe_series.h.
FESeries::Legendre< dim, spacedim >::Legendre | ( | 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 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 191 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 207 of file fe_series_legendre.cc.
|
private |
Number of coefficients in each direction
Definition at line 220 of file fe_series.h.
|
private |
hp::FECollection for which transformation matrices will be calculated.
Definition at line 225 of file fe_series.h.
|
private |
hp::QCollection used in calculation of transformation matrices.
Definition at line 230 of file fe_series.h.
|
private |
Transformation matrices for each FiniteElement.
Definition at line 235 of file fe_series.h.
|
private |
Auxiliary vector to store unrolled coefficients.
Definition at line 240 of file fe_series.h.