v0.9.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, 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...
 
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () 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, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
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, bool symm=true)
 
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 do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
 
virtual MoFEMErrorCode opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
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< intiNdices
 
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...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 

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::ForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
cell_forces.cpp.

Definition at line 528 of file NonLinearElasticElement.hpp.

Constructor & Destructor Documentation

◆ OpRhsPiolaKirchhoff()

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

Definition at line 631 of file NonLinearElasticElement.cpp.

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 637 of file NonLinearElasticElement.cpp.

639  {
641 
642  int nb_dofs = row_data.getIndices().size();
643  int *indices_ptr = &row_data.getIndices()[0];
644  if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
645  iNdices.resize(nb_dofs, false);
646  noalias(iNdices) = row_data.getIndices();
647  indices_ptr = &iNdices[0];
648  VectorDofs &dofs = row_data.getFieldDofs();
649  VectorDofs::iterator dit = dofs.begin();
650  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
651  if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
653  iNdices[ii] = -1;
654  }
655  }
656  }
657  CHKERR VecSetValues(getFEMethod()->snes_f, nb_dofs, indices_ptr, &nf[0],
658  ADD_VALUES);
660 }
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:596
ublas::vector< boost::shared_ptr< const FEDofEntity >, DofsAllocator > VectorDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

◆ doWork()

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

Definition at line 662 of file NonLinearElasticElement.cpp.

664  {
666 
667  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
668  dAta.tEts.end()) {
670  }
671 
672  const int nb_dofs = row_data.getIndices().size();
673  if (nb_dofs == 0)
675  if ((unsigned int)nb_dofs > 3 * row_data.getN().size2()) {
676  SETERRQ(PETSC_COMM_SELF, 1, "data inconsistency");
677  }
678  const int nb_base_functions = row_data.getN().size2();
679  const int nb_gauss_pts = row_data.getN().size1();
680 
681  nf.resize(nb_dofs, false);
682  nf.clear();
683 
684  FTensor::Tensor1<double *, 3> diff_base_functions =
685  row_data.getFTensor1DiffN<3>();
688 
689  for (int gg = 0; gg != nb_gauss_pts; gg++) {
690  double val = getVolume() * getGaussPts()(3, gg);
691  if ((!aLe) && getHoGaussPtsDetJac().size() > 0) {
692  val *= getHoGaussPtsDetJac()[gg]; ///< higher order geometry
693  }
694  MatrixDouble3by3 &stress = commonData.sTress[gg];
696  &stress(0, 0), &stress(0, 1), &stress(0, 2), &stress(1, 0),
697  &stress(1, 1), &stress(1, 2), &stress(2, 0), &stress(2, 1),
698  &stress(2, 2));
699  FTensor::Tensor1<double *, 3> rhs(&nf[0], &nf[1], &nf[2], 3);
700  int bb = 0;
701  for (; bb != nb_dofs / 3; bb++) {
702  rhs(i) += val * t3(i, j) * diff_base_functions(j);
703  ++rhs;
704  ++diff_base_functions;
705  }
706  for (; bb != nb_base_functions; bb++) {
707  ++diff_base_functions;
708  }
709  }
710 
711  CHKERR aSemble(row_side, row_type, row_data);
712 
714 }
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:477
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:508
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:97
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
virtual MoFEMErrorCode aSemble(int row_side, EntityType row_type, DataForcesAndSourcesCore::EntData &row_data)

Member Data Documentation

◆ aLe

bool NonlinearElasticElement::OpRhsPiolaKirchhoff::aLe

Definition at line 534 of file NonLinearElasticElement.hpp.

◆ commonData

CommonData& NonlinearElasticElement::OpRhsPiolaKirchhoff::commonData

Definition at line 532 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpRhsPiolaKirchhoff::dAta

Definition at line 531 of file NonLinearElasticElement.hpp.

◆ fieldDisp

bool NonlinearElasticElement::OpRhsPiolaKirchhoff::fieldDisp

Definition at line 533 of file NonLinearElasticElement.hpp.

◆ iNdices

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

Definition at line 536 of file NonLinearElasticElement.hpp.

◆ nf

VectorDouble NonlinearElasticElement::OpRhsPiolaKirchhoff::nf

Definition at line 540 of file NonLinearElasticElement.hpp.


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