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

Edge finite elementUser is implementing own operator at Gauss points level, by own object derived from EdgeElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to rowOpPtrVector and rowColOpPtrVector. More...

#include <src/finite_elements/EdgeElementForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 default operator for EDGE element More...
 

Public Member Functions

 EdgeElementForcesAndSourcesCore (Interface &m_field)
 
MoFEMErrorCode calculateEdgeDirection ()
 
MoFEMErrorCode setIntegrationPts ()
 
MoFEMErrorCode calculateCoordsAtIntegrationPts ()
 
MoFEMErrorCode calculateBaseFunctionsOnElement ()
 
MoFEMErrorCode calculateHoCoordsAtIntegrationPts ()
 
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...
 
MoFEMErrorCode getDataOrderSpaceAndBase (const std::string &field_name, const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get maximal approximation order on entity More...
 
MoFEMErrorCode getEdgesSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getTrisSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getQuadSense (DataForcesAndSourcesCore &data) const
 
MoFEMErrorCode getEdgesDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getTrisDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getQuadDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getTetDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getPrismDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 
MoFEMErrorCode getEdgesDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTrisDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getQuadDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTetDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getPrismDataOrderSpaceAndBase (DataForcesAndSourcesCore &data, const std::string &field_name) 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 getTypeIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorInt &indices, VectorInt &local_indices) const
 get indices by type (generic function) More...
 
MoFEMErrorCode getTypeIndices (const boost::string_ref field_name, FENumeredDofEntity_multiIndex &dofs, EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get indices by type (generic function) 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 getEdgesRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Edges row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getEdgesColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Edges col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTrisRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tris row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTrisColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tris col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTetsRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tets row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getTetsColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Tets col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getQuadRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Quad row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getQuadColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Quad col indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getPrismRowIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Prism row indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getPrismColIndices (DataForcesAndSourcesCore &data, const std::string &field_name) const
 get Prism col indices from FENumeredDofEntity_multiIndex More...
 
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 getTypeFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, EntityType type, int side_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
 Get field data on entities. More...
 
MoFEMErrorCode getTypeFieldData (const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 
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 getEdgesFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTrisFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getQuadFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getTetsFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) const
 
MoFEMErrorCode getPrismFieldData (DataForcesAndSourcesCore &data, const std::string &field_name) 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...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop 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

DataForcesAndSourcesCore dataH1
 
DerivedDataForcesAndSourcesCore derivedDataH1
 
DataForcesAndSourcesCore dataHcurl
 
DerivedDataForcesAndSourcesCore derivedDataHcurl
 
DataForcesAndSourcesCore dataNoField
 
DataForcesAndSourcesCore dataNoFieldCol
 
std::string meshPositionsFieldName
 
MatrixDouble tangentAtGaussPts
 
OpGetHoTangentOnEdge opGetHoTangentOnEdge
 
OpSetCovariantPiolaTransformOnEdge opCovariantTransform
 
double lEngth
 
int numNodes
 
const EntityHandlecOnn
 
VectorDouble dIrection
 
VectorDouble cOords
 
MatrixDouble gaussPts
 
MatrixDouble coordsAtGaussPts
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
boost::ptr_vector< UserDataOperatoropPtrVector
 
- 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

Edge finite element

User is implementing own operator at Gauss points level, by own object derived from EdgeElementForcesAndSourcesCoreL::UserDataOperator. Arbitrary number of operator added pushing objects to rowOpPtrVector and rowColOpPtrVector.

Examples:
forces_and_sources_testing_edge_element.cpp.

Definition at line 40 of file EdgeElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ EdgeElementForcesAndSourcesCore()

MoFEM::EdgeElementForcesAndSourcesCore::EdgeElementForcesAndSourcesCore ( Interface m_field)

Definition at line 53 of file EdgeElementForcesAndSourcesCore.hpp.

54  : ForcesAndSourcesCore(m_field), dataH1(MBEDGE), derivedDataH1(dataH1),
56  dataNoFieldCol(MBEDGE), meshPositionsFieldName("MESH_NODE_POSITIONS"),
ForcesAndSourcesCore(Interface &m_field)
OpSetCovariantPiolaTransformOnEdge opCovariantTransform

Member Function Documentation

◆ calculateBaseFunctionsOnElement()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateBaseFunctionsOnElement ( )

Definition at line 167 of file EdgeElementForcesAndSourcesCore.cpp.

167  {
169  for (int b = AINSWORTH_LEGENDRE_BASE; b != LASTBASE; b++) {
170  if (dataH1.bAse.test(b)) {
171  switch (ApproximationBaseArray[b]) {
175  if (dataH1.spacesOnEntities[MBVERTEX].test(H1) &&
176  dataH1.basesOnEntities[MBVERTEX].test(b)) {
177  CHKERR EdgePolynomialBase().getValue(
178  gaussPts,
179  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
180  dataH1, H1, ApproximationBaseArray[b], NOBASE)));
181  }
182  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL) &&
183  dataH1.basesOnEntities[MBEDGE].test(b)) {
184  CHKERR EdgePolynomialBase().getValue(
185  gaussPts,
186  boost::shared_ptr<BaseFunctionCtx>(new EntPolynomialBaseCtx(
187  dataHcurl, HCURL, ApproximationBaseArray[b], NOBASE)));
188  }
189  if (dataH1.spacesOnEntities[MBEDGE].test(HDIV)) {
190  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
191  "Make no sense to have Hdiv space on edge (If you need h-div "
192  "space in 2d use h-curl space intead and use orthogonal "
193  "vector");
194  }
195  if (dataH1.spacesOnEntities[MBEDGE].test(L2)) {
196  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
197  "Not yet implemented");
198  }
199  break;
200  default:
201  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
202  "Not yet implemented");
203  }
204  }
205  }
207 }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:181
std::bitset< LASTBASE > basesOnEntities[MBMAXTYPE]
bases on entity types
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:495
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
field with continuous tangents
Definition: definitions.h:180
#define CHKERR
Inline error check.
Definition: definitions.h:614
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:439
continuous field
Definition: definitions.h:179
std::bitset< LASTBASE > bAse
bases on element
field with C-1 continuity
Definition: definitions.h:182

