v0.14.0
Classes | Public Types | Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
MoFEM::EdgeElementForcesAndSourcesCore Struct Reference

Edge finite element. 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 Types

enum  Switches
 
- Public Types inherited from MoFEM::ForcesAndSourcesCore
typedef boost::function< int(int order_row, int order_col, int order_data)> RuleHookFun
 
typedef boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> GaussHookFun
 
- 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

 EdgeElementForcesAndSourcesCore (Interface &m_field)
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore
 ForcesAndSourcesCore (Interface &m_field)
 
boost::ptr_deque< 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...
 
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...
 
auto & getEntData (const FieldSpace space, const EntityType type, const int side)
 Get the entity data. More...
 
auto & getDataOnElementBySpaceArray ()
 Get data on entities and space. More...
 
auto & getDerivedDataOnElementBySpaceArray ()
 Get derived data on entities and space. More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
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
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- 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...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. 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...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. 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
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface reference 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
 
- Public Member Functions inherited from MoFEM::SnesMethod
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Public Attributes

std::string meshPositionsFieldName
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore
InterfacemField
 
RuleHookFun getRuleHook
 Hook to get rule. More...
 
GaussHookFun 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
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. 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
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- 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 number 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...
 
PetscReal ts_dt
 time step size 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...
 

Static Public Attributes

static FTensor::Tensor1< double, 3 > tFaceOrientation
 
- 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)
 

Protected Member Functions

MoFEMErrorCode calculateEdgeDirection ()
 
MoFEMErrorCode setIntegrationPts ()
 
