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

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

#include <src/finite_elements/FaceElementForcesAndSourcesCore.hpp>

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

Classes

struct  UserDataOperator
 default operator for TRI element More...
 

Public Types

enum  Switches { NO_HO_GEOMETRY = 1 << 0, NO_CONTRAVARIANT_TRANSFORM_HDIV = 1 << 1, 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::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

 FaceElementForcesAndSourcesCoreBase (Interface &m_field)
 
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...
 
MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 Get number of DOFs on element. 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 ()
 
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)
 Copy data from other base method to this base method. More...
 
- 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
 Name of the field with geometry. More...
 
- 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
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexrowPtr
 Pointer to finite element rows dofs view. More...
 
boost::shared_ptr< const FENumeredDofEntity_multiIndexcolPtr
 Pointer to finite element columns dofs view. More...
 
boost::shared_ptr< const FEDofEntity_multiIndexdataPtr
 Pointer to finite element data dofs. More...
 
boost::shared_ptr< const FieldEntity_vector_viewrowFieldEntsPtr
 Pointer to finite element field entities row view. More...
 
boost::shared_ptr< const FieldEntity_vector_viewcolFieldEntsPtr
 Pointer to finite element field entities column view. More...
 
boost::shared_ptr< const FieldEntity_multiIndex_spaceType_viewdataFieldEntsPtr
 Pointer to finite element field entities data view. More...
 
- 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
 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
 

Protected Member Functions

MoFEMErrorCode getNumberOfNodes (int &num_nodes) const
 
virtual MoFEMErrorCode calculateAreaAndNormal ()
 Calculate element area and normal of the face. More...
 
virtual MoFEMErrorCode setIntegrationPts ()
 Set integration points. More...
 
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement ()
 Determine approximation space and order of base functions. More...
 
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts ()
 Calculate coordinate at integration points. More...
 
virtual MoFEMErrorCode calculateHoNormal ()
 Calculate normal on curved elements. More...
 
- 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 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...
 
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 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

double aRea
 
int num_nodes
 
const EntityHandleconn
 
VectorDouble nOrmal
 
VectorDouble tangentOne
 
VectorDouble tangentTwo
 
VectorDouble coords
 
MatrixDouble coordsAtGaussPts
 
MatrixDouble hoCoordsAtGaussPts
 
MatrixDouble normalsAtGaussPts
 
MatrixDouble tangentOneAtGaussPts
 
MatrixDouble tangentTwoAtGaussPts
 
OpGetCoordsAndNormalsOnFace opHOCoordsAndNormals
 
OpSetContravariantPiolaTransformOnFace opContravariantTransform
 
OpSetCovariantPiolaTransformOnFace opCovariantTransform
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore
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
 
boost::ptr_vector< UserDataOperatoropPtrVector
 Vector of finite element users data operators. More...
 

Friends

class UserDataOperator
 
class VolumeElementForcesAndSourcesCoreOnSideBase
 

Detailed Description

Face finite element

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

Definition at line 35 of file FaceElementForcesAndSourcesCore.hpp.

Member Enumeration Documentation

◆ Switches

Constructor & Destructor Documentation

◆ FaceElementForcesAndSourcesCoreBase()

MoFEM::FaceElementForcesAndSourcesCoreBase::FaceElementForcesAndSourcesCoreBase ( Interface m_field)

Definition at line 23 of file FaceElementForcesAndSourcesCore.cpp.

25  : ForcesAndSourcesCore(m_field),
26  meshPositionsFieldName("MESH_NODE_POSITIONS"),
34  boost::shared_ptr<BaseFunction>(new TriPolynomialBase());
35 }
ForcesAndSourcesCore(Interface &m_field)
std::string meshPositionsFieldName
Name of the field with geometry.
OpSetContravariantPiolaTransformOnFace opContravariantTransform
auto & getElementPolynomialBase()
Get the Entity Polynomial Base object.

Member Function Documentation

◆ calculateAreaAndNormal()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateAreaAndNormal ( )
protectedvirtual