◆ calculateCoordsAtIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateCoordsAtIntegrationPts ( )

Definition at line 143 of file EdgeElementForcesAndSourcesCore.cpp.

143  {
145  const int nb_gauss_pts = gaussPts.size2();
146  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
147  if (nb_gauss_pts) {
149  &*gaussPts.data().begin());
151  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
152  &coordsAtGaussPts(0, 2));
153  FTensor::Tensor1<double, 3> t_coords_node0(cOords[0], cOords[1], cOords[2]);
154  FTensor::Tensor1<double, 3> t_coords_node1(cOords[3], cOords[4], cOords[5]);
156  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
157  t_coords(i) = N_MBEDGE0(t_ksi) * t_coords_node0(i) +
158  N_MBEDGE1(t_ksi) * t_coords_node1(i);
159  ++t_coords;
160  ++t_ksi;
161  }
162  }
164 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:519
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:80
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:79

◆ calculateEdgeDirection()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateEdgeDirection ( )

Definition at line 78 of file EdgeElementForcesAndSourcesCore.cpp.

78  {
81  CHKERR mField.get_moab().get_connectivity(ent, cOnn, numNodes, true);
82  cOords.resize(numNodes * 3, false);
83  CHKERR mField.get_moab().get_coords(cOnn, numNodes, &*cOords.data().begin());
84  dIrection.resize(3, false);
85  cblas_dcopy(3, &cOords[3], 1, &*dIrection.data().begin(), 1);
86  cblas_daxpy(3, -1., &cOords[0], 1, &*dIrection.data().begin(), 1);
87  lEngth = cblas_dnrm2(3, &*dIrection.data().begin(), 1);
89 }
virtual moab::Interface & get_moab()=0
void cblas_daxpy(const int N, const double alpha, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_daxpy.c:11
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:495
double cblas_dnrm2(const int N, const double *X, const int incX)
Definition: cblas_dnrm2.c:12
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
#define CHKERR
Inline error check.
Definition: definitions.h:614
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:439

◆ calculateHoCoordsAtIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateHoCoordsAtIntegrationPts ( )

Definition at line 210 of file EdgeElementForcesAndSourcesCore.cpp.

210  {
212  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
213  dataPtr->get<FieldName_mi_tag>().end()) {
218  } else {
219  tangentAtGaussPts.resize(0, 3, false);
220  }
222 }
MoFEMErrorCode getEdgesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
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:495
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
#define CHKERR
Inline error check.
Definition: definitions.h:614
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)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:439
MoFEMErrorCode getEdgesDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const

◆ operator()()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::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.

Definition at line 224 of file EdgeElementForcesAndSourcesCore.cpp.

