v0.14.0
ContactPrismElementForcesAndSourcesCore.hpp
Go to the documentation of this file.
1 /** \file ContactPrismElementForcesAndSourcesCore.hpp
2  \brief Implementation of the contact prism element.
3 
4  These elements are used to enforce contact constraints in the interface
5  between two solids.
6 
7 */
8 
9 
10 
11 #ifndef __CONTACTPRISMELEMENTFORCESANDSURCESCORE_HPP__
12 #define __CONTACTPRISMELEMENTFORCESANDSURCESCORE_HPP__
13 
14 namespace MoFEM {
15 
16 struct VolumeElementForcesAndSourcesCoreOnContactPrismSide;
17 
18 /** \brief ContactPrism finite element
19  \ingroup mofem_forces_and_sources_prism_element
20 
21  User is implementing own operator at Gauss points level, by own class
22  derived from ContactPrismElementForcesAndSourcesCoreL::UserDataOperator.
23  Arbitrary number of operator added pushing instances to rowOpPtrVector and
24  rowColOpPtrVector.
25 
26  */
28 
30 
31  /** \brief default operator for Contact Prism element
32  * \ingroup mofem_forces_and_sources_prism_element
33  */
34  struct UserDataOperator;
35 
37 
38  inline const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
40  return dataOnMaster;
41  }
42 
43  inline const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
45  return dataOnSlave;
46  }
47 
49 
51 
52 protected:
53  std::array<double, 2> aRea; ///< Array storing master and slave faces areas
54 
56  normal; ///< vector storing vector normal to master or slave element
58  MatrixDouble coordsAtGaussPtsMaster; ///< matrix storing master Gauss points
59  ///< global coordinates
60  MatrixDouble coordsAtGaussPtsSlave; ///< matrix storing slave Gauss points
61  ///< global coordinates
62 
63  MatrixDouble gaussPtsMaster; ///< matrix storing master Gauss points local
64  ///< coordinates and weights
65  MatrixDouble gaussPtsSlave; ///< matrix storing slave Gauss points local
66  ///< coordinates and weights
67 
71 
72  /**
73  * @brief Entity data on element entity rows fields
74  *
75  *
76  * FIXME: that should be moved to private class data and acessed only by
77  * member function
78  */
79  const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
81  const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE> dataOnSlave;
82 
83  /**
84  * @brief Entity data on element entity columns fields
85  *
86  * FIXME: that should be moved to private class data and acessed only by
87  * member function
88  */
89  const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
91  const std::array<boost::shared_ptr<EntitiesFieldData>, LASTSPACE>
93 
96 
105 
106  MoFEMErrorCode setDefaultGaussPts(const int rule);
107 
108  /**
109  * @brief Iterate user data operators
110  *
111  * @return MoFEMErrorCode
112  */
114 
115  /**
116  * @brief Iterate user data operators
117  *
118  * @return MoFEMErrorCode
119  */
121  FieldApproximationBase m_s_base,
122  EntitiesFieldData &m_s_data);
123 
124  /** \brief function that gets entity field data.
125  *
126  * \param master_data data fot master face
127  * \param slave_data data fot master face
128  * \param field_name field name of interest
129  * \param type_lo lowest dimension entity type to be searched
130  * \param type_hi highest dimension entity type to be searched
131  */
133  EntitiesFieldData &master_data, EntitiesFieldData &slave_data,
134  const std::string &field_name, const EntityType type_lo = MBVERTEX,
135  const EntityType type_hi = MBPOLYHEDRON) const;
136 
137  /** \brief function that gets entity indices.
138  *
139  * \param master_data data fot master face
140  * \param slave_data data fot master face
141  * \param field_name field name of interest
142  * \param dofs MultiIndex container keeping FENumeredDofEntity.
143  * \param type_lo lowest dimension entity type to be searched
144  * \param type_hi highest dimension entity type to be searched
145  */
146  template <typename EXTRACTOR>
148  getEntityIndices(EntitiesFieldData &master_data,
149  EntitiesFieldData &slave_data, const std::string &field_name,
150  FieldEntity_vector_view &ents_field,
151  const EntityType type_lo, const EntityType type_hi,
152  EXTRACTOR &&extractor) const;
153 
154  /** \brief function that gets nodes indices.
155  *
156  * \param field_name field name of interest
157  * \param dofs MultiIndex container keeping FENumeredDofEntity.
158  * \param master_nodes_indices vector containing global master nodes indices
159  * \param master_local_nodes_indices vector containing local master nodes
160  * indices
161  * \param slave_nodes_indices vector containing global master nodes indices
162  * \param slave_local_nodes_indices vector containing local master nodes
163  * indices
164  */
165  template <typename EXTRACTOR>
166  MoFEMErrorCode getNodesIndices(const std::string field_name,
167  FieldEntity_vector_view &ents_field,
168  VectorInt &master_nodes_indices,
169  VectorInt &master_local_nodes_indices,
170  VectorInt &slave_nodes_indices,
171  VectorInt &slave_local_nodes_indices,
172  EXTRACTOR &&extractor) const;
173 
174  /** \brief function that gets nodes field data.
175  *
176  * \param field_name field name of interest
177  * \param dofs MultiIndex container keeping FENumeredDofEntity.
178  * \param master_nodes_data vector containing master nodes data
179  * \param slave_nodes_data vector containing master nodes data
180  * \param master_nodes_dofs vector containing master nodes dofs
181  * \param slave_nodes_dofs vector containing slave nodes dofs
182  * \param master_space approximation energy space at master
183  * \param slave_space approximation energy space at slave
184  * \param master_base base for master face
185  * \param slave_base base for slave face
186  */
187  MoFEMErrorCode getNodesFieldData(const std::string field_name,
188  VectorDouble &master_nodes_data,
189  VectorDouble &slave_nodes_data,
190  VectorDofs &master_nodes_dofs,
191  VectorDofs &slave_nodes_dofs,
192 
193  VectorFieldEntities &master_field_entities,
194  VectorFieldEntities &slave_field_entities,
195 
196  FieldSpace &master_space,
197  FieldSpace &slave_space,
198  FieldApproximationBase &master_base,
199  FieldApproximationBase &slave_base) const;
200 
201 private:
203 };
204 
205 /** \brief default operator for Contact Prism element
206  * \ingroup mofem_forces_and_sources_prism_element
207  */
210 
213 
214  UserDataOperator(const std::string &field_name, const char type)
216 
217  UserDataOperator(const std::string &row_field_name,
218  const std::string &col_field_name, const char type)
219  : ForcesAndSourcesCore::UserDataOperator(row_field_name, col_field_name,
220  type) {}
221 
222  UserDataOperator(const std::string &row_field_name,
223  const std::string &col_field_name, const char type,
224  const char face_type)
225  : ForcesAndSourcesCore::UserDataOperator(row_field_name, col_field_name,
226  type),
227  faceType(face_type) {}
228 
229  UserDataOperator(const std::string &field_name, const char type,
230  const char face_type)
232  faceType(face_type) {}
233 
234  /**
235  * \brief Controls loop over faces and face combination on element
236  *
237  * FACEMASTER is used when column or row data needs to be accessed
238  located at master face
239  * FACESLAVE is used when column or row data needs to be accessed
240  located at slave face
241  * FACEMASTERMASTER is used for accessing simultaneously row and col
242  data located at master face.
243  * FACEMASTERSLAVE is used for accessing simultaneously row data that
244  is located on master face and col data located at slave face.
245  * FACESLAVEMASTER is used for accessing simultaneously row data that
246  is located on slave face and col data located at master face.
247  * FACESLAVESLAVE is used for accessing simultaneously row and col
248  data located at slave face.
249  *
250  */
251  enum FaceType {
252  FACEMASTER = 1 << 0,
253  FACESLAVE = 1 << 1,
255  FACEMASTERSLAVE = 1 << 3,
256  FACESLAVEMASTER = 1 << 4,
257  FACESLAVESLAVE = 1 << 5,
258  FACELAST = 1 << 6
259  };
260 
261  char faceType;
262 
263  /**
264  * \brief Get operator types
265  * @return Return operator type
266  */
267  inline int getFaceType() const;
268 
269  inline boost::shared_ptr<const NumeredEntFiniteElement>
271 
272  /** \brief get face aRea Master
273  */
274  inline double getAreaMaster();
275 
276  /** \brief get face aRea Slave
277  */
278  inline double getAreaSlave();
279 
280  /** \brief get face normal vector to Master face
281  */
283 
284  /** \brief get first face tangent vector to Master face
285  */
287 
288  /** \brief get second face tangent vector to Master face
289  */
291 
292  /** \brief get face normal vector to Slave face
293  */
294  inline VectorAdaptor getNormalSlave();
295 
296  /** \brief get first face tangent vector to Slave face
297  */
299 
300  /** \brief get second face tangent vector to Slave face
301  */
303 
304  /** \brief get Gauss point at Master face
305  */
307 
308  /** \brief get Gauss point at Slave face
309  */
310  inline MatrixDouble &getGaussPtsSlave();
311 
312  /**
313  * @brief Get integration weights for slave side
314  *
315  * \code
316  * auto t_w = getFTensor0IntegrationWeight();
317  * for(int gg = 0; gg!=getGaussPts.size2(); ++gg) {
318  * // integrate something
319  * ++t_w;
320  * }
321  * \endcode
322  *
323  * @return FTensor::Tensor0<FTensor::PackPtr<double *, 1>>
324  */
325  inline auto getFTensor0IntegrationWeightSlave();
326 
327  /**
328  * @brief Get integration weights for master side
329  *
330  * \code
331  * auto t_w = getFTensor0IntegrationWeight();
332  * for(int gg = 0; gg!=getGaussPts.size2(); ++gg) {
333  * // integrate something
334  * ++t_w;
335  * }
336  * \endcode
337  *
338  * @return FTensor::Tensor0<FTensor::PackPtr<double *, 1>>
339  */
341 
342  /** \brief get triangle coordinates
343 
344  Vector has 9 elements, i.e. coordinates on Master face
345 
346  */
347  inline VectorDouble getCoordsMaster();
348 
349  /** \brief get triangle coordinates
350 
351  Vector has 9 elements, i.e. coordinates on Slave face
352 
353  */
354  inline VectorDouble getCoordsSlave();
355 
356  /** \brief get coordinates at Gauss pts on full prism.
357 
358  Matrix has size (nb integration points on master)x(3),
359  i.e. coordinates on face Master
360 
361  */
363 
364  /** \brief get coordinates at Gauss pts on full prism.
365 
366  Matrix has size (nb integration points on slave)x(3),
367  i.e. coordinates on face Slave
368 
369  */
371 
372  /** \brief return pointer to triangle finite element object
373  */
376 
377  /**
378  *
379  * User call this function to loop over elements on the side of face. This
380  * function calls MoFEM::VolumeElementForcesAndSourcesCoreOnContactPrismSide
381  * with is operator to do calculations.
382  *
383  * @param fe_name Name of the element
384  * @param method Finite element object
385  * @param side_type states the side from which side element will work (0
386  * for master 1 for slave)
387  * @return error code
388  */
390  const string fe_name,
392  const int side_type, const EntityHandle ent_for_side);
393 
394 protected:
395  inline ForcesAndSourcesCore *getSidePtrFE() const;
396 };
397 
398 boost::shared_ptr<const NumeredEntFiniteElement>
401  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
403 };
404 
405 inline int
407  return faceType;
408 }
409 
410 inline double
412  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
413 }
414 
415 inline double
417  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[1];
418 }
419 
420 inline VectorAdaptor
422  double *data = &(
423  static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->normal[0]);
424  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
425 }
426 
429  double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
430  ->tangentMasterOne[0]);
431  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
432 }
433 
436  double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
437  ->tangentMasterTwo[0]);
438  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
439 }
440 
441 inline VectorAdaptor
443  double *data = &(
444  static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->normal[3]);
445  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
446 }
447 
450  double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
451  ->tangentSlaveOne[0]);
452  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
453 }
454 
457  double *data = &(static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
458  ->tangentSlaveTwo[0]);
459  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
460 }
461 
462 inline MatrixDouble &
464  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
465  ->gaussPtsMaster;
466 }
467 
468 inline MatrixDouble &
470  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
471  ->gaussPtsSlave;
472 }
473 
477  &(getGaussPtsSlave()(getGaussPtsSlave().size1() - 1, 0)));
478 }
479 
483  &(getGaussPtsMaster()(getGaussPtsMaster().size1() - 1, 0)));
484 }
485 
486 inline VectorDouble
488  double *data = &(
489  static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->coords[0]);
490  return VectorAdaptor(9, ublas::shallow_array_adaptor<double>(9, data));
491 }
492 
493 inline VectorDouble
495  double *data = &(
496  static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)->coords[9]);
497  return VectorAdaptor(9, ublas::shallow_array_adaptor<double>(9, data));
498 }
499 
502  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
504 }
505 
508  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE)
510 }
511 
515  return static_cast<ContactPrismElementForcesAndSourcesCore *>(ptrFE);
516 }
517 
518 } // namespace MoFEM
519 
520 #endif //__CONTACTPRISMELEMENTFORCESANDSURCESCORE_HPP__
UBlasMatrix< double >
MoFEM::ContactPrismElementForcesAndSourcesCore::tangentSlaveOne
VectorDouble tangentSlaveOne
Definition: ContactPrismElementForcesAndSourcesCore.hpp:68
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPtsMaster
MatrixDouble & getGaussPtsMaster()
get Gauss point at Master face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:463
MoFEM::ContactPrismElementForcesAndSourcesCore::tangentMasterOne
VectorDouble tangentMasterOne
Definition: ContactPrismElementForcesAndSourcesCore.hpp:69
MoFEM::ContactPrismElementForcesAndSourcesCore::derivedDataOnSlave
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > derivedDataOnSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:92
MoFEM::ContactPrismElementForcesAndSourcesCore::dataHcurlMaster
EntitiesFieldData & dataHcurlMaster
Definition: ContactPrismElementForcesAndSourcesCore.hpp:99
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACEMASTERSLAVE
@ FACEMASTERSLAVE
Definition: ContactPrismElementForcesAndSourcesCore.hpp:255
MoFEM::ContactPrismElementForcesAndSourcesCore::dataL2Slave
EntitiesFieldData & dataL2Slave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:104
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACESLAVESLAVE
@ FACESLAVESLAVE
Definition: ContactPrismElementForcesAndSourcesCore.hpp:257
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator
UserDataOperator(const std::string &field_name, const char type, const char face_type)
Definition: ContactPrismElementForcesAndSourcesCore.hpp:229
MoFEM::ContactPrismElementForcesAndSourcesCore::opContravariantTransform
OpSetContravariantPiolaTransformOnFace opContravariantTransform
Definition: ContactPrismElementForcesAndSourcesCore.hpp:70
EntityHandle
MoFEM::ContactPrismElementForcesAndSourcesCore::derivedDataOnMaster
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > derivedDataOnMaster
Entity data on element entity columns fields.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:90
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getNormalSlave
VectorAdaptor getNormalSlave()
get face normal vector to Slave face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:442
MoFEM::ContactPrismElementForcesAndSourcesCore::ContactPrismElementForcesAndSourcesCore
ContactPrismElementForcesAndSourcesCore(Interface &m_field)
Definition: ContactPrismElementForcesAndSourcesCore.cpp:10
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::FieldEntity_vector_view
std::vector< boost::weak_ptr< FieldEntity > > FieldEntity_vector_view
Definition: FieldEntsMultiIndices.hpp:478
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator
UserDataOperator(const std::string &row_field_name, const std::string &col_field_name, const char type, const char face_type)
Definition: ContactPrismElementForcesAndSourcesCore.hpp:222
MoFEM::ContactPrismElementForcesAndSourcesCore::tangentMasterTwo
VectorDouble tangentMasterTwo
Definition: ContactPrismElementForcesAndSourcesCore.hpp:69
MoFEM::ContactPrismElementForcesAndSourcesCore::getGaussPtsMasterFromEleSide
MatrixDouble & getGaussPtsMasterFromEleSide()
Definition: ContactPrismElementForcesAndSourcesCore.hpp:48
MoFEM::ContactPrismElementForcesAndSourcesCore::dataHdivMaster
EntitiesFieldData & dataHdivMaster
Definition: ContactPrismElementForcesAndSourcesCore.hpp:101
MoFEM::ContactPrismElementForcesAndSourcesCore::getEntityFieldData
MoFEMErrorCode getEntityFieldData(EntitiesFieldData &master_data, EntitiesFieldData &slave_data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
function that gets entity field data.
Definition: ContactPrismElementForcesAndSourcesCore.cpp:784
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getNumeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Definition: ContactPrismElementForcesAndSourcesCore.hpp:400
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentSlaveOne
VectorAdaptor getTangentSlaveOne()
get first face tangent vector to Slave face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:449
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentMasterOne
VectorAdaptor getTangentMasterOne()
get first face tangent vector to Master face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:428
MoFEM::ContactPrismElementForcesAndSourcesCore::getValueHdivDemkowiczBase
MoFEMErrorCode getValueHdivDemkowiczBase(MatrixDouble &pts, FieldApproximationBase m_s_base, EntitiesFieldData &m_s_data)
Iterate user data operators.
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeightMaster
auto getFTensor0IntegrationWeightMaster()
Get integration weights for master side.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:481
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::ContactPrismElementForcesAndSourcesCore::getDataOnMasterFromEleSide
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > getDataOnMasterFromEleSide()
Definition: ContactPrismElementForcesAndSourcesCore.hpp:39
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getFaceType
int getFaceType() const
Get operator types.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:406
MoFEM::ContactPrismElementForcesAndSourcesCore::dataNoFieldSlave
EntitiesFieldData & dataNoFieldSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:98
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator
UserDataOperator(const std::string &row_field_name, const std::string &col_field_name, const char type)
Definition: ContactPrismElementForcesAndSourcesCore.hpp:217
MoFEM::ContactPrismElementForcesAndSourcesCore::getDataOnSlaveFromEleSide
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > getDataOnSlaveFromEleSide()
Definition: ContactPrismElementForcesAndSourcesCore.hpp:44
MoFEM::ContactPrismElementForcesAndSourcesCore::dataOnSlave
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:81
FieldSpace
FieldSpace
approximation spaces
Definition: definitions.h:82
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::loopSideVolumes
MoFEMErrorCode loopSideVolumes(const string fe_name, VolumeElementForcesAndSourcesCoreOnContactPrismSide &fe_method, const int side_type, const EntityHandle ent_for_side)
Definition: ContactPrismElementForcesAndSourcesCore.cpp:1181
MoFEM::ContactPrismElementForcesAndSourcesCore
ContactPrism finite element.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:27
MoFEM::ContactPrismElementForcesAndSourcesCore::getNodesIndices
MoFEMErrorCode getNodesIndices(const std::string field_name, FieldEntity_vector_view &ents_field, VectorInt &master_nodes_indices, VectorInt &master_local_nodes_indices, VectorInt &slave_nodes_indices, VectorInt &slave_local_nodes_indices, EXTRACTOR &&extractor) const
function that gets nodes indices.
Definition: ContactPrismElementForcesAndSourcesCore.cpp:1081
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACESLAVEMASTER
@ FACESLAVEMASTER
Definition: ContactPrismElementForcesAndSourcesCore.hpp:256
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentMasterTwo
VectorAdaptor getTangentMasterTwo()
get second face tangent vector to Master face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:435
MoFEM::ContactPrismElementForcesAndSourcesCore::setDefaultGaussPts
MoFEMErrorCode setDefaultGaussPts(const int rule)
Definition: ContactPrismElementForcesAndSourcesCore.cpp:92
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::VolumeElementForcesAndSourcesCoreOnContactPrismSide
Base volume element used to integrate on contact surface (could be extended to other volume elements)
Definition: VolumeElementForcesAndSourcesCoreOnContactPrismSide.hpp:23
MoFEM::ContactPrismElementForcesAndSourcesCore::getNodesFieldData
MoFEMErrorCode getNodesFieldData(const std::string field_name, VectorDouble &master_nodes_data, VectorDouble &slave_nodes_data, VectorDofs &master_nodes_dofs, VectorDofs &slave_nodes_dofs, VectorFieldEntities &master_field_entities, VectorFieldEntities &slave_field_entities, FieldSpace &master_space, FieldSpace &slave_space, FieldApproximationBase &master_base, FieldApproximationBase &slave_base) const
function that gets nodes field data.
Definition: ContactPrismElementForcesAndSourcesCore.cpp:880
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPtsSlave
MatrixDouble & getCoordsAtGaussPtsSlave()
get coordinates at Gauss pts on full prism.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:507
MoFEM::ContactPrismElementForcesAndSourcesCore::normal
VectorDouble normal
vector storing vector normal to master or slave element
Definition: ContactPrismElementForcesAndSourcesCore.hpp:56
MoFEM::ContactPrismElementForcesAndSourcesCore::dataL2Master
EntitiesFieldData & dataL2Master
Definition: ContactPrismElementForcesAndSourcesCore.hpp:103
convert.type
type
Definition: convert.py:64
MoFEM::ContactPrismElementForcesAndSourcesCore::tangentSlaveTwo
VectorDouble tangentSlaveTwo
Definition: ContactPrismElementForcesAndSourcesCore.hpp:68
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getAreaMaster
double getAreaMaster()
get face aRea Master
Definition: ContactPrismElementForcesAndSourcesCore.hpp:411
MoFEM::VectorDofs
ublas::vector< FEDofEntity *, DofsAllocator > VectorDofs
Definition: EntitiesFieldData.hpp:23
MoFEM::ContactPrismElementForcesAndSourcesCore::gaussPtsMaster
MatrixDouble gaussPtsMaster
Definition: ContactPrismElementForcesAndSourcesCore.hpp:63
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator
default operator for Contact Prism element
Definition: ContactPrismElementForcesAndSourcesCore.hpp:208
MoFEM::ContactPrismElementForcesAndSourcesCore::coordsAtGaussPtsSlave
MatrixDouble coordsAtGaussPtsSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:60
LASTSPACE
@ LASTSPACE
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:89
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator
UserDataOperator(const std::string &field_name, const char type)
Definition: ContactPrismElementForcesAndSourcesCore.hpp:214
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getTangentSlaveTwo
VectorAdaptor getTangentSlaveTwo()
get second face tangent vector to Slave face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:456
MoFEM::ContactPrismElementForcesAndSourcesCore::dataH1Slave
EntitiesFieldData & dataH1Slave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:95
MoFEM::ContactPrismElementForcesAndSourcesCore::getEntityIndices
MoFEMErrorCode getEntityIndices(EntitiesFieldData &master_data, EntitiesFieldData &slave_data, const std::string &field_name, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
function that gets entity indices.
Definition: ContactPrismElementForcesAndSourcesCore.cpp:987
MoFEM::ContactPrismElementForcesAndSourcesCore::coords
VectorDouble coords
Definition: ContactPrismElementForcesAndSourcesCore.hpp:57
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::UserDataOperator
UserDataOperator(const FieldSpace space)
Definition: ContactPrismElementForcesAndSourcesCore.hpp:211
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getSidePtrFE
ForcesAndSourcesCore * getSidePtrFE() const
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getNormalMaster
VectorAdaptor getNormalMaster()
get face normal vector to Master face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:421
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsSlave
VectorDouble getCoordsSlave()
get triangle coordinates
Definition: ContactPrismElementForcesAndSourcesCore.hpp:494
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getAreaSlave
double getAreaSlave()
get face aRea Slave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:416
MoFEM::ContactPrismElementForcesAndSourcesCore::dataNoFieldMaster
EntitiesFieldData & dataNoFieldMaster
Definition: ContactPrismElementForcesAndSourcesCore.hpp:97
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MoFEM::ContactPrismElementForcesAndSourcesCore::dataHcurlSlave
EntitiesFieldData & dataHcurlSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:100
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getContactPrismElementForcesAndSourcesCore
const ContactPrismElementForcesAndSourcesCore * getContactPrismElementForcesAndSourcesCore()
return pointer to triangle finite element object
Definition: ContactPrismElementForcesAndSourcesCore.hpp:514
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
FTensor::Tensor0
Definition: Tensor0.hpp:16
MoFEM::ContactPrismElementForcesAndSourcesCore::nbGaussPts
int nbGaussPts
Definition: ContactPrismElementForcesAndSourcesCore.hpp:202
MoFEM::OpSetContravariantPiolaTransformOnFace
transform Hdiv base fluxes from reference element to physical triangle
Definition: DataOperators.hpp:416
MoFEM::Types::VectorInt
UBlasVector< int > VectorInt
Definition: Types.hpp:67
MoFEM::ContactPrismElementForcesAndSourcesCore::dataH1Master
EntitiesFieldData & dataH1Master
Definition: ContactPrismElementForcesAndSourcesCore.hpp:94
MoFEM::ContactPrismElementForcesAndSourcesCore::operator()
MoFEMErrorCode operator()()
function is run for every finite element
Definition: ContactPrismElementForcesAndSourcesCore.cpp:140
MoFEM::ContactPrismElementForcesAndSourcesCore::coordsAtGaussPtsMaster
MatrixDouble coordsAtGaussPtsMaster
Definition: ContactPrismElementForcesAndSourcesCore.hpp:58
MoFEM::ContactPrismElementForcesAndSourcesCore::dataOnMaster
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnMaster
Entity data on element entity rows fields.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:80
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsMaster
VectorDouble getCoordsMaster()
get triangle coordinates
Definition: ContactPrismElementForcesAndSourcesCore.hpp:487
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPtsSlave
MatrixDouble & getGaussPtsSlave()
get Gauss point at Slave face
Definition: ContactPrismElementForcesAndSourcesCore.hpp:469
UBlasVector< double >
FieldApproximationBase
FieldApproximationBase
approximation base
Definition: definitions.h:58
MoFEM::ContactPrismElementForcesAndSourcesCore::dataHdivSlave
EntitiesFieldData & dataHdivSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:102
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPtsMaster
MatrixDouble & getCoordsAtGaussPtsMaster()
get coordinates at Gauss pts on full prism.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:501
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACELAST
@ FACELAST
Definition: ContactPrismElementForcesAndSourcesCore.hpp:258
MoFEM::VectorFieldEntities
ublas::vector< FieldEntity *, FieldEntAllocator > VectorFieldEntities
Definition: EntitiesFieldData.hpp:31
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FaceType
FaceType
Definition: ContactPrismElementForcesAndSourcesCore.hpp:251
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACEMASTERMASTER
@ FACEMASTERMASTER
Definition: ContactPrismElementForcesAndSourcesCore.hpp:254
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::faceType
char faceType
Definition: ContactPrismElementForcesAndSourcesCore.hpp:261
MoFEM::EntitiesFieldData
data structure for finite element entity
Definition: EntitiesFieldData.hpp:40
MoFEM::ContactPrismElementForcesAndSourcesCore::loopOverOperators
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
Definition: ContactPrismElementForcesAndSourcesCore.cpp:482
MoFEM::ContactPrismElementForcesAndSourcesCore::aRea
std::array< double, 2 > aRea
Array storing master and slave faces areas.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:53
MoFEM::ForcesAndSourcesCore::UserDataOperator::ptrFE
ForcesAndSourcesCore * ptrFE
Definition: ForcesAndSourcesCore.hpp:985
MoFEM::ContactPrismElementForcesAndSourcesCore::getGaussPtsSlaveFromEleSide
MatrixDouble & getGaussPtsSlaveFromEleSide()
Definition: ContactPrismElementForcesAndSourcesCore.hpp:50
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACESLAVE
@ FACESLAVE
Definition: ContactPrismElementForcesAndSourcesCore.hpp:253
MoFEM::ContactPrismElementForcesAndSourcesCore::gaussPtsSlave
MatrixDouble gaussPtsSlave
Definition: ContactPrismElementForcesAndSourcesCore.hpp:65
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::FACEMASTER
@ FACEMASTER
Definition: ContactPrismElementForcesAndSourcesCore.hpp:252
MoFEM::ContactPrismElementForcesAndSourcesCore::UserDataOperator::getFTensor0IntegrationWeightSlave
auto getFTensor0IntegrationWeightSlave()
Get integration weights for slave side.
Definition: ContactPrismElementForcesAndSourcesCore.hpp:475