v0.13.2
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
OpElasticTools::OpCalculateNeoHookeStressTangent Struct Reference

#include <users_modules/multifield_plasticity/src/ElasticOperators.hpp>

Inheritance diagram for OpElasticTools::OpCalculateNeoHookeStressTangent:
[legend]
Collaboration diagram for OpElasticTools::OpCalculateNeoHookeStressTangent:
[legend]

Public Member Functions

 OpCalculateNeoHookeStressTangent (const std::string field_name, boost::shared_ptr< CommonData > common_data_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Private Attributes

boost::shared_ptr< CommonDatacommonDataPtr
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 878 of file ElasticOperators.hpp.

Constructor & Destructor Documentation

◆ OpCalculateNeoHookeStressTangent()

OpElasticTools::OpCalculateNeoHookeStressTangent::OpCalculateNeoHookeStressTangent ( const std::string  field_name,
boost::shared_ptr< CommonData common_data_ptr 
)

Definition at line 319 of file ElasticOperators.cpp.

322 commonDataPtr(common_data_ptr) {
323 // Operator is only executed for vertices
324 std::fill(&doEntities[MBEDGE], &doEntities[MBMAXTYPE], false);
325}
constexpr auto field_name
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
@ OPROW
operator doWork function is executed on FE rows

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpElasticTools::OpCalculateNeoHookeStressTangent::doWork ( int  side,
EntityType  type,
EntData data 
)
virtual

Operator for linear form, usually to calculate values on right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 326 of file ElasticOperators.cpp.

328 {
330
331 const size_t nb_gauss_pts = data.getN().size1();
332 auto t_stress = getFTensor2SymmetricFromMat<3>(*(commonDataPtr->mStressPtr));
333 commonDataPtr->mPiolaStressPtr->resize(9, nb_gauss_pts, false);
334 auto t_piola = getFTensor2FromMat<3, 3>(*(commonDataPtr->mPiolaStressPtr));
335 auto F = getFTensor2FromMat<3, 3>(*(commonDataPtr->mDeformationGradient));
336 auto t_detF = getFTensor0FromVec(*(commonDataPtr->detF));
337
338 MatrixDouble &dP_dF = *(commonDataPtr->materialTangent);
339 dP_dF.resize(81, nb_gauss_pts, false);
340 auto t_dP_dF = getFTensor4FromMat<3, 3, 3, 3>(dP_dF);
341
342 MatrixDouble &dE_dF = *(commonDataPtr->dE_dF_mat);
343 dE_dF.resize(81, nb_gauss_pts, false);
344 auto t_dE_dF = getFTensor4FromMat<3, 3, 3, 3>(dE_dF);
345 const bool is_lhs = (int)getTSCtx() == 3;
346 constexpr auto t_kd = Kronecker_Delta<int>();
347
348 const double lambda = LAMBDA((*cache).young_modulus, (*cache).poisson);
349 const double mu = MU((*cache).young_modulus, (*cache).poisson);
350
351 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
352
354 CHKERR invertTensor3by3(F, t_detF, invF);
355 const double log_det = log(t_detF);
356
357 if (is_lhs) {
358 // compute tangent
359 const double coef = mu - lambda * log_det;
360 t_dP_dF(i, J, k, L) =
361 invF(J, i) * (invF(L, k) * lambda) + invF(L, i) * (invF(J, k) * coef);
362
363 t_dP_dF(0, 0, 0, 0) += mu;
364 t_dP_dF(0, 1, 0, 1) += mu;
365 t_dP_dF(0, 2, 0, 2) += mu;
366 t_dP_dF(1, 0, 1, 0) += mu;
367 t_dP_dF(1, 1, 1, 1) += mu;
368 t_dP_dF(1, 2, 1, 2) += mu;
369 t_dP_dF(2, 0, 2, 0) += mu;
370 t_dP_dF(2, 1, 2, 1) += mu;
371 t_dP_dF(2, 2, 2, 2) += mu;
372
373 } else {
374
375 // compute only stress
376 const double var = lambda * log_det - mu;
377 t_piola(i, J) = mu * F(i, J) + var * invF(J, i);
378 }
379
380 ++t_detF;
381 ++t_piola;
382 ++t_dP_dF;
383 ++t_stress;
384 ++F;
385 }
386
388}
static Index< 'J', 3 > J
Kronecker Delta class.
#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
#define CHKERR
Inline error check.
Definition: definitions.h:535
#define MU(E, NU)
Definition: fem_tools.h:23
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
@ F
constexpr double lambda
surface tension
constexpr auto t_kd
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'k', 3 > k
MoFEMErrorCode invertTensor3by3(ublas::matrix< T, L, A > &jac_data, ublas::vector< T, A > &det_data, ublas::matrix< T, L, A > &inv_jac_data)
Calculate inverse of tensor rank 2 at integration points.
Definition: Templates.hpp:1363
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
constexpr double mu
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<CommonData> OpElasticTools::OpCalculateNeoHookeStressTangent::commonDataPtr
private

Definition at line 885 of file ElasticOperators.hpp.


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