224  {
226 
227  if (numeredEntFiniteElementPtr->getEntType() != MBEDGE)
229 
233  dataH1.dataOnEntities[MBEDGE][0].getSense() =
234  1; // set sense to 1, this is this entity
235 
236  // Hcurl
237  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
239  dataHcurl.dataOnEntities[MBEDGE][0].getSense() =
240  1; // set sense to 1, this is this entity
241  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
242  }
243 
244  /// Use the some node base
245  dataHcurl.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE) =
246  dataH1.dataOnEntities[MBVERTEX][0].getNSharedPtr(NOBASE);
247 
252 
253  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
254  // cerr << dataHcurl.dataOnEntities[MBEDGE][0].getN(AINSWORTH_LEGENDRE_BASE)
255  // << endl;
257  }
258 
261  std::vector<std::string> last_eval_field_name(2);
262  DataForcesAndSourcesCore *op_data[2];
263  FieldSpace space[2];
264  FieldApproximationBase base[2];
265 
266  boost::ptr_vector<UserDataOperator>::iterator oit, hi_oit;
267  oit = opPtrVector.begin();
268  hi_oit = opPtrVector.end();
269 
270  for (; oit != hi_oit; oit++) {
271 
272  oit->setPtrFE(this);
273 
274  for (int ss = 0; ss != 2; ss++) {
275 
276  std::string field_name = !ss ? oit->rowFieldName : oit->colFieldName;
277  const Field *field_struture = mField.get_field_structure(field_name);
278  BitFieldId data_id = field_struture->getId();
279 
280  if ((oit->getNumeredEntFiniteElementPtr()->getBitFieldIdData() & data_id)
281  .none()) {
282  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
283  "no data field < %s > on finite element < %s >",
284  field_name.c_str(), feName.c_str());
285  }
286 
287  if (oit->getOpType() & types[ss] ||
288  oit->getOpType() & UserDataOperator::OPROWCOL) {
289 
290  space[ss] = field_struture->getSpace();
291  switch (space[ss]) {
292  case NOSPACE:
293  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
294  case H1:
295  op_data[ss] = !ss ? &dataH1 : &derivedDataH1;
296  break;
297  case HCURL:
298  op_data[ss] = !ss ? &dataHcurl : &derivedDataHcurl;
299  break;
300  case HDIV:
301  SETERRQ(
302  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
303  "not make sanes on edge in 3d space (for 1d/2d not implemented)");
304  break;
305  case L2:
306  SETERRQ(
307  PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
308  "not make sanes on edge in 3d space (for 1d/2d not implemented)");
309  break;
310  case NOFIELD:
311  op_data[ss] = !ss ? &dataNoField : &dataNoFieldCol;
312  break;
313  case LASTSPACE:
314  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
315  break;
316  }
317 
318  base[ss] = field_struture->getApproxBase();
319  switch (base[ss]) {
323  break;
324  default:
325  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
326  "unknown or not implemented base");
327  break;
328  }
329 
330  if (last_eval_field_name[ss] != field_name) {
331  switch (space[ss]) {
332  case NOSPACE:
333  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
334  case H1:
335  if (!ss) {
336  CHKERR getRowNodesIndices(*op_data[ss], field_name);
337  } else {
338  CHKERR getColNodesIndices(*op_data[ss], field_name);
339  }
340  CHKERR getNodesFieldData(*op_data[ss], field_name);
341  case HCURL:
342  if (!ss) {
343  CHKERR getEdgesRowIndices(*op_data[ss], field_name);
344  } else {
345  CHKERR getEdgesColIndices(*op_data[ss], field_name);
346  }
347  CHKERR getEdgesDataOrderSpaceAndBase(*op_data[ss], field_name);
348  CHKERR getEdgesFieldData(*op_data[ss], field_name);
349  break;
350  case HDIV:
351  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
352  "not make sanes on edge");
353  break;
354  case L2:
355  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
356  "not make sanes on edge");
357  break;
358  case NOFIELD:
359  if (!getNinTheLoop()) {
360  // NOFIELD data are the same for each element, can be retrieved
361  // only once
362  if (!ss) {
363  CHKERR getNoFieldRowIndices(*op_data[ss], field_name);
364  } else {
365  CHKERR getNoFieldColIndices(*op_data[ss], field_name);
366  }
367  CHKERR getNoFieldFieldData(*op_data[ss], field_name);
368  }
369  break;
370  case LASTSPACE:
371  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "unknown space");
372  break;
373  }
374  last_eval_field_name[ss] = field_name;
375  }
376  }
377  }
378 
379  if (oit->getOpType() & UserDataOperator::OPROW) {
380  try {
381  ierr = oit->opRhs(*op_data[0], oit->doVertices, oit->doEdges, false,
382  false, false, false);
383  CHKERRG(ierr);
384  } catch (std::exception &ex) {
385  std::ostringstream ss;
386  ss << "Operator "
387  << boost::typeindex::type_id_runtime(*oit).pretty_name()
388  << " operator number "
389  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
390  opPtrVector.begin(), oit)
391  << " throw in method: " << ex.what() << " at line " << __LINE__
392  << " in file " << __FILE__;
393  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
394  }
395  }
396 
397  if (oit->getOpType() & UserDataOperator::OPCOL) {
398  try {
399  ierr = oit->opRhs(*op_data[1], oit->doVertices, oit->doEdges, false,
400  false, false, false);
401  CHKERRG(ierr);
402  } catch (std::exception &ex) {
403  std::ostringstream ss;
404  ss << "Operator "
405  << boost::typeindex::type_id_runtime(*oit).pretty_name()
406  << " operator number "
407  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
408  opPtrVector.begin(), oit)
409  << " throw in method: " << ex.what() << " at line " << __LINE__
410  << " in file " << __FILE__;
411  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
412  }
413  }
414 
415  if (oit->getOpType() & UserDataOperator::OPROWCOL) {
416  try {
417  ierr = oit->opLhs(*op_data[0], *op_data[1], oit->sYmm);
418  CHKERRG(ierr);
419  } catch (std::exception &ex) {
420  std::ostringstream ss;
421  ss << "Operator "
422  << boost::typeindex::type_id_runtime(*oit).pretty_name()
423  << " operator number "
424  << std::distance<boost::ptr_vector<UserDataOperator>::iterator>(
425  opPtrVector.begin(), oit)
426  << " throw in method: " << ex.what() << " at line " << __LINE__
427  << " in file " << __FILE__;
428  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
429  }
430  }
431  }
432 
434 }
OpSetCovariantPiolaTransformOnEdge opCovariantTransform
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
boost::ptr_vector< UserDataOperator > opPtrVector
field with continuous normal traction
Definition: definitions.h:181
MoFEMErrorCode getNoFieldColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
MoFEMErrorCode getEdgesRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Edges row indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getColNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesDataOrder(DataForcesAndSourcesCore &data, const FieldSpace space) const
scalar or vector of scalars describe (no true field)
Definition: definitions.h:178
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:495
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:562
std::string feName
std::bitset< BITFIELDID_SIZE > BitFieldId
Definition: Common.hpp:149
int getNinTheLoop() const
get number of evaluated element in the loop
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
MoFEMErrorCode getNoFieldRowIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get col NoField indices
OpType
Controls loop over entities on element.
MoFEMErrorCode getNoFieldFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs) const
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
FieldSpace in [ 0, LASTSPACE )
Definition: definitions.h:183
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
FieldApproximationBase
approximation base
Definition: definitions.h:140
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142
field with continuous tangents
Definition: definitions.h:180
FieldSpace
approximation spaces
Definition: definitions.h:176
#define CHKERR
Inline error check.
Definition: definitions.h:614
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 getEdgesColIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get Edges col indices from FENumeredDofEntity_multiIndex
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:439
continuous field
Definition: definitions.h:179
virtual const Field * get_field_structure(const std::string &name)=0
get field structure
MoFEMErrorCode getRowNodesIndices(DataForcesAndSourcesCore &data, const std::string &field_name) const
get row node indices from FENumeredDofEntity_multiIndex
MoFEMErrorCode getEdgesDataOrderSpaceAndBase(DataForcesAndSourcesCore &data, const std::string &field_name) const
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
field with C-1 continuity
Definition: definitions.h:182

