v0.15.0
Loading...
Searching...
No Matches
OpSimpleRodK Struct Reference

Assemble contribution of SimpleRod element to LHS *. More...

Inheritance diagram for OpSimpleRodK:
[legend]
Collaboration diagram for OpSimpleRodK:
[legend]

Public Member Functions

 OpSimpleRodK (boost::shared_ptr< DataAtIntegrationPtsSimpleRods > &common_data_ptr, BlockOptionDataSimpleRods &data, const std::string field_name)
 
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.
 
- Public Member Functions inherited from MoFEM::EdgeElementForcesAndSourcesCore::UserDataOperator
const EntityHandlegetConn ()
 get element connectivity
 
double getLength ()
 get edge length
 
VectorDoublegetDirection ()
 get edge direction
 
auto getFTensor1Normal ()
 get edge normal NOTE: it should be used only in 2D analysis
 
auto getFTensor1Normal (const FTensor::Tensor1< double, 3 > &vec)
 get ftensor1 edge normal
 
auto getFTensor1NormalsAtGaussPts (const FTensor::Tensor1< double, 3 > &vec)
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points
 
VectorDoublegetCoords ()
 get edge node coordinates
 
MatrixDoublegetTangentAtGaussPts ()
 get tangent vector to edge curve at integration points
 
DEPRECATED MatrixDoublegetTangetAtGaussPts ()
 
const EdgeElementForcesAndSourcesCoregetEdgeFE ()
 get pointer to this finite element
 
FTensor::Tensor1< double, 3 > getFTensor1Direction ()
 
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1Coords ()
 get get coords at gauss points
 
DEPRECATED FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getTensor1Coords ()
 
template<int DIM = 3>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, DIM > getFTensor1TangentAtGaussPts ()
 Get tangent vector at integration points.
 
MoFEMErrorCode loopSideFaces (const string fe_name, FaceElementForcesAndSourcesCoreOnSide &fe_side)
 
template<>
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getFTensor1TangentAtGaussPts ()
 
- 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 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

boost::shared_ptr< DataAtIntegrationPtsSimpleRodscommonDataPtr
 
BlockOptionDataSimpleRodsdAta
 
MatrixDouble locK
 
MatrixDouble transLocK
 
- 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::EdgeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Assemble contribution of SimpleRod element to LHS *.

Definition at line 111 of file SimpleRodElement.cpp.

Constructor & Destructor Documentation

◆ OpSimpleRodK()

OpSimpleRodK::OpSimpleRodK ( boost::shared_ptr< DataAtIntegrationPtsSimpleRods > & common_data_ptr,
BlockOptionDataSimpleRods & data,
const std::string field_name )
inline

Definition at line 119 of file SimpleRodElement.cpp.

123 field_name.c_str(), field_name.c_str(), OPROWCOL),
124 commonDataPtr(common_data_ptr), dAta(data) {
125 sYmm = false;
126 }
constexpr auto field_name
bool sYmm
If true assume that matrix is symmetric structure.
@ OPROWCOL
operator doWork is executed on FE rows &columns
BlockOptionDataSimpleRods & dAta
boost::shared_ptr< DataAtIntegrationPtsSimpleRods > commonDataPtr

Member Function Documentation

◆ doWork()

MoFEMErrorCode OpSimpleRodK::doWork ( int row_side,
int col_side,
EntityType row_type,
EntityType col_type,
EntitiesFieldData::EntData & row_data,
EntitiesFieldData::EntData & col_data )
inlinevirtual

Operator for bi-linear form, usually to calculate values on left hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 128 of file SimpleRodElement.cpp.

131 {
133
134 // check if the edge have associated degrees of freedom
135 const int row_nb_dofs = row_data.getIndices().size();
136 if (!row_nb_dofs)
138
139 const int col_nb_dofs = col_data.getIndices().size();
140 if (!col_nb_dofs)
142
143 if (dAta.eDges.find(getFEEntityHandle()) == dAta.eDges.end()) {
145 }
146
147 CHKERR commonDataPtr->getBlockData(dAta);
148 // size associated to the entity
149 locK.resize(row_nb_dofs, col_nb_dofs, false);
150 locK.clear();
151
152 double tension_stiffness = commonDataPtr->simpleRodYoungModulus *
153 commonDataPtr->simpleRodSectionArea;
154
155 VectorDouble coords;
156 coords = getCoords();
157 double L = getLength();
158 double coeff = tension_stiffness / (L * L * L);
159
160 double x21 = coords(3) - coords(0);
161 double y21 = coords(4) - coords(1);
162 double z21 = coords(5) - coords(2);
163
164 // Calculate element matrix
165 locK(0, 0) = coeff * x21 * x21;
166 locK(0, 1) = coeff * x21 * y21;
167 locK(0, 2) = coeff * x21 * z21;
168 locK(0, 3) = -coeff * x21 * x21;
169 locK(0, 4) = -coeff * x21 * y21;
170 locK(0, 5) = -coeff * x21 * z21;
171
172 locK(1, 0) = locK(0, 1);
173 locK(1, 1) = coeff * y21 * y21;
174 locK(1, 2) = coeff * y21 * z21;
175 locK(1, 3) = -coeff * y21 * x21;
176 locK(1, 4) = -coeff * y21 * y21;
177 locK(1, 5) = -coeff * y21 * z21;
178
179 locK(2, 0) = locK(0, 2);
180 locK(2, 1) = locK(1, 2);
181 locK(2, 2) = coeff * z21 * z21;
182 locK(2, 3) = -coeff * z21 * x21;
183 locK(2, 4) = -coeff * z21 * y21;
184 locK(2, 5) = -coeff * z21 * z21;
185
186 locK(3, 0) = locK(0, 3);
187 locK(3, 1) = locK(1, 3);
188 locK(3, 2) = locK(2, 3);
189 locK(3, 3) = coeff * x21 * x21;
190 locK(3, 4) = coeff * x21 * y21;
191 locK(3, 5) = coeff * x21 * z21;
192
193 locK(4, 0) = locK(0, 4);
194 locK(4, 1) = locK(1, 4);
195 locK(4, 2) = locK(2, 4);
196 locK(4, 3) = locK(3, 4);
197 locK(4, 4) = coeff * y21 * y21;
198 locK(4, 5) = coeff * y21 * z21;
199
200 locK(5, 0) = locK(0, 5);
201 locK(5, 1) = locK(1, 5);
202 locK(5, 2) = locK(2, 5);
203 locK(5, 3) = locK(3, 5);
204 locK(5, 4) = locK(4, 5);
205 locK(5, 5) = coeff * z21 * z21;
206
207 CHKERR MatSetValues(getKSPB(), row_data, col_data, &locK(0, 0), ADD_VALUES);
208
210 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#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.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
MatrixDouble locK

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<DataAtIntegrationPtsSimpleRods> OpSimpleRodK::commonDataPtr

Definition at line 113 of file SimpleRodElement.cpp.

◆ dAta

BlockOptionDataSimpleRods& OpSimpleRodK::dAta

Definition at line 114 of file SimpleRodElement.cpp.

◆ locK

MatrixDouble OpSimpleRodK::locK

Definition at line 116 of file SimpleRodElement.cpp.

◆ transLocK

MatrixDouble OpSimpleRodK::transLocK

Definition at line 117 of file SimpleRodElement.cpp.


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