v0.14.0
Public Member Functions | Public Attributes | List of all members
MetaSpringBC::OpSpringFs Struct Reference

RHS-operator for the spring boundary condition element. More...

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

Inheritance diagram for MetaSpringBC::OpSpringFs:
[legend]
Collaboration diagram for MetaSpringBC::OpSpringFs:
[legend]

Public Member Functions

MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Integrates and assembles to global RHS vector virtual work of springs. More...
 
 OpSpringFs (boost::shared_ptr< MetaSpringBC::DataAtIntegrationPtsSprings > &common_data_ptr, MetaSpringBC::BlockOptionDataSprings &data, const std::string field_name)
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- 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...
 
std::string getFEName () const
 Get name of the element. More...
 
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 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...
 
doublegetMeasure ()
 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, 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. More...
 
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. More...
 
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. More...
 
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. 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 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

VectorDouble nF
 
boost::shared_ptr< MetaSpringBC::DataAtIntegrationPtsSpringscommonDataPtr
 
MetaSpringBC::BlockOptionDataSpringsdAta
 
bool is_spatial_position = true
 
- 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,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
- 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 Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

RHS-operator for the spring boundary condition element.

Integrates virtual work of springs on displacement or spatial positions and assembles components to RHS global vector.

Definition at line 134 of file SpringElement.hpp.

Constructor & Destructor Documentation

◆ OpSpringFs()

MetaSpringBC::OpSpringFs::OpSpringFs ( boost::shared_ptr< MetaSpringBC::DataAtIntegrationPtsSprings > &  common_data_ptr,
MetaSpringBC::BlockOptionDataSprings data,
const std::string  field_name 
)
inline
Parameters
common_data_ptrPointer to the common data for spring element
dataVariable containing data for normal and tangential stiffnesses of springs attached to the current element
field_nameString of field name for spatial positions or displacements for rows

Definition at line 180 of file SpringElement.hpp.

185  field_name.c_str(), OPROW),
186  commonDataPtr(common_data_ptr), dAta(data) {
187  if (field_name.compare(0, 16, "SPATIAL_POSITION") != 0)
188  is_spatial_position = false;
189  }

Member Function Documentation

◆ doWork()

MoFEMErrorCode MetaSpringBC::OpSpringFs::doWork ( int  side,
EntityType  type,
EntData data 
)
virtual

Integrates and assembles to global RHS vector virtual work of springs.

*Computes virtual work of springs on spatial positions or displacements for configurational changes:

\[ f_s({\mathbf x}, {\mathbf X}, \delta \mathbf{x}) = \int\limits_{\partial \Omega }^{} {{\delta \mathbf{x}^T} \cdot \left[ k_{\rm n} \left( {\mathbf N} \otimes {\mathbf N} \right) + k_{\rm t} \left( {\mathbf I} - {\mathbf N} \otimes {\mathbf N} \right) \right] \left( {\mathbf x} - {\mathbf X} \right) \partial \Omega } \]

where \( \delta \mathbf{x} \) is either virtual displacement or variation of spatial positions, \( k_{\rm n} \) is the stiffness of the springs normal to the surface, \( k_{\rm t} \) is the stiffness of the springs tangential to the surface, \( {\mathbf N} \) is the normal to the surface vector based on material positions, \( {\mathbf x} \) is the vector of spatial positions or displacements of the surface with springs and \( {\mathbf X} \) is the vector of material positions that is zero when displacements are considered

Reimplemented from MoFEM::DataOperator.

Definition at line 13 of file SpringElement.cpp.

