v0.13.0
Public Member Functions | Public Attributes | List of all members
NonlinearElasticElement::OpRhsPiolaKirchhoff Struct Reference

#include <users_modules/basic_finite_elements/src/NonLinearElasticElement.hpp>

Inheritance diagram for NonlinearElasticElement::OpRhsPiolaKirchhoff:
[legend]
Collaboration diagram for NonlinearElasticElement::OpRhsPiolaKirchhoff:
[legend]

Public Member Functions

 OpRhsPiolaKirchhoff (const std::string field_name, BlockData &data, CommonData &common_data)
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 
virtual MoFEMErrorCode aSemble (int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
double & getVolume ()
 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...
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () 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...
 
const std::string & getFEName () const
 Get name of the element. More...
 
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 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...
 
double & getMeasure ()
 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...
 

Public Attributes

BlockDatadAta
 
CommonDatacommonData
 
bool fieldDisp
 
bool aLe
 
ublas::vector< int > iNdices
 
VectorDouble nf
 
- 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...
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 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)>
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
cell_forces.cpp.

Definition at line 533 of file NonLinearElasticElement.hpp.

Constructor & Destructor Documentation

◆ OpRhsPiolaKirchhoff()

NonlinearElasticElement::OpRhsPiolaKirchhoff::OpRhsPiolaKirchhoff ( const std::string  field_name,
BlockData data,
CommonData common_data 
)

Definition at line 606 of file NonLinearElasticElement.cpp.

609  field_name, UserDataOperator::OPROW),
610  dAta(data), commonData(common_data), aLe(false) {}
ForcesAndSourcesCore::UserDataOperator UserDataOperator

Member Function Documentation

◆ aSemble()

MoFEMErrorCode NonlinearElasticElement::OpRhsPiolaKirchhoff::aSemble ( int  row_side,
EntityType  row_type,
EntitiesFieldData::EntData row_data 
)
virtual

Reimplemented in Smoother::OpRhsSmoother.

Definition at line 612 of file NonLinearElasticElement.cpp.

614  {
616 
617  int nb_dofs = row_data.getIndices().size();
618  int *indices_ptr = &row_data.getIndices()[0];
619  if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
620  iNdices.resize(nb_dofs, false);
621  noalias(iNdices) = row_data.getIndices();
622  indices_ptr = &iNdices[0];
623  VectorDofs &dofs = row_data.getFieldDofs();
624  VectorDofs::iterator dit = dofs.begin();
625  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
626  if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
628  iNdices[ii] = -1;
629  }
630  }
631  }
632  CHKERR VecSetValues(getFEMethod()->snes_f, nb_dofs, indices_ptr, &nf[0],
633  ADD_VALUES);
635 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
#define CHKERR
Inline error check.
Definition: definitions.h:548
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
const VectorInt & getIndices() const
Get global indices of dofs on entity.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

◆ doWork()

MoFEMErrorCode NonlinearElasticElement::OpRhsPiolaKirchhoff::doWork ( int  row_side,
EntityType  row_type,
EntitiesFieldData::EntData row_data 
)

Definition at line 637 of file NonLinearElasticElement.cpp.

639  {
641 
642  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
643  dAta.tEts.end()) {
645  }
646 
647  const int nb_dofs = row_data.getIndices().size();
648  if (nb_dofs == 0)
650  if ((unsigned int)nb_dofs > 3 * row_data.getN().size2()) {
651  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
652  }
653  const int nb_base_functions = row_data.getN().size2();
654  const int nb_gauss_pts = row_data.getN().size1();
655 
656  nf.resize(nb_dofs, false);
657  nf.clear();
658 
659  auto diff_base_functions = row_data.getFTensor1DiffN<3>();
662 
663  for (int gg = 0; gg != nb_gauss_pts; gg++) {
664  double val = getVolume() * getGaussPts()(3, gg);
665  MatrixDouble3by3 &stress = commonData.sTress[gg];
667  &stress(0, 0), &stress(0, 1), &stress(0, 2), &stress(1, 0),
668  &stress(1, 1), &stress(1, 2), &stress(2, 0), &stress(2, 1),
669  &stress(2, 2));
670  FTensor::Tensor1<double *, 3> rhs(&nf[0], &nf[1], &nf[2], 3);
671  int bb = 0;
672  for (; bb != nb_dofs / 3; bb++) {
673  rhs(i) += val * t3(i, j) * diff_base_functions(j);
674  ++rhs;
675  ++diff_base_functions;
676  }
677  for (; bb != nb_base_functions; bb++) {
678  ++diff_base_functions;
679  }
680  }
681 
682  CHKERR aSemble(row_side, row_type, row_data);
683 
685 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MatrixBoundedArray< double, 9 > MatrixDouble3by3
Definition: Types.hpp:116
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
Range tEts
constrains elements in block set
std::vector< MatrixDouble3by3 > sTress
virtual MoFEMErrorCode aSemble(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
FTensor::Index< 'j', 3 > j
FTensor::Index< 'i', 3 > i

Member Data Documentation

◆ aLe

bool NonlinearElasticElement::OpRhsPiolaKirchhoff::aLe

Definition at line 539 of file NonLinearElasticElement.hpp.

◆ commonData

CommonData& NonlinearElasticElement::OpRhsPiolaKirchhoff::commonData

Definition at line 537 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpRhsPiolaKirchhoff::dAta

Definition at line 536 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpRhsPiolaKirchhoff::fieldDisp

Definition at line 538 of file NonLinearElasticElement.hpp.

◆ iNdices

ublas::vector<int> NonlinearElasticElement::OpRhsPiolaKirchhoff::iNdices

Definition at line 541 of file NonLinearElasticElement.hpp.

◆ nf

VectorDouble NonlinearElasticElement::OpRhsPiolaKirchhoff::nf

Definition at line 545 of file NonLinearElasticElement.hpp.


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