![]() |
Reference documentation for deal.II version Git 74ec1a6c4b 2021-01-19 10:11:04 -0500
|
#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 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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 |
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 253 of file fe_series.h.
using FESeries::Legendre< dim, spacedim >::CoefficientType = double |
Definition at line 256 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 188 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 212 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 271 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 260 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 241 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 228 of file fe_series_legendre.cc.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 136 of file subscriptor.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 156 of file subscriptor.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 301 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 318 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 204 of file subscriptor.cc.
|
inlineinherited |
Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.
Definition at line 310 of file subscriptor.h.
|
private |
Number of coefficients in each direction for each finite element in the registered hp::FECollection.
Definition at line 346 of file fe_series.h.
|
private |
hp::FECollection for which transformation matrices will be calculated.
Definition at line 351 of file fe_series.h.
|
private |
hp::QCollection used in calculation of transformation matrices.
Definition at line 356 of file fe_series.h.
|
private |
Transformation matrices for each FiniteElement.
Definition at line 361 of file fe_series.h.
|
private |
Auxiliary vector to store unrolled coefficients.
Definition at line 366 of file fe_series.h.