Reference documentation for deal.II version 8.5.0
Public Member Functions | Private Member Functions | Private Attributes | List of all members
FESeries::Legendre< dim > Class Template Reference

#include <deal.II/fe/fe_series.h>

Inheritance diagram for FESeries::Legendre< dim >:

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)
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 Subscriptor (const Subscriptor &)
 Subscriptor (Subscriptor &&)
virtual ~Subscriptor ()
Subscriptoroperator= (const Subscriptor &)
Subscriptoroperator= (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 Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, char *arg2, std::string &arg3)
static::ExceptionBase & ExcNoSubscriber (char *arg1, char *arg2)

Detailed Description

template<int dim>
class FESeries::Legendre< dim >

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.

Denis Davydov, 2016.

Definition at line 188 of file fe_series.h.

Constructor & Destructor Documentation

template<int dim>
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

Member Function Documentation

template<int dim>
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

template<int dim>
void FESeries::Legendre< dim >::ensure_existence ( const unsigned int  fe_index)

Ensure that the transformation matrix for FiniteElement index fe_index is calculated. If not, calculate it.

Member Data Documentation

template<int dim>
const unsigned int FESeries::Legendre< dim >::N

Number of coefficients in each direction

Definition at line 215 of file fe_series.h.

template<int dim>
SmartPointer<const hp::FECollection<dim> > FESeries::Legendre< dim >::fe_collection

hp::FECollection for which transformation matrices will be calculated.

Definition at line 220 of file fe_series.h.

template<int dim>
SmartPointer<const hp::QCollection<dim> > FESeries::Legendre< dim >::q_collection

hp::QCollection used in calculation of transformation matrices.

Definition at line 225 of file fe_series.h.

template<int dim>
std::vector<FullMatrix<double> > FESeries::Legendre< dim >::legendre_transform_matrices

Transformation matrices for each FiniteElement.

Definition at line 236 of file fe_series.h.

template<int dim>
std::vector<double> FESeries::Legendre< dim >::unrolled_coefficients

Auxiliary vector to store unrolled coefficients.

Definition at line 241 of file fe_series.h.

The documentation for this class was generated from the following files: