v0.10.0
Classes | Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
MoFEM::EdgeElementForcesAndSourcesCoreBase Struct Reference

Edge finite element. More...

#include <src/finite_elements/EdgeElementForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 default operator for EDGE element More...
 

Public Types

enum  Switches { NO_HO_GEOMETRY = 1 << 0, NO_COVARIANT_TRANSFORM_HCURL = 1 << 2 }
 
- 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::PetscData
enum  DataContext {
  CTX_SET_NONE = 0, CTX_SET_F = 1 << 0, CTX_SET_A = 1 << 1, CTX_SET_B = 1 << 2,
  CTX_SET_X = 1 << 3, CTX_SET_X_T = 1 << 4, CTX_SET_X_TT = 1 << 6, CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- 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
}
 

Public Member Functions

template<int SWITCH>
MoFEMErrorCode OpSwitch ()
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
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 operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
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...
 
const DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side) const
 Get the entity data. More...
 
DataForcesAndSourcesCore::EntDatagetEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (const MOFEMuuid &uuid, UnknownInterface **iface) const
 
 FEMethod ()
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. More...
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
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
 
- Public Member Functions inherited from MoFEM::BasicMethod
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
 PetscData ()
 
virtual ~PetscData ()=default
 
- 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 ()=default
 
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 ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
DEPRECATED MoFEMErrorCode setSnesCtx (SNESContext ctx)
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 
DEPRECATED MoFEMErrorCode setTsCtx (TSContext ctx)
 

Public Attributes

std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
RuleHookFun setRuleHook
 Set function to calculate integration rule. More...
 
MatrixDouble gaussPts
 Matrix of integration points. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 

Protected Member Functions

 EdgeElementForcesAndSourcesCoreBase (Interface &m_field)
 
MoFEMErrorCode calculateEdgeDirection ()
 
MoFEMErrorCode setIntegrationPts ()
 
MoFEMErrorCode calculateCoordsAtIntegrationPts ()
 
MoFEMErrorCode calculateHoCoordsAtIntegrationPts ()
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< DataForcesAndSourcesCore::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (DataForcesAndSourcesCore &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (DataForcesAndSourcesCore &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceTriNodes (DataForcesAndSourcesCore &data) const
 Get nodes on triangles. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (DataForcesAndSourcesCore &data) const
 Get field approximation space and base on entities. More...
 
virtual int getRule (int order_row, int order_col, int order_data)
 another variant of getRule More...
 
virtual MoFEMErrorCode setGaussPts (int order_row, int order_col, int order_data)
 set user specific integration rule More...
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement (const FieldApproximationBase b)
 Calculate base functions. More...
 
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement ()
 Calculate base functions. More...
 
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement ()
 Calculate Bernstein-Bezier base. More...
 
MoFEMErrorCode createDataOnElement ()
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getNodesIndices (const std::string &field_name, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) 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...
 
template<typename EXTRACTOR >
MoFEMErrorCode getEntityIndices (DataForcesAndSourcesCore &data, const std::string &field_name, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) 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, boost::shared_ptr< 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 getNoFieldFieldData (const std::string field_name, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (DataForcesAndSourcesCore &data, const std::string 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 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)
 
virtual MoFEMErrorCode setGaussPts (int order)
 
- 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...
 

Protected Attributes

MatrixDouble tangentAtGaussPts
 
OpGetHoTangentOnEdge opGetHoTangentOnEdge
 
OpSetCovariantPiolaTransformOnEdge opCovariantTransform
 
double lEngth
 
int numNodes
 
const EntityHandlecOnn
 
VectorDouble dIrection
 
VectorDouble cOords
 
MatrixDouble coordsAtGaussPts
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
DataForcesAndSourcesCoredataNoField
 
DataForcesAndSourcesCoredataH1
 
DataForcesAndSourcesCoredataHcurl
 
DataForcesAndSourcesCoredataHdiv
 
DataForcesAndSourcesCoredataL2
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 

Friends

class FaceElementForcesAndSourcesCoreOnSideBase
 

Additional Inherited Members

- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

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
helmholtz.cpp.

Definition at line 44 of file EdgeElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ Switches

Enumerator
NO_HO_GEOMETRY 
NO_COVARIANT_TRANSFORM_HCURL 

Definition at line 146 of file EdgeElementForcesAndSourcesCore.hpp.

146  {
147  NO_HO_GEOMETRY = 1 << 0,
149  };

Constructor & Destructor Documentation

◆ EdgeElementForcesAndSourcesCoreBase()

MoFEM::EdgeElementForcesAndSourcesCoreBase::EdgeElementForcesAndSourcesCoreBase ( Interface m_field)
protected

Definition at line 34 of file EdgeElementForcesAndSourcesCore.cpp.

36  : ForcesAndSourcesCore(m_field),
37  meshPositionsFieldName("MESH_NODE_POSITIONS"),
41  boost::shared_ptr<BaseFunction>(new EdgePolynomialBase());
42 }

Member Function Documentation

◆ calculateCoordsAtIntegrationPts()

MoFEMErrorCode MoFEM::EdgeElementForcesAndSourcesCoreBase::calculateCoordsAtIntegrationPts ( )
protected

Definition at line 110 of file EdgeElementForcesAndSourcesCore.cpp.

110  {
112  const int nb_gauss_pts = gaussPts.size2();
113  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
114  if (nb_gauss_pts) {
116  &*gaussPts.data().begin());
118  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
119  &coordsAtGaussPts(0, 2));
120  FTensor::Tensor1<double, 3> t_coords_node0(cOords[0], cOords[1], cOords[2]);
121  FTensor::Tensor1<double, 3> t_coords_node1(cOords[3], cOords[4], cOords[5]);
123  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
124  t_coords(i) = N_MBEDGE0(t_ksi) * t_coords_node0(i) +
125  N_MBEDGE1(t_ksi) * t_coords_node1(i);
126  ++t_coords;
127  ++t_ksi;
128  }
129  }
131 }

