v0.9.1
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, DataForcesAndSourcesCore::EntData &row_data)
 
virtual MoFEMErrorCode aSemble (int row_side, EntityType row_type, DataForcesAndSourcesCore::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...
 
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...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHoGaussPtsJac ()
 
MatrixDoublegetHoGaussPtsInvJac ()
 
VectorDoublegetHoGaussPtsDetJac ()
 
auto getFTenosr0HoMeasure ()
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 Get coordinates at integration points taking geometry from field. More...
 
auto getFTensor2HoGaussPtsJac ()
 
auto getFTensor2HoGaussPtsInvJac ()
 
VolumeElementForcesAndSourcesCoreBasegetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
MoFEMErrorCode getDivergenceOfHDivBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, VectorDouble &div)
 Get divergence of base functions at integration point. More...
 
MoFEMErrorCode getCurlOfHCurlBaseFunctions (int side, EntityType type, DataForcesAndSourcesCore::EntData &data, int gg, MatrixDouble &curl)
 Get curl of base functions at integration point. More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, 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...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (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 ()
 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
 
DEPRECATED Vec getSnesF () const
 
DEPRECATED Vec getSnesX () const
 
DEPRECATED Mat getSnesA () const
 
DEPRECATED 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
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- 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, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &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
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, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCoreBase::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
cell_forces.cpp.

Definition at line 529 of file NonLinearElasticElement.hpp.

Constructor & Destructor Documentation

◆ OpRhsPiolaKirchhoff()

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

Definition at line 611 of file NonLinearElasticElement.cpp.

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

Member Function Documentation

◆ aSemble()

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

Reimplemented in Smoother::OpRhsSmoother.

Definition at line 617 of file NonLinearElasticElement.cpp.

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

◆ doWork()

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

Definition at line 642 of file NonLinearElasticElement.cpp.

644  {
646 
647  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
648  dAta.tEts.end()) {
650  }
651 
652  const int nb_dofs = row_data.getIndices().size();
653  if (nb_dofs == 0)
655  if ((unsigned int)nb_dofs > 3 * row_data.getN().size2()) {
656  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
657  }
658  const int nb_base_functions = row_data.getN().size2();
659  const int nb_gauss_pts = row_data.getN().size1();
660 
661  nf.resize(nb_dofs, false);
662  nf.clear();
663 
664  FTensor::Tensor1<double *, 3> diff_base_functions =
665  row_data.getFTensor1DiffN<3>();
668 
669  for (int gg = 0; gg != nb_gauss_pts; gg++) {
670  double val = getVolume() * getGaussPts()(3, gg);
671  if ((!aLe) && getHoGaussPtsDetJac().size() > 0) {
672  val *= getHoGaussPtsDetJac()[gg]; ///< higher order geometry
673  }
674  MatrixDouble3by3 &stress = commonData.sTress[gg];
676  &stress(0, 0), &stress(0, 1), &stress(0, 2), &stress(1, 0),
677  &stress(1, 1), &stress(1, 2), &stress(2, 0), &stress(2, 1),
678  &stress(2, 2));
679  FTensor::Tensor1<double *, 3> rhs(&nf[0], &nf[1], &nf[2], 3);
680  int bb = 0;
681  for (; bb != nb_dofs / 3; bb++) {
682  rhs(i) += val * t3(i, j) * diff_base_functions(j);
683  ++rhs;
684  ++diff_base_functions;
685  }
686  for (; bb != nb_base_functions; bb++) {
687  ++diff_base_functions;
688  }
689  }
690 
691  CHKERR aSemble(row_side, row_type, row_data);
692 
694 }
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
std::vector< MatrixDouble3by3 > sTress
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:483
FTensor::Tensor1< double *, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:514
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Range tEts
constrains elements in block set
MatrixBoundedArray< double, 9 > MatrixDouble3by3
Definition: Types.hpp:101
#define CHKERR
Inline error check.
Definition: definitions.h:602
FTensor::Index< 'j', 2 > j
Definition: ContactOps.hpp:27
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:413
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Index< 'i', 2 > i
[Common data]
Definition: ContactOps.hpp:26
virtual MoFEMErrorCode aSemble(int row_side, EntityType row_type, DataForcesAndSourcesCore::EntData &row_data)

Member Data Documentation

◆ aLe

bool NonlinearElasticElement::OpRhsPiolaKirchhoff::aLe

Definition at line 535 of file NonLinearElasticElement.hpp.

◆ commonData

CommonData& NonlinearElasticElement::OpRhsPiolaKirchhoff::commonData

Definition at line 533 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpRhsPiolaKirchhoff::dAta

Definition at line 532 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpRhsPiolaKirchhoff::fieldDisp

Definition at line 534 of file NonLinearElasticElement.hpp.

◆ iNdices

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

Definition at line 537 of file NonLinearElasticElement.hpp.

◆ nf

VectorDouble NonlinearElasticElement::OpRhsPiolaKirchhoff::nf

Definition at line 541 of file NonLinearElasticElement.hpp.


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