v0.8.15
Classes | Public Member Functions | Public Attributes | List of all members
MoFEM::FaceElementForcesAndSourcesCore Struct Reference

Face finite elementUser is implementing own operator at Gauss point level, by own object derived from FaceElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to OpPtrVector. More...

#include <src/finite_elements/FaceElementForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::FaceElementForcesAndSourcesCore:
[legend]
Collaboration diagram for MoFEM::FaceElementForcesAndSourcesCore:
[legend]

Classes

struct  UserDataOperator
 default operator for TRI element More...
 

Public Member Functions

 FaceElementForcesAndSourcesCore (Interface &m_field)
 
virtual MoFEMErrorCode calculateAreaAndNormal ()
 Calculate element area and normal of the face. More...
 
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
DEPRECATED MoFEMErrorCode setIntegartionPts ()
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode calculateHoNormal ()
 Calculate normal on curved elements. More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
virtual ~ForcesAndSourcesCore ()
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 
int getMaxDataOrder () const
 Get max order of approximation for data fields. More...
 
int getMaxRowOrder () const
 Get max order of approximation for field in rows. More...
 
int getMaxColOrder () const
 Get max order of approximation for field in columns. More...
 
MoFEMErrorCode getSense (EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order of approximation on the entity More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (DataForcesAndSourcesCore &data) const
 
template<EntityType type>
MoFEMErrorCode getEntityFieldDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getNodesIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices, VectorInt &local_nodes_indices) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEntityIndices (DataForcesAndSourcesCore &data, const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const std::string &field_name, FENumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get col NoField indices More...
 
MoFEMErrorCode getNodesFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 
MoFEMErrorCode getNoFieldFieldData (DataForcesAndSourcesCore &data, const boost::string_ref field_name) const
 
MoFEMErrorCode getNodesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &data) const
 Get field approximation space and base on entities. More...
 
MoFEMErrorCode getProblemNodesIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, VectorInt &nodes_indices) const
 get indices of nodal indices which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemTypeIndices (const std::string &field_name, const NumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices) const
 get indices by type (generic function) which are declared for problem but not this particular element More...
 
MoFEMErrorCode getProblemNodesRowIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeRowIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
MoFEMErrorCode getProblemNodesColIndices (const std::string &field_name, VectorInt &nodes_indices) const
 
MoFEMErrorCode getProblemTypeColIndices (const std::string &field_name, EntityType type, int side_number, VectorInt &indices) const
 
virtual int getRule (int order)
 set integration rule for finite element More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order)
 It will be removed in the future use other variant. More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
boost::ptr_vector< UserDataOperator > & getOpPtrVector ()
 Use to push back operator for row operator. More...
 
auto & getElementPolynomialBase ()
 Get the Entity Polynomial Base object. More...
 
auto & getUserPolynomialBase ()
 Get the User Polynomial Base object. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calculateBaseFunctionsOnElement ()
 Calculate base functions. More...
 
MoFEMErrorCode createDataOnElement ()
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type, const int side_number) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityType type) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_begin (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
template<class MULTIINDEX >
MULTIINDEX::iterator get_end (const MULTIINDEX &index, const std::string &field_name, const EntityHandle ent) const
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()
 
MoFEMErrorCode setKspCtx (const KSPContext &ctx)
 set operator type More...
 
MoFEMErrorCode setKsp (KSP ksp)
 set solver More...
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (const MOFEMuuid &uuid, bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE , bool VERIFY = false>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, IFACE *&iface) const
 Get interface by uuid and return reference to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()
 
virtual MoFEMErrorCode getLibVersion (Version &version) const
 Get library version. More...
 
virtual const MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version) const
 Get database major version. More...
 
virtual MoFEMErrorCode getInterfaceVersion (Version &version) const
 Get database major version. More...
 
template<>
MoFEMErrorCode getInterface (const MOFEMuuid &uuid, UnknownInterface *&iface) const
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()
 
MoFEMErrorCode setSnesCtx (const SNESContext &ctx)
 Set SNES context. More...
 