◆ calculateEdgeDirection()

MoFEMErrorCode MoFEM::EdgeElementForcesAndSourcesCoreBase::calculateEdgeDirection ( )
protected

Definition at line 44 of file EdgeElementForcesAndSourcesCore.cpp.

44  {
47  CHKERR mField.get_moab().get_connectivity(ent, cOnn, numNodes, true);
48  cOords.resize(numNodes * 3, false);
49  CHKERR mField.get_moab().get_coords(cOnn, numNodes, &*cOords.data().begin());
50  dIrection.resize(3, false);
51  cblas_dcopy(3, &cOords[3], 1, &*dIrection.data().begin(), 1);
52  cblas_daxpy(3, -1., &cOords[0], 1, &*dIrection.data().begin(), 1);
53  lEngth = cblas_dnrm2(3, &*dIrection.data().begin(), 1);
55 }

◆ calculateHoCoordsAtIntegrationPts()

MoFEMErrorCode MoFEM::EdgeElementForcesAndSourcesCoreBase::calculateHoCoordsAtIntegrationPts ( )
protected

Definition at line 134 of file EdgeElementForcesAndSourcesCore.cpp.

134  {
136 
137  auto check_field = [&]() {
138  auto field_it =
139  fieldsPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName);
140  if (field_it != fieldsPtr->get<FieldName_mi_tag>().end()) {
141  if (
142 
143  (numeredEntFiniteElementPtr->getBitFieldIdData() &
144  (*field_it)->getId())
145  .any()
146 
147  )
148  return true;
149  }
150  return false;
151  };
152 
153  if (check_field()) {
157  } else {
158  tangentAtGaussPts.resize(0, 3, false);
159  }
161 }

◆ OpSwitch()

template<int SWITCH>
MoFEMErrorCode MoFEM::EdgeElementForcesAndSourcesCoreBase::OpSwitch

Definition at line 203 of file EdgeElementForcesAndSourcesCore.hpp.

203  {
205 
206  if (numeredEntFiniteElementPtr->getEntType() != MBEDGE)
208 
210 
211  DataForcesAndSourcesCore &data_curl = *dataOnElement[HCURL];
212 
215  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
216  dataH1.dataOnEntities[MBEDGE][0].getSense() =
217  1; // set sense to 1, this is this entity
218 
219  // Hcurl
220  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
221  CHKERR getEntityDataOrder<MBEDGE>(data_curl, HCURL);
222  data_curl.dataOnEntities[MBEDGE][0].getSense() =
223  1; // set sense to 1, this is this entity
224  data_curl.spacesOnEntities[MBEDGE].set(HCURL);
225  }
226 
231  if (!(SWITCH & NO_HO_GEOMETRY))
233 
234  if (!(SWITCH & NO_COVARIANT_TRANSFORM_HCURL))
235  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL))
236  CHKERR opCovariantTransform.opRhs(data_curl);
237 
238  // Iterate over operators
240 
242 }

◆ setIntegrationPts()

MoFEMErrorCode MoFEM::EdgeElementForcesAndSourcesCoreBase::setIntegrationPts ( )
protected

Definition at line 57 of file EdgeElementForcesAndSourcesCore.cpp.

57  {
59  int order_data = getMaxDataOrder();
60  int order_row = getMaxRowOrder();
61  int order_col = getMaxColOrder();
62  int rule = getRule(order_row, order_col, order_data);
63 
64  int nb_gauss_pts;
65  if (rule >= 0) {
66  if (rule < QUAD_1D_TABLE_SIZE) {
67  if (QUAD_1D_TABLE[rule]->dim != 1) {
68  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
69  }
70  if (QUAD_1D_TABLE[rule]->order < rule) {
71  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
72  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order, rule);
73  }
74  nb_gauss_pts = QUAD_1D_TABLE[rule]->npoints;
75  gaussPts.resize(2, nb_gauss_pts, false);
76  cblas_dcopy(nb_gauss_pts, &QUAD_1D_TABLE[rule]->points[1], 2,
77  &gaussPts(0, 0), 1);
78  cblas_dcopy(nb_gauss_pts, QUAD_1D_TABLE[rule]->weights, 1,
79  &gaussPts(1, 0), 1);
80  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 2,
81  false);
82  double *shape_ptr =
83  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
84  cblas_dcopy(2 * nb_gauss_pts, QUAD_1D_TABLE[rule]->points, 1, shape_ptr,
85  1);
86  } else {
87  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
88  "rule > quadrature order %d < %d", rule, QUAD_1D_TABLE_SIZE);
89  nb_gauss_pts = 0;
90  }
91  } else {
92  // If rule is negative, set user defined integration points
93  CHKERR setGaussPts(order_row, order_col, order_data);
94  nb_gauss_pts = gaussPts.size2();
95  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 2,
96  false);
97  if (nb_gauss_pts) {
98  for (int gg = 0; gg != nb_gauss_pts; gg++) {
99  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0) =
100  N_MBEDGE0(gaussPts(0, gg));
101  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 1) =
102  N_MBEDGE1(gaussPts(0, gg));
103  }
104  }
105  }
107 }

Friends And Related Function Documentation

◆ FaceElementForcesAndSourcesCoreOnSideBase

Definition at line 173 of file EdgeElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ cOnn

const EntityHandle* MoFEM::EdgeElementForcesAndSourcesCoreBase::cOnn
protected

Definition at line 163 of file EdgeElementForcesAndSourcesCore.hpp.

◆ cOords

VectorDouble MoFEM::EdgeElementForcesAndSourcesCoreBase::cOords
protected

Definition at line 165 of file EdgeElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCoreBase::coordsAtGaussPts
protected

Definition at line 166 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dIrection

VectorDouble MoFEM::EdgeElementForcesAndSourcesCoreBase::dIrection
protected

Definition at line 164 of file EdgeElementForcesAndSourcesCore.hpp.

◆ lEngth

double MoFEM::EdgeElementForcesAndSourcesCoreBase::lEngth
protected

Definition at line 160 of file EdgeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::EdgeElementForcesAndSourcesCoreBase::meshPositionsFieldName

Definition at line 46 of file EdgeElementForcesAndSourcesCore.hpp.

