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

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

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

Public Member Functions

 OpLhsPiolaKirchhoff_dx (const std::string vel_field, const std::string field_name, BlockData &data, CommonData &common_data)
 
virtual MoFEMErrorCode getJac (EntitiesFieldData::EntData &col_data, int gg)
 Directive of Piola Kirchhoff stress over spatial DOFs. More...
 
virtual MoFEMErrorCode aSemble (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_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
 
int tAg
 
bool aLe
 
ublas::vector< int > rowIndices
 
ublas::vector< int > colIndices
 
MatrixDouble k
 
MatrixDouble trans_k
 
MatrixDouble jac
 
MatrixDouble F
 
- 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 569 of file NonLinearElasticElement.hpp.

Constructor & Destructor Documentation

◆ OpLhsPiolaKirchhoff_dx()

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

Definition at line 734 of file NonLinearElasticElement.cpp.

738  vel_field, field_name, UserDataOperator::OPROWCOL),
739  dAta(data), commonData(common_data), aLe(false) {}
ForcesAndSourcesCore::UserDataOperator UserDataOperator

Member Function Documentation

◆ aSemble()

MoFEMErrorCode NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::aSemble ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)
virtual

Reimplemented in Smoother::OpLhsSmoother, and NonlinearElasticElement::OpLhsPiolaKirchhoff_dX.

Definition at line 829 of file NonLinearElasticElement.cpp.

832  {
834 
835  int nb_row = row_data.getIndices().size();
836  int nb_col = col_data.getIndices().size();
837 
838  int *row_indices_ptr = &row_data.getIndices()[0];
839  int *col_indices_ptr = &col_data.getIndices()[0];
840 
841  if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
842  rowIndices.resize(nb_row, false);
843  noalias(rowIndices) = row_data.getIndices();
844  row_indices_ptr = &rowIndices[0];
845  VectorDofs &dofs = row_data.getFieldDofs();
846  VectorDofs::iterator dit = dofs.begin();
847  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
848  if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
850  rowIndices[ii] = -1;
851  }
852  }
853  }
854 
855  if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
856  colIndices.resize(nb_col, false);
857  noalias(colIndices) = col_data.getIndices();
858  col_indices_ptr = &colIndices[0];
859  VectorDofs &dofs = col_data.getFieldDofs();
860  VectorDofs::iterator dit = dofs.begin();
861  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
862  if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
864  colIndices[ii] = -1;
865  }
866  }
867  }
868 
869  CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr, nb_col,
870  col_indices_ptr, &k(0, 0), ADD_VALUES);
871 
872  // is symmetric
873  if (row_side != col_side || row_type != col_type) {
874 
875  row_indices_ptr = &row_data.getIndices()[0];
876  col_indices_ptr = &col_data.getIndices()[0];
877 
878  if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
879  rowIndices.resize(nb_row, false);
880  noalias(rowIndices) = row_data.getIndices();
881  row_indices_ptr = &rowIndices[0];
882  VectorDofs &dofs = row_data.getFieldDofs();
883  VectorDofs::iterator dit = dofs.begin();
884  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
885  if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
887  rowIndices[ii] = -1;
888  }
889  }
890  }
891 
892  if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
893  colIndices.resize(nb_col, false);
894  noalias(colIndices) = col_data.getIndices();
895  col_indices_ptr = &colIndices[0];
896  VectorDofs &dofs = col_data.getFieldDofs();
897  VectorDofs::iterator dit = dofs.begin();
898  for (int ii = 0; dit != dofs.end(); dit++, ii++) {
899  if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
901  colIndices[ii] = -1;
902  }
903  }
904  }
905 
906  trans_k.resize(nb_col, nb_row, false);
907  noalias(trans_k) = trans(k);
908  CHKERR MatSetValues(getFEMethod()->snes_B, nb_col, col_indices_ptr, nb_row,
909  row_indices_ptr, &trans_k(0, 0), ADD_VALUES);
910  }
911 
913 }
#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::OpLhsPiolaKirchhoff_dx::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData row_data,
EntitiesFieldData::EntData col_data 
)

Definition at line 915 of file NonLinearElasticElement.cpp.