MoFEMErrorCode calculateCoordsAtIntegrationPts ()
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore
MoFEMErrorCode getEntitySense (const EntityType type, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 get sense (orientation) of entity More...
 
MoFEMErrorCode getEntityDataOrder (const EntityType type, const FieldSpace space, boost::ptr_vector< EntitiesFieldData::EntData > &data) const
 Get the entity data order. More...
 
template<EntityType type>
MoFEMErrorCode getEntitySense (EntitiesFieldData &data) const
 Get the entity sense (orientation) More...
 
template<EntityType type>
MoFEMErrorCode getEntityDataOrder (EntitiesFieldData &data, const FieldSpace space) const
 Get the entity data order for given space. More...
 
MoFEMErrorCode getFaceNodes (EntitiesFieldData &data) const
 Get nodes on faces. More...
 
MoFEMErrorCode getSpacesAndBaseOnEntities (EntitiesFieldData &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 (EntityType type)
 Create a entity data on element object. More...
 
MoFEMErrorCode loopOverOperators ()
 Iterate user data operators. More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getNodesIndices (const int bit_number, FieldEntity_vector_view &ents_field, VectorInt &nodes_indices, VectorInt &local_nodes_indices, EXTRACTOR &&extractor) const
 get node indices More...
 
MoFEMErrorCode getRowNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get row node indices from FENumeredDofEntity_multiIndex More...
 
MoFEMErrorCode getColNodesIndices (EntitiesFieldData &data, const int bit_number) const
 get col node indices from FENumeredDofEntity_multiIndex More...
 
template<typename EXTRACTOR >
MoFEMErrorCode getEntityIndices (EntitiesFieldData &data, const int bit_number, FieldEntity_vector_view &ents_field, const EntityType type_lo, const EntityType type_hi, EXTRACTOR &&extractor) const
 
MoFEMErrorCode getEntityRowIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getEntityColIndices (EntitiesFieldData &data, const int bit_number, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
 
MoFEMErrorCode getNoFieldIndices (const int bit_number, boost::shared_ptr< FENumeredDofEntity_multiIndex > dofs, VectorInt &nodes_indices) const
 get NoField indices More...
 
MoFEMErrorCode getNoFieldRowIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getNoFieldColIndices (EntitiesFieldData &data, const int bit_number) const
 get col NoField indices More...
 
MoFEMErrorCode getBitRefLevelOnData ()
 
MoFEMErrorCode getNoFieldFieldData (const int bit_number, VectorDouble &ent_field_data, VectorDofs &ent_field_dofs, VectorFieldEntities &ent_field) const
 Get field data on nodes. More...
 
MoFEMErrorCode getNoFieldFieldData (EntitiesFieldData &data, const int bit_number) const
 
MoFEMErrorCode getNodesFieldData (EntitiesFieldData &data, const int bit_number) const
 Get data on nodes. More...
 
MoFEMErrorCode getEntityFieldData (EntitiesFieldData &data, const int bit_number, 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 Attributes

MatrixDouble tangentAtGaussPts
 
MatrixDouble normalsAtGaussPts
 
doublelEngth
 
int numNodes
 
const EntityHandlecOnn
 
VectorDouble dIrection
 
VectorDouble cOords
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEdataOnElement
 Entity data on element entity rows fields. More...
 
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACEderivedDataOnElement
 Entity data on element entity columns fields. More...
 
EntitiesFieldDatadataNoField
 
EntitiesFieldDatadataH1
 
EntitiesFieldDatadataHcurl
 
EntitiesFieldDatadataHdiv
 
EntitiesFieldDatadataL2
 
boost::ptr_deque< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 
EntityType lastEvaluatedElementEntityType
 Last evaluated type of element entity. More...
 
MatrixDouble coordsAtGaussPts
 coordinated at gauss points More...
 
double elementMeasure
 

Friends

class FaceElementForcesAndSourcesCoreOnSide
 

Additional Inherited Members

- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. 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
continuity_check_on_skeleton_with_simple_2d_for_h1.cpp, continuity_check_on_skeleton_with_simple_2d_for_hcurl.cpp, continuity_check_on_skeleton_with_simple_2d_for_hdiv.cpp, elasticity.cpp, forces_and_sources_testing_edge_element.cpp, hcurl_divergence_operator_2d.cpp, helmholtz.cpp, phase.cpp, and prism_elements_from_surface.cpp.

Definition at line 30 of file EdgeElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ Switches

Definition at line 40 of file EdgeElementForcesAndSourcesCore.hpp.

40 {};

Constructor & Destructor Documentation

◆ EdgeElementForcesAndSourcesCore()

EdgeElementForcesAndSourcesCore::EdgeElementForcesAndSourcesCore ( Interface m_field)

Definition at line 25 of file EdgeElementForcesAndSourcesCore.cpp.

27  : ForcesAndSourcesCore(m_field),
28  meshPositionsFieldName("MESH_NODE_POSITIONS"), lEngth(elementMeasure) {
30  boost::shared_ptr<BaseFunction>(new EdgePolynomialBase());
32  "Problem with creation data on element");
33 }

Member Function Documentation

◆ calculateCoordsAtIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateCoordsAtIntegrationPts ( )
protected

Definition at line 101 of file EdgeElementForcesAndSourcesCore.cpp.

101  {
104  const auto nb_gauss_pts = gaussPts.size2();
105  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
106  tangentAtGaussPts.resize(nb_gauss_pts, 3, false);
107 
108  if (nb_gauss_pts) {
110  &*gaussPts.data().begin());
112  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
113  &coordsAtGaussPts(0, 2));
115  &tangentAtGaussPts(0, 0), &tangentAtGaussPts(0, 1),
116  &tangentAtGaussPts(0, 2));
117 
118  FTensor::Tensor1<double, 3> t_coords_node0(cOords[0], cOords[1], cOords[2]);
119  FTensor::Tensor1<double, 3> t_coords_node1(cOords[3], cOords[4], cOords[5]);
121  t_dir(i) = t_coords_node1(i) - t_coords_node0(i);
122 
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_tangent(i) = t_dir(i);
127  ++t_tangent;
128  ++t_coords;
129  ++t_ksi;
130  }
131  }
133 }

◆ calculateEdgeDirection()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateEdgeDirection ( )
protected

Definition at line 35 of file EdgeElementForcesAndSourcesCore.cpp.

35  {
38  CHKERR mField.get_moab().get_connectivity(ent, cOnn, numNodes, true);
39  cOords.resize(numNodes * 3, false);
40  CHKERR mField.get_moab().get_coords(cOnn, numNodes, &*cOords.data().begin());
41  dIrection.resize(3, false);
42  cblas_dcopy(3, &cOords[3], 1, &*dIrection.data().begin(), 1);
43  cblas_daxpy(3, -1., &cOords[0], 1, &*dIrection.data().begin(), 1);
44  lEngth = cblas_dnrm2(3, &*dIrection.data().begin(), 1);
46 }

◆ 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 150 of file EdgeElementForcesAndSourcesCore.cpp.

