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

Public Attributes

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

Additional Inherited Members

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

Detailed Description

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()

EdgeElementForcesAndSourcesCore::EdgeElementForcesAndSourcesCore ( Interface m_field)

Definition at line 78 of file EdgeElementForcesAndSourcesCore.cpp.

80  : ForcesAndSourcesCore(m_field),
81  meshPositionsFieldName("MESH_NODE_POSITIONS"),
85  boost::shared_ptr<BaseFunction>(new EdgePolynomialBase());
86 }
ForcesAndSourcesCore(Interface &m_field)
OpSetCovariantPiolaTransformOnEdge opCovariantTransform
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

Member Function Documentation

◆ calculateCoordsAtIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateCoordsAtIntegrationPts ( )

Definition at line 154 of file EdgeElementForcesAndSourcesCore.cpp.

154  {
156  const int nb_gauss_pts = gaussPts.size2();
157  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
158  if (nb_gauss_pts) {
160  &*gaussPts.data().begin());
162  &coordsAtGaussPts(0, 0), &coordsAtGaussPts(0, 1),
163  &coordsAtGaussPts(0, 2));
164  FTensor::Tensor1<double, 3> t_coords_node0(cOords[0], cOords[1], cOords[2]);
165  FTensor::Tensor1<double, 3> t_coords_node1(cOords[3], cOords[4], cOords[5]);
167  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
168  t_coords(i) = N_MBEDGE0(t_ksi) * t_coords_node0(i) +
169  N_MBEDGE1(t_ksi) * t_coords_node1(i);
170  ++t_coords;
171  ++t_ksi;
172  }
173  }
175 }
MatrixDouble gaussPts
Matrix of integration points.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:483
#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:490
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:79

◆ calculateEdgeDirection()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateEdgeDirection ( )

Definition at line 88 of file EdgeElementForcesAndSourcesCore.cpp.

88  {
91  CHKERR mField.get_moab().get_connectivity(ent, cOnn, numNodes, true);
92  cOords.resize(numNodes * 3, false);
93  CHKERR mField.get_moab().get_coords(cOnn, numNodes, &*cOords.data().begin());
94  dIrection.resize(3, false);
95  cblas_dcopy(3, &cOords[3], 1, &*dIrection.data().begin(), 1);
96  cblas_daxpy(3, -1., &cOords[0], 1, &*dIrection.data().begin(), 1);
97  lEngth = cblas_dnrm2(3, &*dIrection.data().begin(), 1);
99 }
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:459
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:578
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403

◆ calculateHoCoordsAtIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::calculateHoCoordsAtIntegrationPts ( )

Definition at line 178 of file EdgeElementForcesAndSourcesCore.cpp.

178  {
180 
181  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
182  dataPtr->get<FieldName_mi_tag>().end()) {
186  } else {
187  tangentAtGaussPts.resize(0, 3, false);
188  }
190 }
DataForcesAndSourcesCore & dataH1
MoFEMErrorCode getNodesFieldData(const boost::string_ref field_name, FEDofEntity_multiIndex &dofs, VectorDouble &nodes_data, VectorDofs &nodes_dofs, FieldSpace &space, FieldApproximationBase &base) const
Get field data on nodes.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
#define CHKERR
Inline error check.
Definition: definitions.h:578
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
MoFEMErrorCode getEntityFieldData(DataForcesAndSourcesCore &data, const std::string &field_name, const EntityType type_lo=MBVERTEX, const EntityType type_hi=MBPOLYHEDRON) const
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403

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

192  {
194 
195  if (numeredEntFiniteElementPtr->getEntType() != MBEDGE)
197 
199 
200  DataForcesAndSourcesCore &data_curl = *dataOnElement[HCURL];
201 
204  CHKERR getEntityFieldDataOrder<MBEDGE>(dataH1, H1);
205  dataH1.dataOnEntities[MBEDGE][0].getSense() =
206  1; // set sense to 1, this is this entity
207 
208  // Hcurl
209  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
210  CHKERR getEntityFieldDataOrder<MBEDGE>(data_curl, HCURL);
211  data_curl.dataOnEntities[MBEDGE][0].getSense() =
212  1; // set sense to 1, this is this entity
213  data_curl.spacesOnEntities[MBEDGE].set(HCURL);
214  }
215 
216 
221 
222  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
223  CHKERR opCovariantTransform.opRhs(data_curl);
224  }
225 
226  // Iterate over operators
228 
230 }
OpSetCovariantPiolaTransformOnEdge opCovariantTransform
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
const boost::shared_ptr< DataForcesAndSourcesCore > dataOnElement[LASTSPACE]
Entity data on element entity rows fields.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:490
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
field with continuous tangents
Definition: definitions.h:169
#define CHKERR
Inline error check.
Definition: definitions.h:578
virtual MoFEMErrorCode opRhs(DataForcesAndSourcesCore &data, const bool do_vertices, const bool do_edges, const bool do_quads, const bool do_tris, const bool do_tets, const bool do_prisms, const bool error_if_no_base=true)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
continuous field
Definition: definitions.h:168
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.