Calculate element area and normal of the face.

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

Returns
Error code

Definition at line 86 of file FaceElementForcesAndSourcesCore.cpp.

86  {
89  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
90  coords.resize(num_nodes * 3, false);
91  CHKERR mField.get_moab().get_coords(conn, num_nodes, &*coords.data().begin());
92  double diff_n[6];
93  CHKERR ShapeDiffMBTRI(diff_n);
94 
95  nOrmal.resize(3, false);
96  tangentOne.resize(3, false);
97  tangentTwo.resize(3, false);
99  &coords[0], &coords[1], &coords[2]);
101  &nOrmal[0], &nOrmal[1], &nOrmal[2]);
103  &tangentOne[0], &tangentOne[1], &tangentOne[2]);
105  &tangentTwo[0], &tangentTwo[1], &tangentTwo[2]);
107  &diff_n[1]);
108 
114  t_t1(i) = 0;
115  t_t2(i) = 0;
116  for (int nn = 0; nn != 3; ++nn) {
117  t_t1(i) += t_coords(i) * t_diff(N0);
118  t_t2(i) += t_coords(i) * t_diff(N1);
119  ++t_coords;
120  ++t_diff;
121  }
122  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
123  aRea = sqrt(t_normal(i) * t_normal(i)) / 2.;
124 
126 }
virtual moab::Interface & get_moab()=0
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
PetscErrorCode ShapeDiffMBTRI(double *diffN)
calculate derivatives of shape functions
Definition: fem_tools.c:206
#define CHKERR
Inline error check.
Definition: definitions.h:596
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use

◆ calculateCoordinatesAtGaussPts()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateCoordinatesAtGaussPts ( )
protectedvirtual

Calculate coordinate at integration points.

Returns
Error code

Definition at line 224 of file FaceElementForcesAndSourcesCore.cpp.

224  {
226 
227  double *shape_functions =
228  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
229  const size_t nb_gauss_pts = gaussPts.size2();
230  coordsAtGaussPts.resize(nb_gauss_pts, 3, false);
231  for (int gg = 0; gg != nb_gauss_pts; ++gg)
232  for (int dd = 0; dd != 3; ++dd)
233  coordsAtGaussPts(gg, dd) =
234  cblas_ddot(3, &shape_functions[3 * gg], 1, &coords[dd], 3);
235 
237 }
MatrixDouble gaussPts
Matrix of integration points.
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:501
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:508
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33

◆ calculateHoNormal()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::calculateHoNormal ( )
protectedvirtual

Calculate normal on curved elements.

Geometry is given by other field.

Returns
error code

Definition at line 239 of file FaceElementForcesAndSourcesCore.cpp.

239  {
241  // Check if field for high-order geometry is set and if it is set calculate
242  // higher-order normals and face tangent vectors.
243  if (dataPtr->get<FieldName_mi_tag>().find(meshPositionsFieldName) !=
244  dataPtr->get<FieldName_mi_tag>().end()) {
245 
246  const Field *field_struture =
248  BitFieldId id = field_struture->getId();
249 
250  if ((numeredEntFiniteElementPtr->getBitFieldIdData() & id).none())
251  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_FOUND,
252  "no MESH_NODE_POSITIONS in element data");
253 
254  // Calculate normal for high-order geometry
255 
261 
262  } else {
263  hoCoordsAtGaussPts.resize(0, 0, false);
264  normalsAtGaussPts.resize(0, 0, false);
265  tangentOneAtGaussPts.resize(0, 0, false);
266  tangentTwoAtGaussPts.resize(0, 0, false);
267  }
269 }
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:477
std::string meshPositionsFieldName
Name of the field with geometry.
boost::shared_ptr< const FEDofEntity_multiIndex > dataPtr
Pointer to finite element data dofs.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
std::bitset< BITFIELDID_SIZE > BitFieldId
Field Id.
Definition: Types.hpp:53
#define CHKERR
Inline error check.
Definition: definitions.h:596
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:407
virtual const Field * get_field_structure(const std::string &name)=0
get field structure