MoFEMErrorCode setSnes (SNES snes)
 Set SNES instance. More...
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()
 
MoFEMErrorCode setTsCtx (const TSContext &ctx)
 Set Ts context. More...
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
MoFEMErrorCode setTs (TS _ts)
 Set TS solver. More...
 

Public Attributes

double aRea
 
int num_nodes
 
const EntityHandleconn
 
VectorDouble nOrmal
 
VectorDouble tangentOne
 
VectorDouble tangentTwo
 
VectorDouble coords
 
MatrixDouble coordsAtGaussPts
 
std::string meshPositionsFieldName
 Name of the field with geometry. More...
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble normalsAtGaussPts
 
MatrixDouble tangentOneAtGaussPts
 
MatrixDouble tangentTwoAtGaussPts
 
OpGetCoordsAndNormalsOnFace opHOCoordsAndNormals
 
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform
 
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
 
int nbGaussPts
 Number of integration points. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
const boost::shared_ptr< DataForcesAndSourcesCoredataOnElement [LASTSPACE]
 Entity data on element entity rows fields. More...
 
const boost::shared_ptr< DataForcesAndSourcesCorederivedDataOnElement [LASTSPACE]
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 
int loopSize
 
int rAnk
 
int sIze
 
const RefEntity_multiIndexrefinedEntitiesPtr
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 
const ProblemproblemPtr
 
const Field_multiIndexfieldsPtr
 
const FieldEntity_multiIndexentitiesPtr
 
const DofEntity_multiIndexdofsPtr
 
const FiniteElement_multiIndexfiniteElementsPtr
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 
boost::function< MoFEMErrorCode()> postProcessHook
 
boost::function< MoFEMErrorCode()> operatorHook
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec ksp_f
 the right hand side vector More...
 
Mat ksp_A
 matrix More...
 
Mat ksp_B
 preconditioner matrix More...
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 
Vec snes_x
 
Vec snes_f
 
Mat snes_A
 
Mat snes_B
 
- Public Attributes inherited from MoFEM::TSMethod
TSContext ts_ctx
 
TS ts
 
Vec ts_u
 
Vec ts_u_t
 
Vec ts_F
 
Mat ts_A
 
Mat ts_B
 
PetscInt ts_step
 
PetscReal ts_a
 
PetscReal ts_t
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION, CTX_OPERATORS, CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION, CTX_SNESSETJACOBIAN, CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION, CTX_TSSETRHSJACOBIAN, CTX_TSSETIFUNCTION, CTX_TSSETIJACOBIAN,
  CTX_TSTSMONITORSET, CTX_TSNONE
}
 
- Protected Member Functions inherited from MoFEM::UnknownInterface
boost::typeindex::type_index getClassIdx (const MOFEMuuid &uid) const
 Get type name for interface Id. More...
 
MOFEMuuid getUId (const boost::typeindex::type_index &class_idx) const
 Get interface Id for class name. More...
 

Detailed Description

Face finite element

User is implementing own operator at Gauss point level, by own object derived from FaceElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to OpPtrVector.

Examples:
hello_world.cpp, MagneticElement.hpp, simple_elasticity.cpp, and simple_interface.cpp.

Definition at line 38 of file FaceElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ FaceElementForcesAndSourcesCore()

FaceElementForcesAndSourcesCore::FaceElementForcesAndSourcesCore ( Interface m_field)

Definition at line 79 of file FaceElementForcesAndSourcesCore.cpp.

81  : ForcesAndSourcesCore(m_field),
82  meshPositionsFieldName("MESH_NODE_POSITIONS"),
90  boost::shared_ptr<BaseFunction>(new TriPolynomialBase());
91 }
ForcesAndSourcesCore(Interface &m_field)
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
std::string meshPositionsFieldName
Name of the field with geometry.
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

Member Function Documentation

◆ calculateAreaAndNormal()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateAreaAndNormal ( )
virtual

