v0.15.0
Loading...
Searching...
No Matches
EshelbianPlasticity::CGGUserPolynomialBase Struct Reference

#include "users_modules/eshelbian_plasticity/src/EshelbianPlasticity.hpp"

Inheritance diagram for EshelbianPlasticity::CGGUserPolynomialBase:
[legend]
Collaboration diagram for EshelbianPlasticity::CGGUserPolynomialBase:
[legend]

Public Types

using CachePhi = boost::tuple<int, int, MatrixDouble>
 

Public Member Functions

 CGGUserPolynomialBase (boost::shared_ptr< CachePhi > cache_phi=nullptr)
 
 ~CGGUserPolynomialBase ()=default
 
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, BaseFunctionUnknownInterface **iface) const
 
MoFEMErrorCode getValue (MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
 

Private Member Functions

MoFEMErrorCode getValueHdivForCGGBubble (MatrixDouble &pts)
 

Private Attributes

MatrixDouble shapeFun
 
boost::shared_ptr< CachePhicachePhiPtr
 

Detailed Description

Definition at line 37 of file EshelbianPlasticity.hpp.

Member Typedef Documentation

◆ CachePhi

using EshelbianPlasticity::CGGUserPolynomialBase::CachePhi = boost::tuple<int, int, MatrixDouble>

Definition at line 39 of file EshelbianPlasticity.hpp.

Constructor & Destructor Documentation

◆ CGGUserPolynomialBase()

EshelbianPlasticity::CGGUserPolynomialBase::CGGUserPolynomialBase ( boost::shared_ptr< CachePhi > cache_phi = nullptr)

Definition at line 2197 of file EshelbianPlasticity.cpp.

2199 : TetPolynomialBase(), cachePhiPtr(cache_phi_otr) {}
Calculate base functions on tetrahedral.

◆ ~CGGUserPolynomialBase()

EshelbianPlasticity::CGGUserPolynomialBase::~CGGUserPolynomialBase ( )
default

Member Function Documentation

◆ getValue()

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

Definition at line 2209 of file EshelbianPlasticity.cpp.

2210 {
2212
2213 this->cTx = ctx_ptr->getInterface<EntPolynomialBaseCtx>();
2214
2215 int nb_gauss_pts = pts.size2();
2216 if (!nb_gauss_pts) {
2218 }
2219
2220 if (pts.size1() < 3) {
2221 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2222 "Wrong dimension of pts, should be at least 3 rows with "
2223 "coordinates");
2224 }
2225
2226 const auto base = this->cTx->bAse;
2227 EntitiesFieldData &data = this->cTx->dAta;
2228
2229 switch (this->cTx->sPace) {
2230 case HDIV:
2232 break;
2233 case L2:
2234 data.dataOnEntities[MBVERTEX][0].getN(base).resize(nb_gauss_pts, 4, false);
2236 &*data.dataOnEntities[MBVERTEX][0].getN(base).data().begin(),
2237 &pts(0, 0), &pts(1, 0), &pts(2, 0), nb_gauss_pts);
2238 data.dataOnEntities[MBVERTEX][0].getDiffN(base).resize(4, 3, false);
2239 std::copy(Tools::diffShapeFunMBTET.begin(), Tools::diffShapeFunMBTET.end(),
2240 data.dataOnEntities[MBVERTEX][0].getDiffN(base).data().begin());
2241 this->cTx->basePolynomialsType0 = Legendre_polynomials;
2242 CHKERR getValueL2AinsworthBase(pts);
2243 break;
2244 default:
2245 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Not yet implemented");
2246 }
2247
2249}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
Definition definitions.h:88
@ 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.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PetscErrorCode Legendre_polynomials(int p, double s, double *diff_s, double *L, double *diffL, const int dim)
Calculate Legendre approximation basis.
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
Class used to pass element data to calculate base functions on tet,triangle,edge.
data structure for finite element entity
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
static MoFEMErrorCode shapeFunMBTET(double *shape, const double *ksi, const double *eta, const double *zeta, const double nb)
Calculate shape functions on tetrahedron.
Definition Tools.hpp:747
static constexpr std::array< double, 12 > diffShapeFunMBTET
Definition Tools.hpp:271

