v0.14.0
Classes | Static Public Member Functions | List of all members
MetaSpringBC Struct Reference

Set of functions declaring elements and setting operators to apply spring boundary condition. More...

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

Classes

struct  BlockOptionDataSprings
 
struct  DataAtIntegrationPtsSprings
 
struct  OpCalculateDeformation
 Operator for computing deformation gradients in side volumes. More...
 
struct  OpGetNormalSpEle
 Computes, for material configuration, normal to face that lies on a surface with springs. More...
 
struct  OpGetTangentSpEle
 Computes, for material configuration, tangent vectors to face that lie on a surface with springs. More...
 
struct  OpSpringALEMaterialLhs
 Base class for LHS-operators for pressure element (material configuration) More...
 
struct  OpSpringALEMaterialLhs_dX_dx
 LHS-operator for the pressure element (material configuration) More...
 
struct  OpSpringALEMaterialLhs_dX_dX
 LHS-operator for the pressure element (material configuration) More...
 
struct  OpSpringFs
 RHS-operator for the spring boundary condition element. More...
 
struct  OpSpringFsMaterial
 RHS-operator for the spring boundary condition element for ALE formulation. More...
 
struct  OpSpringKs
 LHS-operator for the springs element. More...
 
struct  OpSpringKs_dX
 
struct  SpringALEMaterialVolOnSideLhs
 Base class for LHS-operators (material) on side volumes. More...
 
struct  SpringALEMaterialVolOnSideLhs_dX_dx
 LHS-operator (material configuration) on the side volume for spring element. More...
 
struct  SpringALEMaterialVolOnSideLhs_dX_dX
 LHS-operator (material configuration) on the side volume. More...
 

Static Public Member Functions

static MoFEMErrorCode addSpringElements (MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 Declare spring element. More...
 
static MoFEMErrorCode addSpringElementsALE (MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions, Range &spring_triangles)
 Declare spring element in ALE. More...
 
static MoFEMErrorCode setSpringOperators (MoFEM::Interface &m_field, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_lhs_ptr, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_rhs_ptr, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS", double stiffness_scale=1.)
 Implementation of spring element. Set operators to calculate LHS and RHS. More...
 
static MoFEMErrorCode setSpringOperatorsMaterial (MoFEM::Interface &m_field, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_lhs_ptr_dx, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_lhs_ptr_dX, boost::shared_ptr< FaceElementForcesAndSourcesCore > fe_spring_rhs_ptr, boost::shared_ptr< DataAtIntegrationPtsSprings > data_at_integration_pts, const std::string field_name, const std::string mesh_nodals_positions, std::string side_fe_name)
 Implementation of spring element. Set operators to calculate LHS and RHS. More...
 

Detailed Description

Set of functions declaring elements and setting operators to apply spring boundary condition.

Definition at line 19 of file SpringElement.hpp.

Member Function Documentation

◆ addSpringElements()

MoFEMErrorCode MetaSpringBC::addSpringElements ( MoFEM::Interface m_field,
const std::string  field_name,
const std::string  mesh_nodals_positions = "MESH_NODE_POSITIONS" 
)
static

Declare spring element.

Search cubit sidesets and blocksets with spring bc and declare surface element

Blockset has to have name “SPRING_BC”. The first two attributes of the blockset are spring stiffness value.

Parameters
m_fieldInterface instance
field_nameField name (e.g. SPATIAL_POSITION)
mesh_nodals_positionsName of field on which ho-geometry is defined
Returns
Error code
Examples
elasticity.cpp, mortar_contact_thermal.cpp, simple_contact.cpp, and simple_contact_thermal.cpp.

Definition at line 1127 of file SpringElement.cpp.

1129  {
1131 
1132  // Define boundary element that operates on rows, columns and data of a
1133  // given field
1134  CHKERR m_field.add_finite_element("SPRING", MF_ZERO);
1138  if (m_field.check_field(mesh_nodals_positions)) {
1140  mesh_nodals_positions);
1141  }
1142  // Add entities to that element, here we add all triangles with SPRING_BC
1143  // from cubit
1145  if (bit->getName().compare(0, 9, "SPRING_BC") == 0) {
1146  CHKERR m_field.add_ents_to_finite_element_by_dim(bit->getMeshset(), 2,
1147  "SPRING");
1148  }
1149  }
1150 
1152 }