Calculate element area and normal of the face.

Note that at that point is assumed that geometry is exclusively defined by corner nodes.

Returns
Error code

Definition at line 143 of file FaceElementForcesAndSourcesCore.cpp.

143  {
146  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
147  coords.resize(num_nodes * 3, false);
148  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
149  double diff_n[6];
150  CHKERR ShapeDiffMBTRI(diff_n);
151 
152  nOrmal.resize(3, false);
153  tangentOne.resize(3, false);
154  tangentTwo.resize(3, false);
156  &coords[0], &coords[1], &coords[2]);
158  &nOrmal[0], &nOrmal[1], &nOrmal[2]);
160  &tangentOne[0], &tangentOne[1], &tangentOne[2]);
162  &tangentTwo[0], &tangentTwo[1], &tangentTwo[2]);
164  &diff_n[1]);
165 
171  t_t1(i) = 0;
172  t_t2(i) = 0;
173  for (int nn = 0; nn != 3; ++nn) {
174  t_t1(i) += t_coords(i) * t_diff(N0);
175  t_t2(i) += t_coords(i) * t_diff(N1);
176  ++t_coords;
177  ++t_diff;
178  }
179  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
180  aRea = sqrt(t_normal(i) * t_normal(i)) / 2.;
181 
183 }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:186
#define CHKERR
Inline error check.
Definition: definitions.h:578
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateCoordinatesAtGaussPts ( )
virtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 282 of file FaceElementForcesAndSourcesCore.cpp.

282  {
284 
285 
286  double *shape_functions =
287  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
288  coordsAtGaussPts.resize(nbGaussPts, 3, false);
289  for (int gg = 0; gg != nbGaussPts; ++gg) {
290  for (int dd = 0; dd != 3; ++dd) {
291  coordsAtGaussPts(gg, dd) =
292  cblas_ddot(3, &shape_functions[3 * gg], 1, &coords[dd], 3);
293  }
294  }
296 }
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T *> &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12

◆ calculateHoNormal()

MoFEMErrorCode FaceElementForcesAndSourcesCore::calculateHoNormal ( )
virtual

Calculate normal on curved elements.

Geometry is given by other field.

Returns
error code

Definition at line 299 of file FaceElementForcesAndSourcesCore.cpp.

299  {
301  // Check if field for high-order geometry is set and if it is set calculate
302  // higher-order normals and face tangent vectors.
303  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
304  dataPtr->get<FieldName_mi_tag>().end()) {
305 
306  const Field *field_struture =
308  BitFieldId id = field_struture->getId();
309 
310  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none()) {
311  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
312  "no MESH_NODE_POSITIONS in element data");
313  }
314 
315  // Calculate normal for high-order geometry
316 
322 
323  } else {
324  hoCoordsAtGaussPts.resize(0, 0, false);
325  normalsAtGaussPts.resize(0, 0, false);
326  tangentOneAtGaussPts.resize(0, 0, false);
327  tangentTwoAtGaussPts.resize(0, 0, false);
328  }
330 }
DataForcesAndSourcesCore & dataH1
MoFEMErrorCode getNodesFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
Get field data on nodes.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
std::bitset< BITFIELDID_SIZE > BitFieldId
Definition: Common.hpp:149
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
std::string meshPositionsFieldName
Name of the field with geometry.
#define CHKERR
Inline error check.
Definition: definitions.h:578
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
virtual const Field * get_field_structure(const std::string &name)=0
get field structure

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode FaceElementForcesAndSourcesCore::getSpaceBaseAndOrderOnElement ( )
virtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 236 of file FaceElementForcesAndSourcesCore.cpp.