◆ getValueHdivForCGGBubble()

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

number of integration points

Definition at line 2252 of file EshelbianPlasticity.cpp.

2252 {
2254
2255 // This should be used only in case USER_BASE is selected
2256 if (this->cTx->bAse != USER_BASE) {
2257 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2258 "Wrong base, should be USER_BASE");
2259 }
2260 // get access to data structures on element
2261 EntitiesFieldData &data = this->cTx->dAta;
2262 // Get approximation order on element. Note that bubble functions are only
2263 // on tetrahedron.
2264 const int order = data.dataOnEntities[MBTET][0].getOrder();
2265 /// number of integration points
2266 const int nb_gauss_pts = pts.size2();
2267
2268 auto check_cache = [this](int order, int nb_gauss_pts) -> bool {
2269 if (cachePhiPtr) {
2270 return cachePhiPtr->get<0>() == order &&
2271 cachePhiPtr->get<1>() == nb_gauss_pts;
2272 } else {
2273 return false;
2274 }
2275 };
2276
2277 if (check_cache(order, nb_gauss_pts)) {
2278 auto &mat = cachePhiPtr->get<2>();
2279 auto &phi = data.dataOnEntities[MBTET][0].getN(USER_BASE);
2280 phi.resize(mat.size1(), mat.size2(), false);
2281 noalias(phi) = mat;
2282
2283 } else {
2284 // calculate shape functions, i.e. barycentric coordinates
2285 shapeFun.resize(nb_gauss_pts, 4, false);
2286 CHKERR ShapeMBTET(&*shapeFun.data().begin(), &pts(0, 0), &pts(1, 0),
2287 &pts(2, 0), nb_gauss_pts);
2288 // derivatives of shape functions
2289 double diff_shape_fun[12];
2290 CHKERR ShapeDiffMBTET(diff_shape_fun);
2291
2292 const int nb_base_functions = NBVOLUMETET_CCG_BUBBLE(order);
2293 // get base functions and set size
2294 MatrixDouble &phi = data.dataOnEntities[MBTET][0].getN(USER_BASE);
2295 phi.resize(nb_gauss_pts, 9 * nb_base_functions, false);
2296 // finally calculate base functions
2298 &phi(0, 0), &phi(0, 1), &phi(0, 2),
2299
2300 &phi(0, 3), &phi(0, 4), &phi(0, 5),
2301
2302 &phi(0, 6), &phi(0, 7), &phi(0, 8));
2303 CHKERR CGG_BubbleBase_MBTET(order, &shapeFun(0, 0), diff_shape_fun, t_phi,
2304 nb_gauss_pts);
2305
2306 if (cachePhiPtr) {
2307 cachePhiPtr->get<0>() = order;
2308 cachePhiPtr->get<1>() = nb_gauss_pts;
2309 cachePhiPtr->get<2>().resize(phi.size1(), phi.size2(), false);
2310 noalias(cachePhiPtr->get<2>()) = phi;
2311 }
2312 }
2313
2315}
#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 ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
constexpr int order
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
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.
static double phi

◆ query_interface()

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

Definition at line 2201 of file EshelbianPlasticity.cpp.

2203 {
2204 *iface = const_cast<CGGUserPolynomialBase *>(this);
2205 return 0;
2206}
CGGUserPolynomialBase(boost::shared_ptr< CachePhi > cache_phi=nullptr)

Member Data Documentation

◆ cachePhiPtr

boost::shared_ptr<CachePhi> EshelbianPlasticity::CGGUserPolynomialBase::cachePhiPtr
private

Definition at line 51 of file EshelbianPlasticity.hpp.

◆ shapeFun

MatrixDouble EshelbianPlasticity::CGGUserPolynomialBase::shapeFun
private

Definition at line 50 of file EshelbianPlasticity.hpp.


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