v0.8.4
Public Member Functions | List of all members
MoFEM::LegendrePolynomial Struct Reference

Calculating Legendre base functions. More...

#include <src/approximation/LegendrePolynomial.hpp>

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

Public Member Functions

MoFEMErrorCode query_interface (const MOFEMuuid &uuid, MoFEM::UnknownInterface **iface) const
 
 LegendrePolynomial ()
 
 ~LegendrePolynomial ()
 
MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- Public Member Functions inherited from MoFEM::BaseFunction
 BaseFunction ()
 
 ~BaseFunction ()
 
virtual MoFEMErrorCode getValue (MatrixDouble &pts_x, MatrixDouble &pts_t, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 

Additional Inherited Members

- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

Calculating Legendre base functions.

Definition at line 67 of file LegendrePolynomial.hpp.

Constructor & Destructor Documentation

◆ LegendrePolynomial()

MoFEM::LegendrePolynomial::LegendrePolynomial ( )

Definition at line 71 of file LegendrePolynomial.hpp.

71 {}

◆ ~LegendrePolynomial()

MoFEM::LegendrePolynomial::~LegendrePolynomial ( )

Definition at line 72 of file LegendrePolynomial.hpp.

72 {}

Member Function Documentation

◆ getValue()

MoFEMErrorCode MoFEM::LegendrePolynomial::getValue ( MatrixDouble pts,
boost::shared_ptr< BaseFunctionCtx ctx_ptr 
)
virtual

Reimplemented from MoFEM::BaseFunction.

Reimplemented in MoFEM::KernelLobattoPolynomial, and MoFEM::LobattoPolynomial.

Definition at line 65 of file LegendrePolynomial.cpp.

68  {
69 
72  ierr = ctx_ptr->query_interface(IDD_LEGENDRE_BASE_FUNCTION,&iface); CHKERRG(ierr);
73  LegendrePolynomialCtx *ctx = reinterpret_cast<LegendrePolynomialCtx*>(iface);
74  ctx->baseFunPtr->resize(pts.size2(),ctx->P+1,false);
75  ctx->baseDiffFunPtr->resize(pts.size2(),ctx->dIm*(ctx->P+1),false);
76  double *l = NULL;
77  double *diff_l = NULL;
78  for(unsigned int gg = 0;gg<pts.size2();gg++) {
79  if(ctx->baseFunPtr) l = &((*ctx->baseFunPtr)(gg,0));
80  if(ctx->baseDiffFunPtr) diff_l = &((*ctx->baseDiffFunPtr)(gg,0));
81  ierr = (ctx->basePolynomialsType0)(ctx->P,pts(0,gg),ctx->diffS,l,diff_l,ctx->dIm); CHKERRG(ierr);
82  }
84 }
boost::shared_ptr< MatrixDouble > baseFunPtr
boost::shared_ptr< MatrixDouble > baseDiffFunPtr
PetscErrorCode(* basePolynomialsType0)(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
base class for all interface classes
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
static const MOFEMuuid IDD_LEGENDRE_BASE_FUNCTION
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
Class used to give arguments to Legendre base functions.

◆ query_interface()

MoFEMErrorCode MoFEM::LegendrePolynomial::query_interface ( const MOFEMuuid uuid,
MoFEM::UnknownInterface **  iface 
) const
virtual

Reimplemented from MoFEM::BaseFunction.

Reimplemented in MoFEM::KernelLobattoPolynomial, and MoFEM::LobattoPolynomial.

Definition at line 49 of file LegendrePolynomial.cpp.

51  {
52 
54  *iface = NULL;
55  if(uuid == IDD_LEGENDRE_BASE_FUNCTION) {
56  *iface = const_cast<LegendrePolynomial*>(this);
58  } else {
59  SETERRQ(PETSC_COMM_WORLD,MOFEM_DATA_INCONSISTENCY,"wrong interference");
60  }
63 }
Calculating Legendre base functions.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
static const MOFEMuuid IDD_LEGENDRE_BASE_FUNCTION
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
MoFEMErrorCode query_interface(const MOFEMuuid &uuid, MoFEM::UnknownInterface **iface) const

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