DUBeat  1.0.1
High-order discontinuous Galerkin methods and applications to cardiac electrophysiology
DUBValues< dim > Class Template Reference

Class representing the Dubiner basis functions definitions. More...

#include <DUBValues.hpp>

Public Member Functions

 DUBValues (const unsigned int degree)
 Constructor. More...
 
 DUBValues (DUBValues< dim > &DUBValues)=default
 Default copy constructor. More...
 
 DUBValues (const DUBValues< dim > &DUBValues)=default
 Default const copy constructor. More...
 
 DUBValues (DUBValues< dim > &&DUBValues)=default
 Default move constructor. More...
 
unsigned int get_dofs_per_cell () const
 Return the number of degrees of freedom per element. More...
 
unsigned int get_dofs_per_cell (const unsigned int degree) const
 Same as get_dofs_per_cell() but with variable space degree. More...
 
double shape_value (const unsigned int function_no, const dealii::Point< dim > &quadrature_point) const
 Evaluation of the Dubiner basis functions. More...
 
dealii::Tensor< 1, dim > shape_grad (const unsigned int function_no, const dealii::Point< dim > &quadrature_point) const
 Evaluation of the Dubiner basis functions gradients. More...
 
virtual ~DUBValues ()=default
 Default destructor. More...
 
double shape_value (const unsigned int function_no, const dealii::Point< 2 > &quadrature_point) const
 Evaluation of the function basis in two dimensions. More...
 
double shape_value (const unsigned int function_no, const dealii::Point< 3 > &quadrature_point) const
 Evaluation of the function basis in three dimensions. More...
 
dealii::Tensor< 1, 2 > shape_grad (const unsigned int function_no, const dealii::Point< 2 > &quadrature_point) const
 Evaluation of the gradient of the function basis in two dimensions. More...
 
dealii::Tensor< 1, 3 > shape_grad (const unsigned int function_no, const dealii::Point< 3 > &quadrature_point) const
 Evaluation of the gradient of the function basis in three dimensions. More...
 

Public Attributes

unsigned int dofs_per_cell
 Number of degrees of freedom. More...
 

Protected Member Functions

double eval_jacobi_polynomial (const unsigned int n, const unsigned int alpha, const unsigned int beta, const double eval_point) const
 Evaluation of the \((n,\alpha,\beta)\)-Jacobi polynomial (used to evaluate Dubiner basis). More...
 
std::array< unsigned int, dim > fun_coeff_conversion (const unsigned int n) const
 Conversion from the FEM basis taxonomy ( \(1^{st}\) basis function, \(2^{nd}\) basis function \(\ldots\)) to the Dubiner basis taxonomy (where every basis function is indexed by a tuple of indexes \((i,j)\) in 2D and three indexes \((i,j,k)\) in 3D). More...
 
std::array< unsigned int, 2 > fun_coeff_conversion (const unsigned int n) const
 Specialization in two dimensions, i.e. More...
 
std::array< unsigned int, 3 > fun_coeff_conversion (const unsigned int n) const
 Specialization in three dimensions, i.e. More...
 

Protected Attributes

const unsigned int poly_degree
 Polynomial degree used. More...
 
const double tol = 1e-10
 Default tolerance. More...
 

Detailed Description

template<unsigned int dim>
class DUBValues< dim >

Class representing the Dubiner basis functions definitions.

The main utility is the evaluations of these functions and their gradients on the reference cell. Instead, for their evaluations on the real cells see VolumeHandlerDG and FaceHandlerDG.

Definition at line 49 of file DUBValues.hpp.

Constructor & Destructor Documentation

◆ DUBValues() [1/4]

template<unsigned int dim>
DUBValues< dim >::DUBValues ( const unsigned int  degree)
inline

Constructor.

Definition at line 72 of file DUBValues.hpp.

◆ DUBValues() [2/4]

template<unsigned int dim>
DUBValues< dim >::DUBValues ( DUBValues< dim > &  DUBValues)
default

Default copy constructor.

◆ DUBValues() [3/4]

template<unsigned int dim>
DUBValues< dim >::DUBValues ( const DUBValues< dim > &  DUBValues)
default

Default const copy constructor.

◆ DUBValues() [4/4]

template<unsigned int dim>
DUBValues< dim >::DUBValues ( DUBValues< dim > &&  DUBValues)
default

Default move constructor.

◆ ~DUBValues()

template<unsigned int dim>
virtual DUBValues< dim >::~DUBValues ( )
virtualdefault

Default destructor.

Member Function Documentation

◆ eval_jacobi_polynomial()