◆ setIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::setIntegrationPts ( )

Definition at line 91 of file EdgeElementForcesAndSourcesCore.cpp.

91  {
93  int order_data = getMaxDataOrder();
94  int order_row = getMaxRowOrder();
95  int order_col = getMaxColOrder();
96  int rule = getRule(order_row, order_col, order_data);
97  int nb_gauss_pts;
98  if (rule >= 0) {
99  if (rule < QUAD_1D_TABLE_SIZE) {
100  if (QUAD_1D_TABLE[rule]->dim != 1) {
101  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
102  }
103  if (QUAD_1D_TABLE[rule]->order < rule) {
104  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
105  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order, rule);
106  }
107  nb_gauss_pts = QUAD_1D_TABLE[rule]->npoints;
108  gaussPts.resize(2, nb_gauss_pts, false);
109  cblas_dcopy(nb_gauss_pts, &QUAD_1D_TABLE[rule]->points[1], 2,
110  &gaussPts(0, 0), 1);
111  cblas_dcopy(nb_gauss_pts, QUAD_1D_TABLE[rule]->weights, 1,
112  &gaussPts(1, 0), 1);
113  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 2,
114  false);
115  double *shape_ptr =
116  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
117  cblas_dcopy(2 * nb_gauss_pts, QUAD_1D_TABLE[rule]->points, 1, shape_ptr,
118  1);
119  } else {
120  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
121  "rule > quadrature order %d < %d", rule, QUAD_1D_TABLE_SIZE);
122  nb_gauss_pts = 0;
123  }
124  } else {
125  // If rule is negative, set user defined integration points
126  CHKERR setGaussPts(order_row, order_col, order_data);
127  nb_gauss_pts = gaussPts.size2();
128  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 2,
129  false);
130  if (nb_gauss_pts) {
131  for (int gg = 0; gg != nb_gauss_pts; gg++) {
132  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0) =
133  N_MBEDGE0(gaussPts(0, gg));
134  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 1) =
135  N_MBEDGE1(gaussPts(0, gg));
136  }
137  }
138  }
140 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
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:495
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:80
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
virtual MoFEMErrorCode setGaussPts(int order)
It will be removed in the future use other variant.
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:79
#define CHKERR
Inline error check.
Definition: definitions.h:614
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:439
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.
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163