◆ addSpringElementsALE()

MoFEMErrorCode MetaSpringBC::addSpringElementsALE ( MoFEM::Interface m_field,
const std::string  field_name,
const std::string  mesh_nodals_positions,
Range spring_triangles 
)
static

Declare spring element in ALE.

Search cubit sidesets and blocksets with spring bc and declare surface element

Blockset has to have name “SPRING_BC”. The first two attributes of the blockset are spring stiffness value.

Parameters
m_fieldInterface instance
field_nameField name (e.g. SPATIAL_POSITION)
mesh_nodals_positionsName of field on which ho-geometry is defined
Returns
Error code

Definition at line 1154 of file SpringElement.cpp.

1156  {
1158 
1159  // Define boundary element that operates on rows, columns and data of a
1160  // given field
1161  CHKERR m_field.add_finite_element("SPRING_ALE", MF_ZERO);
1162  CHKERR m_field.modify_finite_element_add_field_row("SPRING_ALE", field_name);
1163  CHKERR m_field.modify_finite_element_add_field_col("SPRING_ALE", field_name);
1165  CHKERR m_field.modify_finite_element_add_field_row("SPRING_ALE",
1166  mesh_nodals_positions);
1167  CHKERR m_field.modify_finite_element_add_field_col("SPRING_ALE",
1168  mesh_nodals_positions);
1169  CHKERR m_field.modify_finite_element_add_field_data("SPRING_ALE",
1170  mesh_nodals_positions);
1171 
1172  CHKERR m_field.add_ents_to_finite_element_by_dim(spring_triangles, 2,
1173  "SPRING_ALE");
1174 
1176 }

◆ setSpringOperators()

MoFEMErrorCode MetaSpringBC::setSpringOperators ( MoFEM::Interface m_field,
boost::shared_ptr< FaceElementForcesAndSourcesCore fe_spring_lhs_ptr,
boost::shared_ptr< FaceElementForcesAndSourcesCore fe_spring_rhs_ptr,
const std::string  field_name,
const std::string  mesh_nodals_positions = "MESH_NODE_POSITIONS",
double  stiffness_scale = 1. 
)
static

Implementation of spring element. Set operators to calculate LHS and RHS.

Parameters
m_fieldInterface instance
fe_spring_lhs_ptrPointer to the FE instance for LHS
fe_spring_rhs_ptrPointer to the FE instance for RHS
field_nameField name (e.g. SPATIAL_POSITION)
mesh_nodals_positionsName of field on which ho-geometry is defined
Returns
Error code
Examples
elasticity.cpp, mortar_contact_thermal.cpp, simple_contact.cpp, and simple_contact_thermal.cpp.

Definition at line 1178 of file SpringElement.cpp.

1182  {
1184 
1185  // Push operators to instances for springs
1186  // loop over blocks
1187  boost::shared_ptr<MetaSpringBC::DataAtIntegrationPtsSprings> commonDataPtr =
1188  boost::make_shared<MetaSpringBC::DataAtIntegrationPtsSprings>(m_field, stiffness_scale);
1189  CHKERR commonDataPtr->getParameters();
1190 
1191  if (m_field.check_field(mesh_nodals_positions)) {
1192  CHKERR addHOOpsFace3D(mesh_nodals_positions, *fe_spring_lhs_ptr, false,
1193  false);
1194  CHKERR addHOOpsFace3D(mesh_nodals_positions, *fe_spring_rhs_ptr, false,
1195  false);
1196  }
1197 
1198  for (auto &sitSpring : commonDataPtr->mapSpring) {
1199 
1200  fe_spring_lhs_ptr->getOpPtrVector().push_back(
1201  new OpGetTangentSpEle(mesh_nodals_positions, commonDataPtr));
1202  fe_spring_lhs_ptr->getOpPtrVector().push_back(
1203  new OpGetNormalSpEle(mesh_nodals_positions, commonDataPtr));
1204  fe_spring_lhs_ptr->getOpPtrVector().push_back(
1205  new OpSpringKs(commonDataPtr, sitSpring.second, field_name));
1206 
1207  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1208  new OpCalculateVectorFieldValues<3>(field_name, commonDataPtr->xAtPts));
1209  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1210  new OpCalculateVectorFieldValues<3>(mesh_nodals_positions,
1211  commonDataPtr->xInitAtPts));
1212  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1213  new OpGetTangentSpEle(mesh_nodals_positions, commonDataPtr));
1214  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1215  new OpGetNormalSpEle(mesh_nodals_positions, commonDataPtr));
1216  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1217  new OpSpringFs(commonDataPtr, sitSpring.second, field_name));
1218  }
1219 
1221 }