150  {
152 
153  if (numeredEntFiniteElementPtr->getEntType() != MBEDGE)
155 
158  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
159  dataH1.dataOnEntities[MBEDGE][0].getSense() =
160  1; // set sense to 1, this is this entity
161 
162  // L2
163  if (dataH1.spacesOnEntities[MBEDGE].test(L2)) {
164  auto &data_l2 = *dataOnElement[L2];
165  CHKERR getEntityDataOrder<MBEDGE>(data_l2, L2);
166  data_l2.dataOnEntities[MBEDGE][0].getSense() =
167  1; // set sense to 1, this is this entity
168  data_l2.spacesOnEntities[MBEDGE].set(L2);
169  }
170 
171  // Hcurl
172  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
173  auto &data_curl = *dataOnElement[HCURL];
174  CHKERR getEntityDataOrder<MBEDGE>(data_curl, HCURL);
175  data_curl.dataOnEntities[MBEDGE][0].getSense() =
176  1; // set sense to 1, this is this entity
177  data_curl.spacesOnEntities[MBEDGE].set(HCURL);
178  }
179 
184 
185  // Iterate over operators
187 
189 }

◆ setIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::setIntegrationPts ( )
protected

Definition at line 48 of file EdgeElementForcesAndSourcesCore.cpp.

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

Friends And Related Function Documentation

◆ FaceElementForcesAndSourcesCoreOnSide

Definition at line 62 of file EdgeElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ cOnn

const EntityHandle* MoFEM::EdgeElementForcesAndSourcesCore::cOnn
protected

Definition at line 54 of file EdgeElementForcesAndSourcesCore.hpp.

◆ cOords

VectorDouble MoFEM::EdgeElementForcesAndSourcesCore::cOords
protected

Definition at line 56 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dIrection

VectorDouble MoFEM::EdgeElementForcesAndSourcesCore::dIrection
protected

Definition at line 55 of file EdgeElementForcesAndSourcesCore.hpp.

◆ lEngth

double& MoFEM::EdgeElementForcesAndSourcesCore::lEngth
protected

Definition at line 51 of file EdgeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::EdgeElementForcesAndSourcesCore::meshPositionsFieldName

Definition at line 33 of file EdgeElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::normalsAtGaussPts
protected

Definition at line 49 of file EdgeElementForcesAndSourcesCore.hpp.

◆ numNodes

int MoFEM::EdgeElementForcesAndSourcesCore::numNodes
protected

Definition at line 53 of file EdgeElementForcesAndSourcesCore.hpp.

◆ tangentAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::tangentAtGaussPts
protected

Definition at line 48 of file EdgeElementForcesAndSourcesCore.hpp.

◆ tFaceOrientation

FTensor::Tensor1< double, 3 > EdgeElementForcesAndSourcesCore::tFaceOrientation
static
Initial value:
{
0., 0., 1.}