236  {
238  // Get spaces order/base and sense of entities.
239 
241 
242  // H1
243  if (dataH1.spacesOnEntities[MBEDGE].test(H1)) {
244  CHKERR getEntitySense<MBEDGE>(dataH1);
245  CHKERR getEntityFieldDataOrder<MBEDGE>(dataH1, H1);
246  }
247  if (dataH1.spacesOnEntities[MBTRI].test(H1)) {
248  CHKERR getEntitySense<MBTRI>(dataH1);
249  CHKERR getEntityFieldDataOrder<MBTRI>(dataH1, H1);
250  }
251 
252  // Hcurl
253  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
254  CHKERR getEntitySense<MBEDGE>(dataHcurl);
255  CHKERR getEntityFieldDataOrder<MBEDGE>(dataHcurl, HCURL);
256  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
257  }
258  if (dataH1.spacesOnEntities[MBTRI].test(HCURL)) {
259  CHKERR getEntitySense<MBTRI>(dataHcurl);
260  CHKERR getEntityFieldDataOrder<MBTRI>(dataHcurl, HCURL);
261  dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
262  }
263 
264  // Hdiv
265  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
266  CHKERR getEntitySense<MBTRI>(dataHdiv);
267  CHKERR getEntityFieldDataOrder<MBTRI>(dataHdiv, HDIV);
268  dataHdiv.spacesOnEntities[MBTRI].set(HDIV);
269  }
270 
271  // L2
272  if (dataH1.spacesOnEntities[MBTRI].test(L2)) {
273  CHKERR getEntitySense<MBTRI>(dataL2);
274  CHKERR getEntityFieldDataOrder<MBTRI>(dataL2, L2);
275  dataL2.spacesOnEntities[MBTRI].set(L2);
276  }
277 
279 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:170
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:169
#define CHKERR
Inline error check.
Definition: definitions.h:578
DataForcesAndSourcesCore & dataHdiv
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
continuous field
Definition: definitions.h:168
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
field with C-1 continuity
Definition: definitions.h:171

◆ operator()()

MoFEMErrorCode FaceElementForcesAndSourcesCore::operator() ( )
virtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from MoFEM::ForcesAndSourcesCore.

Reimplemented in NeummanForcesSurfaceComplexForLazy::MyTriangleSpatialFE.

Definition at line 332 of file FaceElementForcesAndSourcesCore.cpp.

332  {
334 
335  if (numeredEntFiniteElementPtr->getEntType() != MBTRI)
338 
339  // Calculate normal and tangent vectors for face geometry given by 3 nodes.
342 
344  if (nbGaussPts == 0)
346 
347 
348  DataForcesAndSourcesCore &data_curl = *dataOnElement[HCURL];
349  DataForcesAndSourcesCore &data_div = *dataOnElement[HDIV];
350 
351  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
353  &*dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
354 
355  /// Use the some node base
359 
360  // Apply Piola transform to HDiv and HCurl spaces, uses previously calculated
361  // faces normal and tangent vectors.
362  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
364  }
365  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
366  CHKERR opCovariantTransform.opRhs(data_curl);
367  }
368 
369  // Iterate over operators
371 
373 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
field with continuous normal traction
Definition: definitions.h:170
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
DataForcesAndSourcesCore & dataH1
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
const boost::shared_ptr< DataForcesAndSourcesCore > dataOnElement[LASTSPACE]
Entity data on element entity rows fields.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
OpSetCovariantPiolaTransformOnTriangle opCovariantTransform
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:186
virtual MoFEMErrorCode calculateHoNormal()
Calculate normal on curved elements.
field with continuous tangents
Definition: definitions.h:169
#define CHKERR
Inline error check.
Definition: definitions.h:578
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
OpSetContravariantPiolaTransformOnTriangle opContravariantTransform
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.

◆ setIntegartionPts()

DEPRECATED MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCore::setIntegartionPts ( )
Deprecated:
method with spelling mistake, use setIntegrationPts

Definition at line 414 of file FaceElementForcesAndSourcesCore.hpp.

414 { return setIntegrationPts(); }
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.

◆ setIntegrationPts()

MoFEMErrorCode FaceElementForcesAndSourcesCore::setIntegrationPts ( )
virtual

Set integration points.

Returns
Error code

Definition at line 185 of file FaceElementForcesAndSourcesCore.cpp.