◆ setSpringOperatorsMaterial()

MoFEMErrorCode MetaSpringBC::setSpringOperatorsMaterial ( MoFEM::Interface m_field,
boost::shared_ptr< FaceElementForcesAndSourcesCore fe_spring_lhs_ptr_dx,
boost::shared_ptr< FaceElementForcesAndSourcesCore fe_spring_lhs_ptr_dX,
boost::shared_ptr< FaceElementForcesAndSourcesCore fe_spring_rhs_ptr,
boost::shared_ptr< DataAtIntegrationPtsSprings data_at_integration_pts,
const std::string  field_name,
const std::string  mesh_nodals_positions,
std::string  side_fe_name 
)
static

Implementation of spring element. Set operators to calculate LHS and RHS.

Parameters
m_fieldInterface instance
fe_spring_lhs_ptrPointer to the FE instance for LHS
fe_spring_rhs_ptrPointer to the FE instance for RHS
field_nameField name (e.g. SPATIAL_POSITION)
mesh_nodals_positionsName of field on which ho-geometry is defined
Returns
Error code

Definition at line 1223 of file SpringElement.cpp.

1231  {
1233 
1234  // Push operators to instances for springs
1235  // loop over blocks
1236  CHKERR data_at_integration_pts->getParameters();
1237 
1238  if (m_field.check_field(mesh_nodals_positions)) {
1239  CHKERR addHOOpsFace3D(mesh_nodals_positions, *fe_spring_lhs_ptr_dx, false,
1240  false);
1241  CHKERR addHOOpsFace3D(mesh_nodals_positions, *fe_spring_lhs_ptr_dX, false,
1242  false);
1243  CHKERR addHOOpsFace3D(mesh_nodals_positions, *fe_spring_rhs_ptr, false,
1244  false);
1245  }
1246 
1247  for (auto &sitSpring : data_at_integration_pts->mapSpring) {
1248 
1249  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> feMatSideRhs =
1250  boost::make_shared<VolumeElementForcesAndSourcesCoreOnSide>(m_field);
1251 
1252  feMatSideRhs->getOpPtrVector().push_back(
1254  mesh_nodals_positions, data_at_integration_pts->HMat));
1255  feMatSideRhs->getOpPtrVector().push_back(
1257  field_name, data_at_integration_pts->hMat));
1258 
1259  feMatSideRhs->getOpPtrVector().push_back(new OpCalculateDeformation(
1260  mesh_nodals_positions, data_at_integration_pts));
1261 
1262  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1263  new OpGetTangentSpEle(mesh_nodals_positions, data_at_integration_pts));
1264  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1265  new OpGetNormalSpEle(mesh_nodals_positions, data_at_integration_pts));
1266 
1267  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1269  data_at_integration_pts->xAtPts));
1270  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1272  mesh_nodals_positions, data_at_integration_pts->xInitAtPts));
1273 
1274  fe_spring_rhs_ptr->getOpPtrVector().push_back(
1275  new OpSpringFsMaterial(mesh_nodals_positions, data_at_integration_pts,
1276  feMatSideRhs, side_fe_name, sitSpring.second));
1277 
1278  fe_spring_lhs_ptr_dx->getOpPtrVector().push_back(
1280  data_at_integration_pts->xAtPts));
1281  fe_spring_lhs_ptr_dx->getOpPtrVector().push_back(
1283  mesh_nodals_positions, data_at_integration_pts->xInitAtPts));
1284  fe_spring_lhs_ptr_dx->getOpPtrVector().push_back(
1285  new OpGetTangentSpEle(mesh_nodals_positions, data_at_integration_pts));
1286  fe_spring_lhs_ptr_dx->getOpPtrVector().push_back(
1287  new OpGetNormalSpEle(mesh_nodals_positions, data_at_integration_pts));
1288 
1289  fe_spring_lhs_ptr_dx->getOpPtrVector().push_back(
1290  new OpSpringKs_dX(data_at_integration_pts, sitSpring.second, field_name,
1291  mesh_nodals_positions));
1292 
1293  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> feMatSideLhs_dx =
1294  boost::make_shared<VolumeElementForcesAndSourcesCoreOnSide>(m_field);
1295 
1296  feMatSideLhs_dx->getOpPtrVector().push_back(
1298  mesh_nodals_positions, data_at_integration_pts->HMat));
1299  feMatSideLhs_dx->getOpPtrVector().push_back(
1301  field_name, data_at_integration_pts->hMat));
1302 
1303  feMatSideLhs_dx->getOpPtrVector().push_back(new OpCalculateDeformation(
1304  mesh_nodals_positions, data_at_integration_pts));
1305 
1306  feMatSideLhs_dx->getOpPtrVector().push_back(
1307  new SpringALEMaterialVolOnSideLhs_dX_dx(
1308  mesh_nodals_positions, field_name, data_at_integration_pts));
1309 
1310  fe_spring_lhs_ptr_dX->getOpPtrVector().push_back(
1312  data_at_integration_pts->xAtPts));
1313  fe_spring_lhs_ptr_dX->getOpPtrVector().push_back(
1315  mesh_nodals_positions, data_at_integration_pts->xInitAtPts));
1316  fe_spring_lhs_ptr_dX->getOpPtrVector().push_back(
1317  new OpGetTangentSpEle(mesh_nodals_positions, data_at_integration_pts));
1318  fe_spring_lhs_ptr_dX->getOpPtrVector().push_back(
1319  new OpGetNormalSpEle(mesh_nodals_positions, data_at_integration_pts));
1320 
1321  fe_spring_lhs_ptr_dX->getOpPtrVector().push_back(
1322  new OpSpringALEMaterialLhs_dX_dx(mesh_nodals_positions, field_name,
1323  data_at_integration_pts,
1324  feMatSideLhs_dx, side_fe_name));
1325 
1326  boost::shared_ptr<VolumeElementForcesAndSourcesCoreOnSide> feMatSideLhs_dX =
1327  boost::make_shared<VolumeElementForcesAndSourcesCoreOnSide>(m_field);
1328 
1329  feMatSideLhs_dX->getOpPtrVector().push_back(
1331  mesh_nodals_positions, data_at_integration_pts->HMat));
1332  feMatSideLhs_dX->getOpPtrVector().push_back(
1334  field_name, data_at_integration_pts->hMat));
1335 
1336  feMatSideLhs_dX->getOpPtrVector().push_back(new OpCalculateDeformation(
1337  mesh_nodals_positions, data_at_integration_pts));
1338 
1339  feMatSideLhs_dX->getOpPtrVector().push_back(
1340  new SpringALEMaterialVolOnSideLhs_dX_dX(mesh_nodals_positions,
1341  mesh_nodals_positions,
1342  data_at_integration_pts));
1343 
1344  fe_spring_lhs_ptr_dX->getOpPtrVector().push_back(
1345  new OpSpringALEMaterialLhs_dX_dX(
1346  mesh_nodals_positions, mesh_nodals_positions,
1347  data_at_integration_pts, feMatSideLhs_dX, side_fe_name));
1348  }
1349 
1351 }

The documentation for this struct was generated from the following files:
MoFEM::addHOOpsFace3D
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
Definition: HODataOperators.hpp:699
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
MoFEM::OpCalculateVectorFieldValues
Get values at integration pts for tensor filed rank 1, i.e. vector field.
Definition: UserDataOperators.hpp:466
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
bit
auto bit
set bit
Definition: hanging_node_approx.cpp:75
MoFEM::CoreInterface::check_field
virtual bool check_field(const std::string &name) const =0
check if field is in database
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1535
MF_ZERO
@ MF_ZERO
Definition: definitions.h:111
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359