v0.14.0
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 ~BaseFunctionUnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::UnknownInterface
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
 
boost::tuple< int, int, MatrixDouble > cachePhi = {0, 0, MatrixDouble()}
 

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 752 of file EshelbianPlasticity.cpp.

Constructor & Destructor Documentation

◆ CGGUserPolynomialBase()

EshelbianPlasticity::CGGUserPolynomialBase::CGGUserPolynomialBase ( )
inline

Definition at line 754 of file EshelbianPlasticity.cpp.

754 {}

◆ ~CGGUserPolynomialBase()

EshelbianPlasticity::CGGUserPolynomialBase::~CGGUserPolynomialBase ( )
inline

Definition at line 755 of file EshelbianPlasticity.cpp.

755 {}

Member Function Documentation

◆ getValue()

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

Definition at line 763 of file EshelbianPlasticity.cpp.

764  {
766 
767  cTx = ctx_ptr->getInterface<EntPolynomialBaseCtx>();
768 
769  int nb_gauss_pts = pts.size2();
770  if (!nb_gauss_pts) {
772  }
773 
774  if (pts.size1() < 3) {
775  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
776  "Wrong dimension of pts, should be at least 3 rows with "
777  "coordinates");
778  }
779 
780  switch (cTx->sPace) {
781  case HDIV:
783  break;
784  default:
785  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
786  }
787 
789  }

◆ getValueHdivForCGGBubble()

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

number of integration points

Definition at line 797 of file EshelbianPlasticity.cpp.

797  {
799 
800  // This should be used only in case USER_BASE is selected
801  if (cTx->bAse != USER_BASE) {
802  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
803  "Wrong base, should be USER_BASE");
804  }
805  // get access to data structures on element
806  EntitiesFieldData &data = cTx->dAta;
807  // Get approximation order on element. Note that bubble functions are only
808  // on tetrahedron.
809  const int order = data.dataOnEntities[MBTET][0].getOrder();
810  /// number of integration points
811  const int nb_gauss_pts = pts.size2();
812 
813  if(cachePhi.get<0>() == order && cachePhi.get<1>() == nb_gauss_pts) {
814  auto &mat = cachePhi.get<2>();
815  auto &phi = data.dataOnEntities[MBTET][0].getN(USER_BASE);
816  phi.resize(mat.size1(), mat.size2(), false);
817  noalias(phi) = mat;
818  } else {
819  // calculate shape functions, i.e. barycentric coordinates
820  shapeFun.resize(nb_gauss_pts, 4, false);
821  CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
822  &pts(2, 0), nb_gauss_pts);
823  // direvatives of shape functions
824  double diff_shape_fun[12];
825  CHKERR ShapeDiffMBTET(diff_shape_fun);
826 
827  const int nb_base_functions = NBVOLUMETET_CCG_BUBBLE(order);
828  // get base functions and set size
829  MatrixDouble &phi = data.dataOnEntities[MBTET][0].getN(USER_BASE);
830  phi.resize(nb_gauss_pts, 9 * nb_base_functions, false);
831  // finally calculate base functions
833  &phi(0, 0), &phi(0, 1), &phi(0, 2),
834 
835  &phi(0, 3), &phi(0, 4), &phi(0, 5),
836 
837  &phi(0, 6), &phi(0, 7), &phi(0, 8));
838  CHKERR CGG_BubbleBase_MBTET(order, &shapeFun(0, 0), diff_shape_fun, t_phi,
839  nb_gauss_pts);
840 
841  cachePhi.get<0>() = order;
842  cachePhi.get<1>() = nb_gauss_pts;
843  cachePhi.get<2>().resize(phi.size1(), phi.size2(), false);
844  noalias(cachePhi.get<2>()) = phi;
845  }
846 
848  }

◆ query_interface()

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

Definition at line 757 of file EshelbianPlasticity.cpp.

758  {
759  *iface = const_cast<CGGUserPolynomialBase *>(this);
760  return 0;
761  }

Member Data Documentation

◆ cachePhi

boost::tuple<int, int, MatrixDouble> EshelbianPlasticity::CGGUserPolynomialBase::cachePhi = {0, 0, MatrixDouble()}
private

Definition at line 795 of file EshelbianPlasticity.cpp.

◆ cTx

EntPolynomialBaseCtx* EshelbianPlasticity::CGGUserPolynomialBase::cTx
private

Definition at line 792 of file EshelbianPlasticity.cpp.

◆ shapeFun

MatrixDouble EshelbianPlasticity::CGGUserPolynomialBase::shapeFun
private

Definition at line 794 of file EshelbianPlasticity.cpp.


The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
EshelbianPlasticity::CGGUserPolynomialBase::CGGUserPolynomialBase
CGGUserPolynomialBase()
Definition: EshelbianPlasticity.cpp:754
MoFEM::EntPolynomialBaseCtx
Class used to pass element data to calculate base functions on tet,triangle,edge.
Definition: EntPolynomialBaseCtx.hpp:22
phi
static double phi
Definition: poisson_2d_dis_galerkin.cpp:29
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::EntPolynomialBaseCtx::dAta
EntitiesFieldData & dAta
Definition: EntPolynomialBaseCtx.hpp:36
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
EshelbianPlasticity::CGGUserPolynomialBase::shapeFun
MatrixDouble shapeFun
Definition: EshelbianPlasticity.cpp:794
USER_BASE
@ USER_BASE
user implemented approximation base
Definition: definitions.h:68
order
constexpr int order
Definition: dg_projection.cpp:18
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::EntPolynomialBaseCtx::sPace
const FieldSpace sPace
Definition: EntPolynomialBaseCtx.hpp:37
EshelbianPlasticity::CGGUserPolynomialBase::getValueHdivForCGGBubble
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
Definition: EshelbianPlasticity.cpp:797
EshelbianPlasticity::CGG_BubbleBase_MBTET
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.
Definition: CGGTonsorialBubbleBase.cpp:20
EshelbianPlasticity::CGGUserPolynomialBase::cTx
EntPolynomialBaseCtx * cTx
Definition: EshelbianPlasticity.cpp:792
MoFEM::EntPolynomialBaseCtx::bAse
const FieldApproximationBase bAse
Definition: EntPolynomialBaseCtx.hpp:38
NBVOLUMETET_CCG_BUBBLE
#define NBVOLUMETET_CCG_BUBBLE(P)
Bubble function for CGG H div space.
Definition: CGGTonsorialBubbleBase.hpp:19
ShapeDiffMBTET
PetscErrorCode ShapeDiffMBTET(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:319
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:56
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
HDIV
@ HDIV
field with continuous normal traction
Definition: definitions.h:87
MOFEM_NOT_IMPLEMENTED
@ MOFEM_NOT_IMPLEMENTED
Definition: definitions.h:32
EshelbianPlasticity::CGGUserPolynomialBase::cachePhi
boost::tuple< int, int, MatrixDouble > cachePhi
Definition: EshelbianPlasticity.cpp:795
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
ShapeMBTET
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