918  {
920 
921  int nb_row = row_data.getIndices().size();
922  int nb_col = col_data.getIndices().size();
923  if (nb_row == 0)
925  if (nb_col == 0)
927 
928  if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
929  dAta.tEts.end()) {
931  }
932 
933  // const int nb_base_functions = row_data.getN().size2();
934  const int nb_gauss_pts = row_data.getN().size1();
935 
939 
940  k.resize(nb_row, nb_col, false);
941  k.clear();
942  jac.resize(9, nb_col, false);
943 
944  for (int gg = 0; gg != nb_gauss_pts; gg++) {
945  CHKERR getJac(col_data, gg);
946  double val = getVolume() * getGaussPts()(3, gg);
948  &jac(3 * 0 + 0, 0), &jac(3 * 0 + 0, 1), &jac(3 * 0 + 0, 2),
949  &jac(3 * 0 + 1, 0), &jac(3 * 0 + 1, 1), &jac(3 * 0 + 1, 2),
950  &jac(3 * 0 + 2, 0), &jac(3 * 0 + 2, 1), &jac(3 * 0 + 2, 2),
951  &jac(3 * 1 + 0, 0), &jac(3 * 1 + 0, 1), &jac(3 * 1 + 0, 2),
952  &jac(3 * 1 + 1, 0), &jac(3 * 1 + 1, 1), &jac(3 * 1 + 1, 2),
953  &jac(3 * 1 + 2, 0), &jac(3 * 1 + 2, 1), &jac(3 * 1 + 2, 2),
954  &jac(3 * 2 + 0, 0), &jac(3 * 2 + 0, 1), &jac(3 * 2 + 0, 2),
955  &jac(3 * 2 + 1, 0), &jac(3 * 2 + 1, 1), &jac(3 * 2 + 1, 2),
956  &jac(3 * 2 + 2, 0), &jac(3 * 2 + 2, 1), &jac(3 * 2 + 2, 2));
957  for (int cc = 0; cc != nb_col / 3; cc++) {
958  auto diff_base_functions = row_data.getFTensor1DiffN<3>(gg, 0);
960  &k(0, 3 * cc + 0), &k(0, 3 * cc + 1), &k(0, 3 * cc + 2),
961  &k(1, 3 * cc + 0), &k(1, 3 * cc + 1), &k(1, 3 * cc + 2),
962  &k(2, 3 * cc + 0), &k(2, 3 * cc + 1), &k(2, 3 * cc + 2), 3 * nb_col);
963  for (int rr = 0; rr != nb_row / 3; rr++) {
964  lhs(i, j) += val * t3_1(i, m, j) * diff_base_functions(m);
965  ++diff_base_functions;
966  ++lhs;
967  }
968  ++t3_1;
969  }
970  }
971 
972  CHKERR aSemble(row_side, col_side, row_type, col_type, row_data, col_data);
973 
975 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
FTensor::Index< 'm', 3 > m
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
FTensor::Tensor2< double *, 3, 3 > & getJac()
get element Jacobian
Range tEts
constrains elements in block set
virtual MoFEMErrorCode aSemble(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
FTensor::Index< 'j', 3 > j
FTensor::Index< 'i', 3 > i

◆ getJac()

MoFEMErrorCode NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::getJac ( EntitiesFieldData::EntData col_data,
int  gg 
)
virtual

Directive of Piola Kirchhoff stress over spatial DOFs.

This project derivative \(\frac{\partial P}{\partial F}\), that is

\[ \frac{\partial P}{\partial x_\textrm{DOF}} = \frac{\partial P}{\partial F}\frac{\partial F}{\partial x_\textrm{DOF}}, \]

where second therm \(\frac{\partial F}{\partial x_\textrm{DOF}}\) is derivative of shape function

Reimplemented in NonlinearElasticElement::OpLhsEshelby_dX, NonlinearElasticElement::OpLhsEshelby_dx, and NonlinearElasticElement::OpLhsPiolaKirchhoff_dX.

Definition at line 824 of file NonLinearElasticElement.cpp.

825  {
826  return get_jac<0>(col_data, gg, commonData.jacStress[gg], jac);
827 }
std::vector< MatrixDouble > jacStress
this is simply material tangent operator

Member Data Documentation

◆ aLe

bool NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::aLe

Definition at line 575 of file NonLinearElasticElement.hpp.

◆ colIndices

ublas::vector<int> NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::colIndices

Definition at line 578 of file NonLinearElasticElement.hpp.

◆ commonData

CommonData& NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::commonData

Definition at line 573 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::dAta

Definition at line 572 of file NonLinearElasticElement.hpp.

◆ F

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::F

Definition at line 584 of file NonLinearElasticElement.hpp.

◆ jac

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::jac

Definition at line 584 of file NonLinearElasticElement.hpp.

◆ k

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::k

Definition at line 584 of file NonLinearElasticElement.hpp.

◆ rowIndices

ublas::vector<int> NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::rowIndices

Definition at line 577 of file NonLinearElasticElement.hpp.

◆ tAg

int NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::tAg

Definition at line 574 of file NonLinearElasticElement.hpp.

◆ trans_k

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::trans_k

Definition at line 584 of file NonLinearElasticElement.hpp.


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