◆ setIntegrationPts()

MoFEMErrorCode EdgeElementForcesAndSourcesCore::setIntegrationPts ( )

Definition at line 101 of file EdgeElementForcesAndSourcesCore.cpp.

101  {
103  int order_data = getMaxDataOrder();
104  int order_row = getMaxRowOrder();
105  int order_col = getMaxColOrder();
106  int rule = getRule(order_row, order_col, order_data);
107 
108  int nb_gauss_pts;
109  if (rule >= 0) {
110  if (rule < QUAD_1D_TABLE_SIZE) {
111  if (QUAD_1D_TABLE[rule]->dim != 1) {
112  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
113  }
114  if (QUAD_1D_TABLE[rule]->order < rule) {
115  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
116  "wrong order %d != %d", QUAD_1D_TABLE[rule]->order, rule);
117  }
118  nb_gauss_pts = QUAD_1D_TABLE[rule]->npoints;
119  gaussPts.resize(2, nb_gauss_pts, false);
120  cblas_dcopy(nb_gauss_pts, &QUAD_1D_TABLE[rule]->points[1], 2,
121  &gaussPts(0, 0), 1);
122  cblas_dcopy(nb_gauss_pts, QUAD_1D_TABLE[rule]->weights, 1,
123  &gaussPts(1, 0), 1);
124  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 2,
125  false);
126  double *shape_ptr =
127  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
128  cblas_dcopy(2 * nb_gauss_pts, QUAD_1D_TABLE[rule]->points, 1, shape_ptr,
129  1);
130  } else {
131  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
132  "rule > quadrature order %d < %d", rule, QUAD_1D_TABLE_SIZE);
133  nb_gauss_pts = 0;
134  }
135  } else {
136  // If rule is negative, set user defined integration points
137  CHKERR setGaussPts(order_row, order_col, order_data);
138  nb_gauss_pts = gaussPts.size2();
139  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 2,
140  false);
141  if (nb_gauss_pts) {
142  for (int gg = 0; gg != nb_gauss_pts; gg++) {
143  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 0) =
144  N_MBEDGE0(gaussPts(0, gg));
145  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE)(gg, 1) =
146  N_MBEDGE1(gaussPts(0, gg));
147  }
148  }
149  }
151 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
MatrixDouble gaussPts
Matrix of integration points.
static QUAD *const QUAD_1D_TABLE[]
Definition: quad.h:164
DataForcesAndSourcesCore & dataH1
int npoints
Definition: quad.h:29
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:459
#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:578
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:403
virtual int getRule(int order)
set integration rule for finite element
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.
#define QUAD_1D_TABLE_SIZE
Definition: quad.h:163

Member Data Documentation

◆ cOnn

const EntityHandle* MoFEM::EdgeElementForcesAndSourcesCore::cOnn

Definition at line 53 of file EdgeElementForcesAndSourcesCore.hpp.

◆ cOords

VectorDouble MoFEM::EdgeElementForcesAndSourcesCore::cOords

Definition at line 55 of file EdgeElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::coordsAtGaussPts

Definition at line 56 of file EdgeElementForcesAndSourcesCore.hpp.

◆ dIrection

VectorDouble MoFEM::EdgeElementForcesAndSourcesCore::dIrection

Definition at line 54 of file EdgeElementForcesAndSourcesCore.hpp.

◆ lEngth

double MoFEM::EdgeElementForcesAndSourcesCore::lEngth

Definition at line 50 of file EdgeElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::EdgeElementForcesAndSourcesCore::meshPositionsFieldName

Definition at line 42 of file EdgeElementForcesAndSourcesCore.hpp.

◆ numNodes

int MoFEM::EdgeElementForcesAndSourcesCore::numNodes

Definition at line 52 of file EdgeElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnEdge MoFEM::EdgeElementForcesAndSourcesCore::opCovariantTransform

Definition at line 46 of file EdgeElementForcesAndSourcesCore.hpp.

◆ opGetHoTangentOnEdge

OpGetHoTangentOnEdge MoFEM::EdgeElementForcesAndSourcesCore::opGetHoTangentOnEdge

Definition at line 45 of file EdgeElementForcesAndSourcesCore.hpp.

◆ tangentAtGaussPts

MatrixDouble MoFEM::EdgeElementForcesAndSourcesCore::tangentAtGaussPts

Definition at line 44 of file EdgeElementForcesAndSourcesCore.hpp.


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