Definition at line 44 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:460
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::ForcesAndSourcesCore::calBernsteinBezierBaseFunctionsOnElement
MoFEMErrorCode calBernsteinBezierBaseFunctionsOnElement()
Calculate Bernstein-Bezier base.
Definition: ForcesAndSourcesCore.cpp:1132
MoFEM::EdgeElementForcesAndSourcesCore::setIntegrationPts
MoFEMErrorCode setIntegrationPts()
Definition: EdgeElementForcesAndSourcesCore.cpp:48
MoFEM::ForcesAndSourcesCore::getMaxRowOrder
int getMaxRowOrder() const
Get max order of approximation for field in rows.
Definition: ForcesAndSourcesCore.cpp:131
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
NOBASE
@ NOBASE
Definition: definitions.h:59
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
L2
@ L2
field with C-1 continuity
Definition: definitions.h:88
MoFEM::EdgeElementForcesAndSourcesCore::tangentAtGaussPts
MatrixDouble tangentAtGaussPts
Definition: EdgeElementForcesAndSourcesCore.hpp:48
MoFEM::ForcesAndSourcesCore::setGaussPts
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: ForcesAndSourcesCore.cpp:1961
MoFEM::ForcesAndSourcesCore::getMaxDataOrder
int getMaxDataOrder() const
Get max order of approximation for data fields.
Definition: ForcesAndSourcesCore.cpp:120
N_MBEDGE0
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:105
MoFEM::EdgeElementForcesAndSourcesCore::cOnn
const EntityHandle * cOnn
Definition: EdgeElementForcesAndSourcesCore.hpp:54
MoFEM::EdgeElementForcesAndSourcesCore::calculateEdgeDirection
MoFEMErrorCode calculateEdgeDirection()
Definition: EdgeElementForcesAndSourcesCore.cpp:35
MoFEM::ForcesAndSourcesCore::getElementPolynomialBase
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.
Definition: ForcesAndSourcesCore.hpp:90
MoFEM::ForcesAndSourcesCore::getMaxColOrder
int getMaxColOrder() const
Get max order of approximation for field in columns.
Definition: ForcesAndSourcesCore.cpp:135
MoFEM::ForcesAndSourcesCore::calHierarchicalBaseFunctionsOnElement
MoFEMErrorCode calHierarchicalBaseFunctionsOnElement()
Calculate base functions.
Definition: ForcesAndSourcesCore.cpp:1112
MoFEM::EntitiesFieldData::spacesOnEntities
std::array< std::bitset< LASTSPACE >, MBMAXTYPE > spacesOnEntities
spaces on entity types
Definition: EntitiesFieldData.hpp:49
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::EdgeElementForcesAndSourcesCore::lEngth
double & lEngth
Definition: EdgeElementForcesAndSourcesCore.hpp:51
MoFEM::EdgeElementForcesAndSourcesCore::numNodes
int numNodes
Definition: EdgeElementForcesAndSourcesCore.hpp:53
QUAD_1D_TABLE
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
MoFEM::ForcesAndSourcesCore::elementMeasure
double elementMeasure
Definition: ForcesAndSourcesCore.hpp:545
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
N_MBEDGE1
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:106
MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
coordinated at gauss points
Definition: ForcesAndSourcesCore.hpp:544
MoFEM::ForcesAndSourcesCore::dataOnElement
const std::array< boost::shared_ptr< EntitiesFieldData >, LASTSPACE > dataOnElement
Entity data on element entity rows fields.
Definition: ForcesAndSourcesCore.hpp:461
MoFEM::ForcesAndSourcesCore::createDataOnElement
MoFEMErrorCode createDataOnElement(EntityType type)
Create a entity data on element object.
Definition: ForcesAndSourcesCore.cpp:1337
MoFEM::ForcesAndSourcesCore::ForcesAndSourcesCore
ForcesAndSourcesCore(Interface &m_field)
Definition: ForcesAndSourcesCore.cpp:40
MoFEM::EdgeElementForcesAndSourcesCore::cOords
VectorDouble cOords
Definition: EdgeElementForcesAndSourcesCore.hpp:56
QUAD_::npoints
int npoints
Definition: quad.h:29
MoFEM::ForcesAndSourcesCore::getRule
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: ForcesAndSourcesCore.cpp:1955
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::ForcesAndSourcesCore::getSpacesAndBaseOnEntities
MoFEMErrorCode getSpacesAndBaseOnEntities(EntitiesFieldData &data) const
Get field approximation space and base on entities.
Definition: ForcesAndSourcesCore.cpp:990
FTensor::Index< 'i', 3 >
MoFEM::ForcesAndSourcesCore::dataH1
EntitiesFieldData & dataH1
Definition: ForcesAndSourcesCore.hpp:471
QUAD_1D_TABLE_SIZE
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163
FTensor::Tensor0
Definition: Tensor0.hpp:16
HCURL
@ HCURL
field with continuous tangents
Definition: definitions.h:86
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::EdgeElementForcesAndSourcesCore::meshPositionsFieldName
std::string meshPositionsFieldName
Definition: EdgeElementForcesAndSourcesCore.hpp:33
MoFEM::EntitiesFieldData::dataOnEntities
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
Definition: EntitiesFieldData.hpp:57
MoFEM::ForcesAndSourcesCore::loopOverOperators
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
Definition: ForcesAndSourcesCore.cpp:1379
MoFEM::EdgeElementForcesAndSourcesCore::calculateCoordsAtIntegrationPts
MoFEMErrorCode calculateCoordsAtIntegrationPts()
Definition: EdgeElementForcesAndSourcesCore.cpp:101
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::EdgeElementForcesAndSourcesCore::dIrection
VectorDouble dIrection
Definition: EdgeElementForcesAndSourcesCore.hpp:55
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:24