◆ numNodes

int MoFEM::EdgeElementForcesAndSourcesCoreBase::numNodes
protected

Definition at line 162 of file EdgeElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnEdge MoFEM::EdgeElementForcesAndSourcesCoreBase::opCovariantTransform
protected

Definition at line 158 of file EdgeElementForcesAndSourcesCore.hpp.

◆ opGetHoTangentOnEdge

OpGetHoTangentOnEdge MoFEM::EdgeElementForcesAndSourcesCoreBase::opGetHoTangentOnEdge
protected

Definition at line 157 of file EdgeElementForcesAndSourcesCore.hpp.

◆ tangentAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCoreBase::tangentAtGaussPts
protected

Definition at line 156 of file EdgeElementForcesAndSourcesCore.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:516
cblas_daxpy
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
EntityHandle
H1
@ H1
continuous field
Definition: definitions.h:177
MoFEM::EdgeElementForcesAndSourcesCoreBase::cOords
VectorDouble cOords
Definition: EdgeElementForcesAndSourcesCore.hpp:165
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:451
MoFEM::EdgeElementForcesAndSourcesCoreBase::calculateHoCoordsAtIntegrationPts
MoFEMErrorCode calculateHoCoordsAtIntegrationPts()
Definition: EdgeElementForcesAndSourcesCore.cpp:134
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:36
MoFEM::EdgeElementForcesAndSourcesCoreBase::numNodes
int numNodes
Definition: EdgeElementForcesAndSourcesCore.hpp:162
NOBASE
@ NOBASE
Definition: definitions.h:151
MoFEM::ForcesAndSourcesCore::setGaussPts
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: ForcesAndSourcesCore.cpp:1597
MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
Definition: ForcesAndSourcesCore.cpp:1087
MoFEM::BasicMethod::fieldsPtr
const Field_multiIndex * fieldsPtr
raw pointer to fields container
Definition: LoopMethods.hpp:269
FTensor::Tensor1
Definition: Tensor1_value.hpp:9
N_MBEDGE0
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:81
dim
const int dim
Definition: elec_phys_2D.cpp:12
QUAD_::npoints
int npoints
Definition: quad.h:29
MoFEM::EdgeElementForcesAndSourcesCoreBase::dIrection
VectorDouble dIrection
Definition: EdgeElementForcesAndSourcesCore.hpp:164
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::ForcesAndSourcesCore::getMaxRowOrder
int getMaxRowOrder() const
Get max order of approximation for field in rows.
Definition: ForcesAndSourcesCore.cpp:138
MoFEM::EdgeElementForcesAndSourcesCoreBase::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
Definition: EdgeElementForcesAndSourcesCore.hpp:166
QUAD_1D_TABLE
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:604
N_MBEDGE1
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:82
MoFEM::ForcesAndSourcesCore::dataH1
DataForcesAndSourcesCore & dataH1
Definition: ForcesAndSourcesCore.hpp:807
MoFEM::ForcesAndSourcesCore::getRule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: ForcesAndSourcesCore.cpp:1591
MoFEM::EdgeElementForcesAndSourcesCoreBase::opCovariantTransform
OpSetCovariantPiolaTransformOnEdge opCovariantTransform
Definition: EdgeElementForcesAndSourcesCore.hpp:158
MoFEM::ForcesAndSourcesCore::loopOverOperators
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
Definition: ForcesAndSourcesCore.cpp:1304
MoFEM::ForcesAndSourcesCore::getEntityFieldData
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
Definition: ForcesAndSourcesCore.cpp:722
FTensor::Index< 'i', 3 >
MoFEM::EdgeElementForcesAndSourcesCoreBase::tangentAtGaussPts
MatrixDouble tangentAtGaussPts
Definition: EdgeElementForcesAndSourcesCore.hpp:156
StdRDOperators::order
const int order
Definition: rd_stdOperators.hpp:28
MoFEM::ForcesAndSourcesCore::getNodesFieldData
MoFEMErrorCode getNodesFieldData(DataForcesAndSourcesCore &data, const std::string &field_name) const
Get data on nodes.
Definition: ForcesAndSourcesCore.cpp:602
MoFEM::EdgeElementForcesAndSourcesCoreBase::meshPositionsFieldName
std::string meshPositionsFieldName
Definition: EdgeElementForcesAndSourcesCore.hpp:46
MoFEM::ForcesAndSourcesCore::getElementPolynomialBase
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
Definition: ForcesAndSourcesCore.hpp:432
FTensor::Tensor0
Definition: Tensor0.hpp:17
MoFEM::ForcesAndSourcesCore::calHierarchicalBaseFunctionsOnElement
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
Definition: ForcesAndSourcesCore.cpp:1067
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:389
MoFEM::EdgeElementForcesAndSourcesCoreBase::lEngth
double lEngth
Definition: EdgeElementForcesAndSourcesCore.hpp:160
MoFEM::DataOperator::opRhs
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool error_if_no_base=false)
Definition: DataOperators.cpp:171
MoFEM::ForcesAndSourcesCore::dataOnElement
const std::array< boost::shared_ptr< DataForcesAndSourcesCore >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
Definition: ForcesAndSourcesCore.hpp:797
QUAD_1D_TABLE_SIZE
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163
MoFEM::DataForcesAndSourcesCore::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: DataStructures.hpp:1040
MoFEM::EdgeElementForcesAndSourcesCoreBase::NO_COVARIANT_TRANSFORM_HCURL
@ NO_COVARIANT_TRANSFORM_HCURL
Definition: EdgeElementForcesAndSourcesCore.hpp:148
MoFEM::ForcesAndSourcesCore::getMaxDataOrder
int getMaxDataOrder() const
Get max order of approximation for data fields.
Definition: ForcesAndSourcesCore.cpp:127
MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
Definition: ForcesAndSourcesCore.cpp:971
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:178
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:123
MoFEM::EdgeElementForcesAndSourcesCoreBase::calculateEdgeDirection
MoFEMErrorCode calculateEdgeDirection()
Definition: EdgeElementForcesAndSourcesCore.cpp:44
MoFEM::ForcesAndSourcesCore::getMaxColOrder
int getMaxColOrder() const
Get max order of approximation for field in columns.
Definition: ForcesAndSourcesCore.cpp:142
MoFEM::DataForcesAndSourcesCore::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: DataStructures.hpp:1034
MoFEM::EdgeElementForcesAndSourcesCoreBase::NO_HO_GEOMETRY
@ NO_HO_GEOMETRY
Definition: EdgeElementForcesAndSourcesCore.hpp:147
MoFEM::EdgeElementForcesAndSourcesCoreBase::cOnn
const EntityHandle * cOnn
Definition: EdgeElementForcesAndSourcesCore.hpp:163
cblas_dnrm2
double cblas_dnrm2(const int N, const double *X, const int incX)
Definition: cblas_dnrm2.c:12
MoFEM::EdgeElementForcesAndSourcesCoreBase::calculateCoordsAtIntegrationPts
MoFEMErrorCode calculateCoordsAtIntegrationPts()
Definition: EdgeElementForcesAndSourcesCore.cpp:110
MoFEM::ForcesAndSourcesCore::createDataOnElement
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
Definition: ForcesAndSourcesCore.cpp:1256
MoFEM::EdgeElementForcesAndSourcesCoreBase::opGetHoTangentOnEdge
OpGetHoTangentOnEdge opGetHoTangentOnEdge
Definition: EdgeElementForcesAndSourcesCore.hpp:157
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:509
MoFEM::EdgeElementForcesAndSourcesCoreBase::setIntegrationPts
MoFEMErrorCode setIntegrationPts()
Definition: EdgeElementForcesAndSourcesCore.cpp:57
i
FTensor::Index< 'i', 3 > i
Definition: matrix_function.cpp:18
MoFEM::ForcesAndSourcesCore::ForcesAndSourcesCore
ForcesAndSourcesCore(Interface &m_field)
Definition: ForcesAndSourcesCore.cpp:54
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:485
cblas_dcopy
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:415