template<unsigned int dim>
double DUBValues< dim >::eval_jacobi_polynomial ( const unsigned int  n,
const unsigned int  alpha,
const unsigned int  beta,
const double  eval_point 
) const
protected

Evaluation of the \((n,\alpha,\beta)\)-Jacobi polynomial (used to evaluate Dubiner basis).

Definition at line 176 of file DUBValues.hpp.

◆ fun_coeff_conversion() [1/3]

template<unsigned int dim>
std::array<unsigned int, dim> DUBValues< dim >::fun_coeff_conversion ( const unsigned int  n) const
protected

Conversion from the FEM basis taxonomy ( \(1^{st}\) basis function, \(2^{nd}\) basis function \(\ldots\)) to the Dubiner basis taxonomy (where every basis function is indexed by a tuple of indexes \((i,j)\) in 2D and three indexes \((i,j,k)\) in 3D).

It is needed to order sequentially all the Dubiner basis functions of a certain element.

◆ fun_coeff_conversion() [2/3]

std::array< unsigned int, 2 > DUBValues< 2 >::fun_coeff_conversion ( const unsigned int  n) const
protected

Specialization in two dimensions, i.e.

from an integer \(n\) returns a tuple of indexes \((i,j)\) such that the \(n\)-th Dubiner function corresponds to the function with indexes \((i,j)\).

Definition at line 121 of file DUBValues.hpp.

◆ fun_coeff_conversion() [3/3]

std::array< unsigned int, 3 > DUBValues< 3 >::fun_coeff_conversion ( const unsigned int  n) const
protected

Specialization in three dimensions, i.e.

from an integer \(n\) returns three indexes \((i,j,k)\) such that the \(n\)-th Dubiner function corresponds to the function with indexes \((i,j,k)\).

Definition at line 143 of file DUBValues.hpp.

◆ get_dofs_per_cell() [1/2]

template<unsigned int dim>
unsigned int DUBValues< dim >::get_dofs_per_cell

Return the number of degrees of freedom per element.

Definition at line 223 of file DUBValues.hpp.

◆ get_dofs_per_cell() [2/2]

template<unsigned int dim>
unsigned int DUBValues< dim >::get_dofs_per_cell ( const unsigned int  degree) const

Same as get_dofs_per_cell() but with variable space degree.

Definition at line 243 of file DUBValues.hpp.

◆ shape_grad() [1/3]

dealii::Tensor< 1, 2 > DUBValues< 2 >::shape_grad ( const unsigned int  function_no,
const dealii::Point< 2 > &  quadrature_point 
) const

Evaluation of the gradient of the function basis in two dimensions.

Definition at line 341 of file DUBValues.hpp.

◆ shape_grad() [2/3]

dealii::Tensor< 1, 3 > DUBValues< 3 >::shape_grad ( const unsigned int  function_no,
const dealii::Point< 3 > &  quadrature_point 
) const

Evaluation of the gradient of the function basis in three dimensions.

Definition at line 411 of file DUBValues.hpp.

◆ shape_grad() [3/3]

template<unsigned int dim>
dealii::Tensor<1, dim> DUBValues< dim >::shape_grad ( const unsigned int  function_no,
const dealii::Point< dim > &  quadrature_point 
) const

Evaluation of the Dubiner basis functions gradients.

◆ shape_value() [1/3]

double DUBValues< 2 >::shape_value ( const unsigned int  function_no,
const dealii::Point< 2 > &  quadrature_point 
) const

Evaluation of the function basis in two dimensions.

Definition at line 264 of file DUBValues.hpp.

◆ shape_value() [2/3]

double DUBValues< 3 >::shape_value ( const unsigned int  function_no,
const dealii::Point< 3 > &  quadrature_point 
) const

Evaluation of the function basis in three dimensions.

Definition at line 299 of file DUBValues.hpp.

◆ shape_value() [3/3]

template<unsigned int dim>
double DUBValues< dim >::shape_value ( const unsigned int  function_no,
const dealii::Point< dim > &  quadrature_point 
) const

Evaluation of the Dubiner basis functions.

Member Data Documentation

◆ dofs_per_cell

template<unsigned int dim>
unsigned int DUBValues< dim >::dofs_per_cell

Number of degrees of freedom.

Definition at line 83 of file DUBValues.hpp.

◆ poly_degree

template<unsigned int dim>
const unsigned int DUBValues< dim >::poly_degree
protected

Polynomial degree used.

Definition at line 53 of file DUBValues.hpp.

◆ tol

template<unsigned int dim>
const double DUBValues< dim >::tol = 1e-10
protected

Default tolerance.

Definition at line 56 of file DUBValues.hpp.


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