v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
EshelbianPlasticity::CGGUserPolynomialBase Struct Reference
Inheritance diagram for EshelbianPlasticity::CGGUserPolynomialBase:
[legend]
Collaboration diagram for EshelbianPlasticity::CGGUserPolynomialBase:
[legend]

Public Member Functions

 CGGUserPolynomialBase ()
 
 ~CGGUserPolynomialBase ()
 
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, BaseFunctionUnknownInterface **iface) const
 
MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- Public Member Functions inherited from MoFEM::BaseFunction
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, MoFEM::UnknownInterface **iface) const
 
virtual MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
virtual MoFEMErrorCode getValue (MatrixDouble &pts_x, MatrixDouble &pts_t, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 
- Public Member Functions inherited from MoFEM::BaseFunctionUnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
virtual ~BaseFunctionUnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register 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 ()=default
 

Private Member Functions

MoFEMErrorCode getValueHdivForCGGBubble (MatrixDouble &pts)
 

Private Attributes

EntPolynomialBaseCtxcTx
 
MatrixDouble shapeFun
 

Additional Inherited Members

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

Detailed Description

Definition at line 483 of file EshelbianPlasticity.cpp.

Constructor & Destructor Documentation

◆ CGGUserPolynomialBase()

EshelbianPlasticity::CGGUserPolynomialBase::CGGUserPolynomialBase ( )
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 485 of file EshelbianPlasticity.cpp.

485{}

◆ ~CGGUserPolynomialBase()

EshelbianPlasticity::CGGUserPolynomialBase::~CGGUserPolynomialBase ( )
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 486 of file EshelbianPlasticity.cpp.

486{}

Member Function Documentation

◆ getValue()

MoFEMErrorCode EshelbianPlasticity::CGGUserPolynomialBase::getValue ( MatrixDouble &  pts,
boost::shared_ptr< BaseFunctionCtx ctx_ptr 
)
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 494 of file EshelbianPlasticity.cpp.

495 {
497
499
500 int nb_gauss_pts = pts.size2();
501 if (!nb_gauss_pts) {
503 }
504
505 if (pts.size1() < 3) {
506 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
507 "Wrong dimension of pts, should be at least 3 rows with "
508 "coordinates");
509 }
510
511 switch (cTx->sPace) {
512 case HDIV:
514 break;
515 default:
516 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
517 }
518
520 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
Class used to pass element data to calculate base functions on tet,triangle,edge.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ getValueHdivForCGGBubble()

MoFEMErrorCode EshelbianPlasticity::CGGUserPolynomialBase::getValueHdivForCGGBubble ( MatrixDouble &  pts)
inlineprivate

number of integration points

Examples
EshelbianPlasticity.cpp.

Definition at line 527 of file EshelbianPlasticity.cpp.

527 {
529
530 // This should be used only in case USER_BASE is selected
531 if (cTx->bAse != USER_BASE) {
532 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
533 "Wrong base, should be USER_BASE");
534 }
535 // get access to data structures on element
536 EntitiesFieldData &data = cTx->dAta;
537 // Get approximation order on element. Note that bubble functions are only
538 // on tetrahedron.
539 const int order = data.dataOnEntities[MBTET][0].getOrder();
540 /// number of integration points
541 const int nb_gauss_pts = pts.size2();
542 // number of base functions
543
544 // calculate shape functions, i.e. barycentric coordinates
545 shapeFun.resize(nb_gauss_pts, 4, false);
546 CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
547 &pts(2, 0), nb_gauss_pts);
548 // direvatives of shape functions
549 double diff_shape_fun[12];
550 CHKERR ShapeDiffMBTET(diff_shape_fun);
551
552 const int nb_base_functions = NBVOLUMETET_CCG_BUBBLE(order);
553 // get base functions and set size
554 MatrixDouble &phi = data.dataOnEntities[MBTET][0].getN(USER_BASE);
555 phi.resize(nb_gauss_pts, 9 * nb_base_functions, false);
556 // finally calculate base functions
558 &phi(0, 0), &phi(0, 1), &phi(0, 2),
559
560 &phi(0, 3), &phi(0, 4), &phi(0, 5),
561
562 &phi(0, 6), &phi(0, 7), &phi(0, 8));
563 CHKERR CGG_BubbleBase_MBTET(order, &shapeFun(0, 0), diff_shape_fun, t_phi,
564 nb_gauss_pts);
565
567 }
#define NBVOLUMETET_CCG_BUBBLE(P)
Bubble function for CGG H div space.
@ USER_BASE
user implemented approximation base
Definition: definitions.h:68
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:319
PetscErrorCode ShapeMBTET(double *N, const double *G_X, const double *G_Y, const double *G_Z, int DIM)
calculate shape functions
Definition: fem_tools.c:306
static double phi
MoFEMErrorCode CGG_BubbleBase_MBTET(const int p, const double *N, const double *diffN, FTensor::Tensor2< FTensor::PackPtr< double *, 9 >, 3, 3 > &phi, const int gdim)
Calculate CGGT tonsorial bubble base.
const FieldApproximationBase bAse
data structure for finite element entity
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities

◆ query_interface()

MoFEMErrorCode EshelbianPlasticity::CGGUserPolynomialBase::query_interface ( boost::typeindex::type_index  type_index,
BaseFunctionUnknownInterface **  iface 
) const
inline
Examples
EshelbianPlasticity.cpp.

Definition at line 488 of file EshelbianPlasticity.cpp.

489 {
490 *iface = const_cast<CGGUserPolynomialBase *>(this);
491 return 0;
492 }

Member Data Documentation

◆ cTx

EntPolynomialBaseCtx* EshelbianPlasticity::CGGUserPolynomialBase::cTx
private
Examples
EshelbianPlasticity.cpp.

Definition at line 523 of file EshelbianPlasticity.cpp.

◆ shapeFun

MatrixDouble EshelbianPlasticity::CGGUserPolynomialBase::shapeFun
private
Examples
EshelbianPlasticity.cpp.

Definition at line 525 of file EshelbianPlasticity.cpp.


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