◆ getNumberOfNodes()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::getNumberOfNodes ( int num_nodes) const
protected

◆ getSpaceBaseAndOrderOnElement()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::getSpaceBaseAndOrderOnElement ( )
protectedvirtual

Determine approximation space and order of base functions.

Returns
Error code

Definition at line 178 of file FaceElementForcesAndSourcesCore.cpp.

178  {
180  // Get spaces order/base and sense of entities.
181 
183 
184  // H1
185  if (dataH1.spacesOnEntities[MBEDGE].test(H1)) {
186  CHKERR getEntitySense<MBEDGE>(dataH1);
187  CHKERR getEntityDataOrder<MBEDGE>(dataH1, H1);
188  }
189  if (dataH1.spacesOnEntities[MBTRI].test(H1)) {
190  CHKERR getEntitySense<MBTRI>(dataH1);
191  CHKERR getEntityDataOrder<MBTRI>(dataH1, H1);
192  }
193 
194  // Hcurl
195  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL)) {
196  CHKERR getEntitySense<MBEDGE>(dataHcurl);
197  CHKERR getEntityDataOrder<MBEDGE>(dataHcurl, HCURL);
198  dataHcurl.spacesOnEntities[MBEDGE].set(HCURL);
199  }
200  if (dataH1.spacesOnEntities[MBTRI].test(HCURL)) {
201  CHKERR getEntitySense<MBTRI>(dataHcurl);
202  CHKERR getEntityDataOrder<MBTRI>(dataHcurl, HCURL);
203  dataHcurl.spacesOnEntities[MBTRI].set(HCURL);
204  }
205 
206  // Hdiv
207  if (dataH1.spacesOnEntities[MBTRI].test(HDIV)) {
208  CHKERR getEntitySense<MBTRI>(dataHdiv);
209  CHKERR getEntityDataOrder<MBTRI>(dataHdiv, HDIV);
210  dataHdiv.spacesOnEntities[MBTRI].set(HDIV);
211  }
212 
213  // L2
214  if (dataH1.spacesOnEntities[MBTRI].test(L2)) {
215  CHKERR getEntitySense<MBTRI>(dataL2);
216  CHKERR getEntityDataOrder<MBTRI>(dataL2, L2);
217  dataL2.spacesOnEntities[MBTRI].set(L2);
218  }
219 
221 }
DataForcesAndSourcesCore & dataL2
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
field with continuous normal traction
Definition: definitions.h:173
DataForcesAndSourcesCore & dataH1
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
DataForcesAndSourcesCore & dataHcurl
field with continuous tangents
Definition: definitions.h:172
#define CHKERR
Inline error check.
Definition: definitions.h:596
DataForcesAndSourcesCore & dataHdiv
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
continuous field
Definition: definitions.h:171
MoFEMErrorCode getSpacesAndBaseOnEntities(DataForcesAndSourcesCore &data) const
Get field approximation space and base on entities.
field with C-1 continuity
Definition: definitions.h:174

◆ OpSwitch()

template<int SWITCH>
MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::OpSwitch ( )

Definition at line 322 of file FaceElementForcesAndSourcesCore.hpp.

