v0.15.0
Loading...
Searching...
No Matches
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.
 
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::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes
 
const EntityHandlegetConn ()
 get element connectivity
 
double getVolume () const
 element volume (linear geometry)
 
doublegetVolume ()
 element volume (linear geometry)
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian
 
VectorDoublegetCoords ()
 nodal coordinates
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object
 
- 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.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
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
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- 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.
 
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.
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

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.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 

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...
 
using AdjCache
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType
 
using DoWorkRhsHookFunType
 
- 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)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 556 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 721 of file NonLinearElasticElement.cpp.

724 : VolumeElementForcesAndSourcesCore::UserDataOperator(
725 vel_field, field_name, UserDataOperator::OPROWCOL),
726 dAta(data), commonData(common_data), aLe(false) {}
constexpr auto field_name

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 NonlinearElasticElement::OpLhsPiolaKirchhoff_dX, Smoother::OpLhsSmoother, and Smoother::OpLhsSmoother.

Definition at line 817 of file NonLinearElasticElement.cpp.

820 {
822
823 int nb_row = row_data.getIndices().size();
824 int nb_col = col_data.getIndices().size();
825
826 int *row_indices_ptr = &row_data.getIndices()[0];
827 int *col_indices_ptr = &col_data.getIndices()[0];
828
829 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
830 rowIndices.resize(nb_row, false);
831 noalias(rowIndices) = row_data.getIndices();
832 row_indices_ptr = &rowIndices[0];
833 VectorDofs &dofs = row_data.getFieldDofs();
834 VectorDofs::iterator dit = dofs.begin();
835 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
836 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
838 rowIndices[ii] = -1;
839 }
840 }
841 }
842
843 if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
844 colIndices.resize(nb_col, false);
845 noalias(colIndices) = col_data.getIndices();
846 col_indices_ptr = &colIndices[0];
847 VectorDofs &dofs = col_data.getFieldDofs();
848 VectorDofs::iterator dit = dofs.begin();
849 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
850 if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
852 colIndices[ii] = -1;
853 }
854 }
855 }
856
857 CHKERR MatSetValues(getFEMethod()->snes_B, nb_row, row_indices_ptr, nb_col,
858 col_indices_ptr, &k(0, 0), ADD_VALUES);
859
860 // is symmetric
861 if (row_side != col_side || row_type != col_type) {
862
863 row_indices_ptr = &row_data.getIndices()[0];
864 col_indices_ptr = &col_data.getIndices()[0];
865
866 if (!dAta.forcesOnlyOnEntitiesCol.empty()) {
867 rowIndices.resize(nb_row, false);
868 noalias(rowIndices) = row_data.getIndices();
869 row_indices_ptr = &rowIndices[0];
870 VectorDofs &dofs = row_data.getFieldDofs();
871 VectorDofs::iterator dit = dofs.begin();
872 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
873 if (dAta.forcesOnlyOnEntitiesCol.find((*dit)->getEnt()) ==
875 rowIndices[ii] = -1;
876 }
877 }
878 }
879
880 if (!dAta.forcesOnlyOnEntitiesRow.empty()) {
881 colIndices.resize(nb_col, false);
882 noalias(colIndices) = col_data.getIndices();
883 col_indices_ptr = &colIndices[0];
884 VectorDofs &dofs = col_data.getFieldDofs();
885 VectorDofs::iterator dit = dofs.begin();
886 for (int ii = 0; dit != dofs.end(); dit++, ii++) {
887 if (dAta.forcesOnlyOnEntitiesRow.find((*dit)->getEnt()) ==
889 colIndices[ii] = -1;
890 }
891 }
892 }
893
894 trans_k.resize(nb_col, nb_row, false);
895 noalias(trans_k) = trans(k);
896 CHKERR MatSetValues(getFEMethod()->snes_B, nb_col, col_indices_ptr, nb_row,
897 row_indices_ptr, &trans_k(0, 0), ADD_VALUES);
898 }
899
901}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
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 903 of file NonLinearElasticElement.cpp.