185  {
187  // Set integration points
188  int order_data = getMaxDataOrder();
189  int order_row = getMaxRowOrder();
190  int order_col = getMaxColOrder();
191  int rule = getRule(order_row, order_col, order_data);
192 
193  if (rule >= 0) {
194  if (rule < QUAD_2D_TABLE_SIZE) {
195  if (QUAD_2D_TABLE[rule]->dim != 2) {
196  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
197  }
198  if (QUAD_2D_TABLE[rule]->order < rule) {
199  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
200  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
201  }
203  gaussPts.resize(3, nbGaussPts, false);
204  cblas_dcopy(nbGaussPts, &QUAD_2D_TABLE[rule]->points[1], 3,
205  &gaussPts(0, 0), 1);
206  cblas_dcopy(nbGaussPts, &QUAD_2D_TABLE[rule]->points[2], 3,
207  &gaussPts(1, 0), 1);
208  cblas_dcopy(nbGaussPts, QUAD_2D_TABLE[rule]->weights, 1, &gaussPts(2, 0),
209  1);
210  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 3,
211  false);
212  double *shape_ptr =
213  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
214  cblas_dcopy(3 * nbGaussPts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr, 1);
215  } else {
216  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
217  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
218  nbGaussPts = 0;
219  }
220  } else {
221  // If rule is negative, set user defined integration points
222  CHKERR setGaussPts(order_row, order_col, order_data);
223  nbGaussPts = gaussPts.size2();
224  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nbGaussPts, 3,
225  false);
226  if (nbGaussPts) {
228  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
229  &gaussPts(0, 0), &gaussPts(1, 0), nbGaussPts);
230  }
231  }
233 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
MatrixDouble gaussPts
Matrix of integration points.
DataForcesAndSourcesCore & dataH1
int npoints
Definition: quad.h:29
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
virtual MoFEMErrorCode setGaussPts(int order)
It will be removed in the future use other variant.
#define CHKERR
Inline error check.
Definition: definitions.h:578
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:175
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
virtual int getRule(int order)
set integration rule for finite element
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.

Member Data Documentation

◆ aRea

double MoFEM::FaceElementForcesAndSourcesCore::aRea

Definition at line 40 of file FaceElementForcesAndSourcesCore.hpp.

◆ conn

const EntityHandle* MoFEM::FaceElementForcesAndSourcesCore::conn

Definition at line 43 of file FaceElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::coords

Definition at line 45 of file FaceElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::coordsAtGaussPts

Definition at line 46 of file FaceElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::hoCoordsAtGaussPts

Definition at line 50 of file FaceElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::FaceElementForcesAndSourcesCore::meshPositionsFieldName

Name of the field with geometry.

Definition at line 48 of file FaceElementForcesAndSourcesCore.hpp.

◆ nbGaussPts

int MoFEM::FaceElementForcesAndSourcesCore::nbGaussPts

Number of integration points.

Definition at line 405 of file FaceElementForcesAndSourcesCore.hpp.

◆ nOrmal

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::nOrmal

Definition at line 44 of file FaceElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPts

Definition at line 51 of file FaceElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::FaceElementForcesAndSourcesCore::num_nodes

Definition at line 41 of file FaceElementForcesAndSourcesCore.hpp.

◆ opContravariantTransform

OpSetContravariantPiolaTransformOnTriangle MoFEM::FaceElementForcesAndSourcesCore::opContravariantTransform

Definition at line 55 of file FaceElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnTriangle MoFEM::FaceElementForcesAndSourcesCore::opCovariantTransform

Definition at line 56 of file FaceElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnFace MoFEM::FaceElementForcesAndSourcesCore::opHOCoordsAndNormals

Definition at line 54 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOne

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOne

Definition at line 44 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOneAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentOneAtGaussPts

Definition at line 52 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwo

VectorDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwo

Definition at line 44 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwoAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCore::tangentTwoAtGaussPts

Definition at line 53 of file FaceElementForcesAndSourcesCore.hpp.


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