Member Data Documentation

◆ cOnn

const EntityHandle* MoFEM::EdgeElementForcesAndSourcesCore::cOnn

Definition at line 63 of file EdgeElementForcesAndSourcesCore.hpp.

◆ cOords

VectorDouble MoFEM::EdgeElementForcesAndSourcesCore::cOords

Definition at line 65 of file EdgeElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::coordsAtGaussPts

Definition at line 67 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dataH1

DataForcesAndSourcesCore MoFEM::EdgeElementForcesAndSourcesCore::dataH1

Definition at line 42 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dataHcurl

DataForcesAndSourcesCore MoFEM::EdgeElementForcesAndSourcesCore::dataHcurl

Definition at line 44 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dataNoField

DataForcesAndSourcesCore MoFEM::EdgeElementForcesAndSourcesCore::dataNoField

Definition at line 46 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dataNoFieldCol

DataForcesAndSourcesCore MoFEM::EdgeElementForcesAndSourcesCore::dataNoFieldCol

Definition at line 46 of file EdgeElementForcesAndSourcesCore.hpp.

◆ derivedDataH1

DerivedDataForcesAndSourcesCore MoFEM::EdgeElementForcesAndSourcesCore::derivedDataH1

Definition at line 43 of file EdgeElementForcesAndSourcesCore.hpp.

◆ derivedDataHcurl

DerivedDataForcesAndSourcesCore MoFEM::EdgeElementForcesAndSourcesCore::derivedDataHcurl

Definition at line 45 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dIrection

VectorDouble MoFEM::EdgeElementForcesAndSourcesCore::dIrection

Definition at line 64 of file EdgeElementForcesAndSourcesCore.hpp.

◆ gaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::gaussPts

Definition at line 66 of file EdgeElementForcesAndSourcesCore.hpp.

◆ lEngth

double MoFEM::EdgeElementForcesAndSourcesCore::lEngth

Definition at line 60 of file EdgeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::EdgeElementForcesAndSourcesCore::meshPositionsFieldName

Definition at line 47 of file EdgeElementForcesAndSourcesCore.hpp.

◆ numNodes

int MoFEM::EdgeElementForcesAndSourcesCore::numNodes

Definition at line 62 of file EdgeElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnEdge MoFEM::EdgeElementForcesAndSourcesCore::opCovariantTransform

Definition at line 51 of file EdgeElementForcesAndSourcesCore.hpp.

◆ opGetHoTangentOnEdge

OpGetHoTangentOnEdge MoFEM::EdgeElementForcesAndSourcesCore::opGetHoTangentOnEdge

Definition at line 50 of file EdgeElementForcesAndSourcesCore.hpp.

◆ tangentAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::tangentAtGaussPts

Definition at line 49 of file EdgeElementForcesAndSourcesCore.hpp.


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