v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MoFEM::LegendrePolynomialCtx Struct Reference

Context class for Legendre basis function arguments. More...

#include "src/approximation/LegendrePolynomial.hpp"

Inheritance diagram for MoFEM::LegendrePolynomialCtx:
[legend]
Collaboration diagram for MoFEM::LegendrePolynomialCtx:
[legend]

Public Member Functions

MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 LegendrePolynomialCtx (int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
 Constructor.
 
 ~LegendrePolynomialCtx ()
 
- Public Member Functions inherited from MoFEM::BaseFunctionUnknownInterface
virtual ~BaseFunctionUnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference to pointer of interface.
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface.
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface.
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface.
 
virtual ~UnknownInterface ()=default
 

Public Attributes

int P
 Approximation order.
 
doublediffS
 Derivatives of s parameter.
 
int dIm
 Dimension.
 
boost::shared_ptr< MatrixDoublebaseFunPtr
 Pointer to base function matrix.
 
boost::shared_ptr< MatrixDoublebaseDiffFunPtr
 Pointer to base function derivatives matrix.
 
PetscErrorCode(* basePolynomialsType0 )(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
 Function pointer for polynomial evaluation.
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version.
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version.
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version.
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version.
 

Detailed Description

Context class for Legendre basis function arguments.

This class stores parameters and data required for evaluating Legendre polynomial basis functions and their derivatives.

Definition at line 90 of file LegendrePolynomial.hpp.

Constructor & Destructor Documentation

◆ LegendrePolynomialCtx()

MoFEM::LegendrePolynomialCtx::LegendrePolynomialCtx ( int  p,
double diff_s,
int  dim,
boost::shared_ptr< MatrixDouble base_fun_ptr,
boost::shared_ptr< MatrixDouble base_diff_fun_ptr 
)
inline

Constructor.

Parameters
papproximation order
diff_sderivatives of s parameter
dimdimension
base_fun_ptrpointer to base function matrix
base_diff_fun_ptrpointer to base function derivatives matrix

Definition at line 114 of file LegendrePolynomial.hpp.

117 : P(p), diffS(diff_s), dIm(dim), baseFunPtr(base_fun_ptr),
118 baseDiffFunPtr(base_diff_fun_ptr),
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, int dim)
Calculate Legendre approximation basis.
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Function pointer for polynomial evaluation.
double * diffS
Derivatives of s parameter.
boost::shared_ptr< MatrixDouble > baseDiffFunPtr
Pointer to base function derivatives matrix.
boost::shared_ptr< MatrixDouble > baseFunPtr
Pointer to base function matrix.

◆ ~LegendrePolynomialCtx()

MoFEM::LegendrePolynomialCtx::~LegendrePolynomialCtx ( )
inline

Definition at line 120 of file LegendrePolynomial.hpp.

120{}

Member Function Documentation

◆ query_interface()

MoFEMErrorCode MoFEM::LegendrePolynomialCtx::query_interface ( boost::typeindex::type_index  type_index,
UnknownInterface **  iface 
) const
virtual

Reimplemented from MoFEM::BaseFunctionCtx.

Reimplemented in MoFEM::LobattoPolynomialCtx.

Definition at line 8 of file LegendrePolynomial.cpp.

9 {
10 *iface = const_cast<LegendrePolynomialCtx *>(this);
11 return 0;
12}
LegendrePolynomialCtx(int p, double *diff_s, int dim, boost::shared_ptr< MatrixDouble > base_fun_ptr, boost::shared_ptr< MatrixDouble > base_diff_fun_ptr)
Constructor.

Member Data Documentation

◆ baseDiffFunPtr

boost::shared_ptr<MatrixDouble> MoFEM::LegendrePolynomialCtx::baseDiffFunPtr

Pointer to base function derivatives matrix.

Definition at line 99 of file LegendrePolynomial.hpp.

◆ baseFunPtr

boost::shared_ptr<MatrixDouble> MoFEM::LegendrePolynomialCtx::baseFunPtr

Pointer to base function matrix.

Definition at line 98 of file LegendrePolynomial.hpp.

◆ basePolynomialsType0

PetscErrorCode(* MoFEM::LegendrePolynomialCtx::basePolynomialsType0) (int p, double s, double *diff_s, double *L, double *diffL, const int dim)

Function pointer for polynomial evaluation.

Definition at line 102 of file LegendrePolynomial.hpp.

◆ diffS

double* MoFEM::LegendrePolynomialCtx::diffS

Derivatives of s parameter.

Definition at line 96 of file LegendrePolynomial.hpp.

◆ dIm

int MoFEM::LegendrePolynomialCtx::dIm

Dimension.

Definition at line 97 of file LegendrePolynomial.hpp.

◆ P

int MoFEM::LegendrePolynomialCtx::P

Approximation order.

Definition at line 95 of file LegendrePolynomial.hpp.


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