![]() |
Reference documentation for deal.II version 8.5.1
|
#include <deal.II/fe/fe_series.h>
Public Member Functions | |
Legendre (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, double > &legendre_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 | |
const unsigned int | N |
SmartPointer< const hp::FECollection< dim > > | fe_collection |
SmartPointer< const hp::QCollection< dim > > | q_collection |
std::vector< FullMatrix< double > > | legendre_transform_matrices |
std::vector< 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 series of Legendre functions on a reference element.
Legendre functions are solutions to Legendre's differential equation
and can be expressed using Rodrigues' formula
These polynomials are orthogonal with respect to the inner product on the interval
and are complete. A family of -orthogonal polynomials on
can be constructed via
An arbitrary scalar FE field on the reference element can be expanded in the complete orthogonal basis as
From the orthogonality property of the basis, it follows that
This class calculates coefficients using
-dimensional Legendre polynomials constructed from
using tensor product rule.
Definition at line 188 of file fe_series.h.
FESeries::Legendre< dim >::Legendre | ( | 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 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 245 of file fe_series.cc.
void FESeries::Legendre< dim >::calculate | ( | const ::Vector< double > & | local_dof_values, |
const unsigned int | cell_active_fe_index, | ||
Table< dim, double > & | 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 259 of file fe_series.cc.
|
private |
Ensure that the transformation matrix for FiniteElement index fe_index
is calculated. If not, calculate it.
|
private |
Number of coefficients in each direction
Definition at line 215 of file fe_series.h.
|
private |
hp::FECollection for which transformation matrices will be calculated.
Definition at line 220 of file fe_series.h.
|
private |
hp::QCollection used in calculation of transformation matrices.
Definition at line 225 of file fe_series.h.
|
private |
Transformation matrices for each FiniteElement.
Definition at line 236 of file fe_series.h.
|
private |
Auxiliary vector to store unrolled coefficients.
Definition at line 241 of file fe_series.h.