14  {
15 
17  // check that the faces have associated degrees of freedom
18  const int nb_dofs = data.getIndices().size();
19  if (nb_dofs == 0)
21  if (dAta.tRis.find(getFEEntityHandle()) == dAta.tRis.end()) {
23  }
24  CHKERR commonDataPtr->getBlockData(dAta);
25 
26  // size of force vector associated to the entity
27  // set equal to the number of degrees of freedom of associated with the
28  // entity
29  nF.resize(nb_dofs, false);
30  nF.clear();
31 
32  // get number of Gauss points
33  const int nb_gauss_pts = data.getN().size1();
34 
35  // get integration weights
36  auto t_w = getFTensor0IntegrationWeight();
37 
38  // FTensor indices
42  const double normal_stiffness = commonDataPtr->springStiffnessNormal;
43  const double tangent_stiffness = commonDataPtr->springStiffnessTangent;
44 
45  // create a 3d vector to be used as the normal to the face with length equal
46  // to the face area
48  FTensor::Tensor2<double, 3, 3> t_normal_projection;
49  FTensor::Tensor2<double, 3, 3> t_tangent_projection;
50 
51  // Extract solution at Gauss points
52  auto t_solution_at_gauss_point =
53  getFTensor1FromMat<3>(*commonDataPtr->xAtPts);
54  auto t_init_solution_at_gauss_point =
55  getFTensor1FromMat<3>(*commonDataPtr->xInitAtPts);
56  FTensor::Tensor1<double, 3> t_displacement_at_gauss_point;
57 
58  auto t_normal_1 = getFTensor1FromMat<3>(commonDataPtr->normalVector);
59  // loop over all Gauss points of the face
60  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
61  t_normal(i) = t_normal_1(i);
62 
63  const double normal_length = std::sqrt(t_normal(k) * t_normal(k));
64  t_normal(i) = t_normal(i) / normal_length;
65  t_normal_projection(i, j) = t_normal(i) * t_normal(j);
66  t_tangent_projection(i, j) = t_normal_projection(i, j);
67  t_tangent_projection(0, 0) -= 1;
68  t_tangent_projection(1, 1) -= 1;
69  t_tangent_projection(2, 2) -= 1;
70  // Calculate the displacement at the Gauss point
71  if (is_spatial_position) { // "SPATIAL_POSITION"
72  t_displacement_at_gauss_point(i) =
73  t_solution_at_gauss_point(i) - t_init_solution_at_gauss_point(i);
74  } else { // e.g. "DISPLACEMENT" or "U"
75  t_displacement_at_gauss_point(i) = t_solution_at_gauss_point(i);
76  }
77 
78  double w = t_w * normal_length; // area was constant
79  if (getNumeredEntFiniteElementPtr()->getEntType() == MBTRI) {
80  w *= 0.5;
81  }
82 
83  auto t_base_func = data.getFTensor0N(gg, 0);
84 
85  // create a vector t_nf whose pointer points an array of 3 pointers
86  // pointing to nF memory location of components
88  &nF[2]);
89 
90  for (int rr = 0; rr != nb_dofs / 3; ++rr) { // loop over the nodes
91 
92  t_nf(i) += w * t_base_func * normal_stiffness *
93  t_normal_projection(i, j) * t_displacement_at_gauss_point(j);
94  t_nf(i) -= w * t_base_func * tangent_stiffness *
95  t_tangent_projection(i, j) * t_displacement_at_gauss_point(j);
96 
97  // move to next base function
98  ++t_base_func;
99  // move the pointer to next element of t_nf
100  ++t_nf;
101  }
102  // move to next integration weight
103  ++t_w;
104  // move to the solutions at the next Gauss point
105  ++t_solution_at_gauss_point;
106  ++t_init_solution_at_gauss_point;
107  ++t_normal_1;
108  }
109  // add computed values of spring in the global right hand side vector
110  Vec f = getFEMethod()->ksp_f != PETSC_NULL ? getFEMethod()->ksp_f
111  : getFEMethod()->snes_f;
112  CHKERR VecSetValues(f, nb_dofs, &data.getIndices()[0], &nF[0], ADD_VALUES);
114 }

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<MetaSpringBC::DataAtIntegrationPtsSprings> MetaSpringBC::OpSpringFs::commonDataPtr

Definition at line 139 of file SpringElement.hpp.

◆ dAta

MetaSpringBC::BlockOptionDataSprings& MetaSpringBC::OpSpringFs::dAta

Definition at line 140 of file SpringElement.hpp.

◆ is_spatial_position

bool MetaSpringBC::OpSpringFs::is_spatial_position = true

Definition at line 141 of file SpringElement.hpp.

◆ nF

VectorDouble MetaSpringBC::OpSpringFs::nF

Definition at line 137 of file SpringElement.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MetaSpringBC::BlockOptionDataSprings::tRis
Range tRis
Definition: SpringElement.hpp:26
PlasticOps::w
double w(double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y)
Definition: PlasticOpsGeneric.hpp:276
FTensor::Tensor1< double, 3 >
MetaSpringBC::OpSpringFs::commonDataPtr
boost::shared_ptr< MetaSpringBC::DataAtIntegrationPtsSprings > commonDataPtr
Definition: SpringElement.hpp:139
MetaSpringBC::OpSpringFs::nF
VectorDouble nF
Definition: SpringElement.hpp:137
MoFEM::VecSetValues
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
Definition: EntitiesFieldData.hpp:1576
FTensor::Tensor2< double, 3, 3 >
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeight
auto getFTensor0IntegrationWeight()
Get integration weights.
Definition: ForcesAndSourcesCore.hpp:1239
MoFEM::EntitiesFieldData::EntData::getFTensor0N
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
Definition: EntitiesFieldData.hpp:1489
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
friend class UserDataOperator
Definition: FaceElementForcesAndSourcesCore.hpp:86
MetaSpringBC::OpSpringFs::dAta
MetaSpringBC::BlockOptionDataSprings & dAta
Definition: SpringElement.hpp:140
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
MoFEM::ForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
Definition: ForcesAndSourcesCore.hpp:999
MetaSpringBC::OpSpringFs::is_spatial_position
bool is_spatial_position
Definition: SpringElement.hpp:141
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', 3 >
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEMethod
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
Definition: ForcesAndSourcesCore.hpp:1041
MoFEM::SnesMethod::snes_f
Vec & snes_f
residual
Definition: LoopMethods.hpp:122
MoFEM::ForcesAndSourcesCore::UserDataOperator::getFEEntityHandle
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
Definition: ForcesAndSourcesCore.hpp:1003
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1305
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
MoFEM::KspMethod::ksp_f
Vec & ksp_f
Definition: LoopMethods.hpp:89
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MoFEM::ForcesAndSourcesCore::UserDataOperator::OPROW
@ OPROW
operator doWork function is executed on FE rows
Definition: ForcesAndSourcesCore.hpp:567