322  {
324 
325  if (numeredEntFiniteElementPtr->getEntType() != MBTRI)
328 
329  // Calculate normal and tangent vectors for face geometry given by 3 nodes.
332 
334  if (gaussPts.size2() == 0)
336 
337  DataForcesAndSourcesCore &data_curl = *dataOnElement[HCURL];
338  DataForcesAndSourcesCore &data_div = *dataOnElement[HDIV];
339 
340  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).resize(3, 2, false);
341  std::copy(
343  dataH1.dataOnEntities[MBVERTEX][0].getDiffN(NOBASE).data().begin());
344 
345  /// Use the some node base
348  if (!(NO_HO_GEOMETRY & SWITCH))
350 
351  // Apply Piola transform to HDiv and HCurl spaces, uses previously
352  // calculated faces normal and tangent vectors.
353  if (!(NO_CONTRAVARIANT_TRANSFORM_HDIV & SWITCH))
354  if (dataH1.spacesOnEntities[MBTRI].test(HDIV))
356 
357  if (!(NO_COVARIANT_TRANSFORM_HCURL & SWITCH))
358  if (dataH1.spacesOnEntities[MBEDGE].test(HCURL))
359  CHKERR opCovariantTransform.opRhs(data_curl);
360 
361  // Iterate over operators
363 
365  }
std::bitset< LASTSPACE > spacesOnEntities[MBMAXTYPE]
spaces on entity types
MoFEMErrorCode calculateBaseFunctionsOnElement()
Calculate base functions.
MatrixDouble gaussPts
Matrix of integration points.
field with continuous normal traction
Definition: definitions.h:173
MoFEMErrorCode createDataOnElement()
Create a entity data on element object.
DataForcesAndSourcesCore & dataH1
virtual MoFEMErrorCode getSpaceBaseAndOrderOnElement()
Determine approximation space and order of base functions.
virtual MoFEMErrorCode setIntegrationPts()
Set integration points.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
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:508
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
OpSetContravariantPiolaTransformOnFace opContravariantTransform
field with continuous tangents
Definition: definitions.h:172
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:83
#define CHKERR
Inline error check.
Definition: definitions.h:596
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)
virtual MoFEMErrorCode calculateCoordinatesAtGaussPts()
Calculate coordinate at integration points.
virtual MoFEMErrorCode calculateHoNormal()
Calculate normal on curved elements.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual MoFEMErrorCode calculateAreaAndNormal()
Calculate element area and normal of the face.

◆ setIntegrationPts()

MoFEMErrorCode MoFEM::FaceElementForcesAndSourcesCoreBase::setIntegrationPts ( )
protectedvirtual

Set integration points.

Returns
Error code

Definition at line 128 of file FaceElementForcesAndSourcesCore.cpp.

128  {
130  // Set integration points
131  int order_data = getMaxDataOrder();
132  int order_row = getMaxRowOrder();
133  int order_col = getMaxColOrder();
134  int rule = getRule(order_row, order_col, order_data);
135 
136  if (rule >= 0) {
137  if (rule < QUAD_2D_TABLE_SIZE) {
138  if (QUAD_2D_TABLE[rule]->dim != 2) {
139  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "wrong dimension");
140  }
141  if (QUAD_2D_TABLE[rule]->order < rule) {
142  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
143  "wrong order %d != %d", QUAD_2D_TABLE[rule]->order, rule);
144  }
145  const size_t nb_gauss_pts = QUAD_2D_TABLE[rule]->npoints;
146  gaussPts.resize(3, nb_gauss_pts, false);
147  cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[1], 3,
148  &gaussPts(0, 0), 1);
149  cblas_dcopy(nb_gauss_pts, &QUAD_2D_TABLE[rule]->points[2], 3,
150  &gaussPts(1, 0), 1);
151  cblas_dcopy(nb_gauss_pts, QUAD_2D_TABLE[rule]->weights, 1, &gaussPts(2, 0),
152  1);
153  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
154  false);
155  double *shape_ptr =
156  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin();
157  cblas_dcopy(3 * nb_gauss_pts, QUAD_2D_TABLE[rule]->points, 1, shape_ptr, 1);
158  } else {
159  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
160  "rule > quadrature order %d < %d", rule, QUAD_2D_TABLE_SIZE);
161  }
162  } else {
163  // If rule is negative, set user defined integration points
164  CHKERR setGaussPts(order_row, order_col, order_data);
165  const size_t nb_gauss_pts = gaussPts.size2();
166  dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).resize(nb_gauss_pts, 3,
167  false);
168  if (nb_gauss_pts) {
170  &*dataH1.dataOnEntities[MBVERTEX][0].getN(NOBASE).data().begin(),
171  &gaussPts(0, 0), &gaussPts(1, 0), nb_gauss_pts);
172  }
173  }
175 }
int getMaxRowOrder() const
Get max order of approximation for field in rows.
MatrixDouble gaussPts
Matrix of integration points.
DataForcesAndSourcesCore & dataH1
int npoints
Definition: quad.h:29
void cblas_dcopy(const int N, const double *X, const int incX, double *Y, const int incY)
Definition: cblas_dcopy.c:11
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
virtual MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
boost::ptr_vector< EntData > dataOnEntities[MBMAXTYPE]
#define QUAD_2D_TABLE_SIZE
Definition: quad.h:174
static QUAD *const QUAD_2D_TABLE[]
Definition: quad.h:175
#define CHKERR
Inline error check.
Definition: definitions.h:596
PetscErrorCode ShapeMBTRI(double *N, const double *X, const double *Y, const int G_DIM)
calculate shape functions on triangle
Definition: fem_tools.c:194
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407
virtual int getRule(int order_row, int order_col, int order_data)
another variant of getRule
int getMaxColOrder() const
Get max order of approximation for field in columns.
int getMaxDataOrder() const
Get max order of approximation for data fields.
const int order
approximation order