906 {
908
909 int nb_row = row_data.getIndices().size();
910 int nb_col = col_data.getIndices().size();
911 if (nb_row == 0)
913 if (nb_col == 0)
915
916 if (dAta.tEts.find(getNumeredEntFiniteElementPtr()->getEnt()) ==
917 dAta.tEts.end()) {
919 }
920
921 // const int nb_base_functions = row_data.getN().size2();
922 const int nb_gauss_pts = row_data.getN().size1();
923
924 FTensor::Index<'i', 3> i;
925 FTensor::Index<'j', 3> j;
926 FTensor::Index<'m', 3> m;
927
928 k.resize(nb_row, nb_col, false);
929 k.clear();
930 jac.resize(9, nb_col, false);
931
932 for (int gg = 0; gg != nb_gauss_pts; gg++) {
933 CHKERR getJac(col_data, gg);
934 double val = getVolume() * getGaussPts()(3, gg);
936 &jac(3 * 0 + 0, 0), &jac(3 * 0 + 0, 1), &jac(3 * 0 + 0, 2),
937 &jac(3 * 0 + 1, 0), &jac(3 * 0 + 1, 1), &jac(3 * 0 + 1, 2),
938 &jac(3 * 0 + 2, 0), &jac(3 * 0 + 2, 1), &jac(3 * 0 + 2, 2),
939 &jac(3 * 1 + 0, 0), &jac(3 * 1 + 0, 1), &jac(3 * 1 + 0, 2),
940 &jac(3 * 1 + 1, 0), &jac(3 * 1 + 1, 1), &jac(3 * 1 + 1, 2),
941 &jac(3 * 1 + 2, 0), &jac(3 * 1 + 2, 1), &jac(3 * 1 + 2, 2),
942 &jac(3 * 2 + 0, 0), &jac(3 * 2 + 0, 1), &jac(3 * 2 + 0, 2),
943 &jac(3 * 2 + 1, 0), &jac(3 * 2 + 1, 1), &jac(3 * 2 + 1, 2),
944 &jac(3 * 2 + 2, 0), &jac(3 * 2 + 2, 1), &jac(3 * 2 + 2, 2));
945 for (int cc = 0; cc != nb_col / 3; cc++) {
946 auto diff_base_functions = row_data.getFTensor1DiffN<3>(gg, 0);
948 &k(0, 3 * cc + 0), &k(0, 3 * cc + 1), &k(0, 3 * cc + 2),
949 &k(1, 3 * cc + 0), &k(1, 3 * cc + 1), &k(1, 3 * cc + 2),
950 &k(2, 3 * cc + 0), &k(2, 3 * cc + 1), &k(2, 3 * cc + 2), 3 * nb_col);
951 for (int rr = 0; rr != nb_row / 3; rr++) {
952 lhs(i, j) += val * t3_1(i, m, j) * diff_base_functions(m);
953 ++diff_base_functions;
954 ++lhs;
955 }
956 ++t3_1;
957 }
958 }
959
960 CHKERR aSemble(row_side, col_side, row_type, col_type, row_data, col_data);
961
963}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
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)

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

813 {
814 return get_jac<0>(col_data, gg, commonData.jacStress[gg], jac);
815}
static MoFEMErrorCode get_jac(EntitiesFieldData::EntData &col_data, int gg, MatrixDouble &jac_stress, MatrixDouble &jac)
std::vector< MatrixDouble > jacStress
this is simply material tangent operator

Member Data Documentation

◆ aLe

bool NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::aLe

Definition at line 562 of file NonLinearElasticElement.hpp.

◆ colIndices

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

Definition at line 565 of file NonLinearElasticElement.hpp.

◆ commonData

CommonData& NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::commonData

Definition at line 560 of file NonLinearElasticElement.hpp.

◆ dAta

BlockData& NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::dAta

Definition at line 559 of file NonLinearElasticElement.hpp.

◆ F

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::F

Definition at line 571 of file NonLinearElasticElement.hpp.

◆ jac

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::jac

Definition at line 571 of file NonLinearElasticElement.hpp.

◆ k

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::k

Definition at line 571 of file NonLinearElasticElement.hpp.

◆ rowIndices

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

Definition at line 564 of file NonLinearElasticElement.hpp.

◆ tAg

int NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::tAg

Definition at line 561 of file NonLinearElasticElement.hpp.

◆ trans_k

MatrixDouble NonlinearElasticElement::OpLhsPiolaKirchhoff_dx::trans_k

Definition at line 571 of file NonLinearElasticElement.hpp.


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