|
| v0.14.0
|
#include <users_modules/eshelbian_plasticit/src/EshelbianPlasticity.hpp>
|
using | CachePhi = boost::tuple< int, int, MatrixDouble > |
|
Definition at line 21 of file EshelbianPlasticity.hpp.
◆ CachePhi
◆ CGGUserPolynomialBase()
EshelbianPlasticity::CGGUserPolynomialBase::CGGUserPolynomialBase |
( |
boost::shared_ptr< CachePhi > |
cache_phi = nullptr | ) |
|
◆ ~CGGUserPolynomialBase()
EshelbianPlasticity::CGGUserPolynomialBase::~CGGUserPolynomialBase |
( |
| ) |
|
|
default |
◆ getValue()
Definition at line 1380 of file EshelbianPlasticity.cpp.
1386 int nb_gauss_pts = pts.size2();
1387 if (!nb_gauss_pts) {
1391 if (pts.size1() < 3) {
1393 "Wrong dimension of pts, should be at least 3 rows with "
1397 const auto base = this->cTx->bAse;
1400 switch (this->cTx->sPace) {
1405 data.
dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 4,
false);
1406 CHKERR Tools::shapeFunMBTET(
1408 &pts(0, 0), &pts(1, 0), &pts(2, 0), nb_gauss_pts);
1409 data.
dataOnEntities[MBVERTEX][0].getDiffN(base).resize(4, 3,
false);
1410 std::copy(Tools::diffShapeFunMBTET.begin(), Tools::diffShapeFunMBTET.end(),
1413 CHKERR getValueL2AinsworthBase(pts);
◆ getValueHdivForCGGBubble()
MoFEMErrorCode EshelbianPlasticity::CGGUserPolynomialBase::getValueHdivForCGGBubble |
( |
MatrixDouble & |
pts | ) |
|
|
private |
number of integration points
Definition at line 1423 of file EshelbianPlasticity.cpp.
1429 "Wrong base, should be USER_BASE");
1437 const int nb_gauss_pts = pts.size2();
1439 auto check_cache = [
this](
int order,
int nb_gauss_pts) ->
bool {
1448 if (check_cache(
order, nb_gauss_pts)) {
1451 phi.resize(mat.size1(), mat.size2(),
false);
1456 shapeFun.resize(nb_gauss_pts, 4,
false);
1458 &pts(2, 0), nb_gauss_pts);
1460 double diff_shape_fun[12];
1466 phi.resize(nb_gauss_pts, 9 * nb_base_functions,
false);
◆ query_interface()
◆ cachePhiPtr
boost::shared_ptr<CachePhi> EshelbianPlasticity::CGGUserPolynomialBase::cachePhiPtr |
|
private |
◆ shapeFun
MatrixDouble EshelbianPlasticity::CGGUserPolynomialBase::shapeFun |
|
private |
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Class used to pass element data to calculate base functions on tet,triangle,edge.
@ L2
field with C-1 continuity
UBlasMatrix< double > MatrixDouble
@ USER_BASE
user implemented approximation base
#define CHKERR
Inline error check.
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
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.
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
boost::shared_ptr< CachePhi > cachePhiPtr
#define NBVOLUMETET_CCG_BUBBLE(P)
Bubble function for CGG H div space.
@ MOFEM_DATA_INCONSISTENCY
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
data structure for finite element entity
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
CGGUserPolynomialBase(boost::shared_ptr< CachePhi > cache_phi=nullptr)
Calculate base functions on tetrahedral.