Friends And Related Function Documentation

◆ UserDataOperator

friend class UserDataOperator
friend

Definition at line 422 of file FaceElementForcesAndSourcesCore.hpp.

◆ VolumeElementForcesAndSourcesCoreOnSideBase

Definition at line 423 of file FaceElementForcesAndSourcesCore.hpp.

Member Data Documentation

◆ aRea

double MoFEM::FaceElementForcesAndSourcesCoreBase::aRea
protected

Definition at line 407 of file FaceElementForcesAndSourcesCore.hpp.

◆ conn

const EntityHandle* MoFEM::FaceElementForcesAndSourcesCoreBase::conn
protected

Definition at line 409 of file FaceElementForcesAndSourcesCore.hpp.

◆ coords

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::coords
protected

Definition at line 411 of file FaceElementForcesAndSourcesCore.hpp.

◆ coordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::coordsAtGaussPts
protected

Definition at line 412 of file FaceElementForcesAndSourcesCore.hpp.

◆ hoCoordsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::hoCoordsAtGaussPts
protected

Definition at line 414 of file FaceElementForcesAndSourcesCore.hpp.

◆ meshPositionsFieldName

std::string MoFEM::FaceElementForcesAndSourcesCoreBase::meshPositionsFieldName

Name of the field with geometry.

Definition at line 37 of file FaceElementForcesAndSourcesCore.hpp.

◆ nOrmal

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::nOrmal
protected

Definition at line 410 of file FaceElementForcesAndSourcesCore.hpp.

◆ normalsAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::normalsAtGaussPts
protected

Definition at line 415 of file FaceElementForcesAndSourcesCore.hpp.

◆ num_nodes

int MoFEM::FaceElementForcesAndSourcesCoreBase::num_nodes
protected

Definition at line 408 of file FaceElementForcesAndSourcesCore.hpp.

◆ opContravariantTransform

OpSetContravariantPiolaTransformOnFace MoFEM::FaceElementForcesAndSourcesCoreBase::opContravariantTransform
protected

Definition at line 419 of file FaceElementForcesAndSourcesCore.hpp.

◆ opCovariantTransform

OpSetCovariantPiolaTransformOnFace MoFEM::FaceElementForcesAndSourcesCoreBase::opCovariantTransform
protected

Definition at line 420 of file FaceElementForcesAndSourcesCore.hpp.

◆ opHOCoordsAndNormals

OpGetCoordsAndNormalsOnFace MoFEM::FaceElementForcesAndSourcesCoreBase::opHOCoordsAndNormals
protected

Definition at line 418 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOne

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentOne
protected

Definition at line 410 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentOneAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentOneAtGaussPts
protected

Definition at line 416 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwo

VectorDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentTwo
protected

Definition at line 410 of file FaceElementForcesAndSourcesCore.hpp.

◆ tangentTwoAtGaussPts

MatrixDouble MoFEM::FaceElementForcesAndSourcesCoreBase::tangentTwoAtGaussPts
protected

Definition at line 417 of file FaceElementForcesAndSourcesCore.hpp.


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