v0.8.23
Public Member Functions | List of all members
MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator Struct Reference

default operator for EDGE element More...

#include <src/finite_elements/EdgeElementForcesAndSourcesCore.hpp>

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

Public Member Functions

const EntityHandlegetConn ()
 get element connectivity More...
 
double getLength ()
 get edge length More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetDirection ()
 get edge direction More...
 
VectorDoublegetCoords ()
 get edge node coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 get coordinate at integration point More...
 
auto getFTensor1CoordsAtGaussPts ()
 get coordinates at Gauss pts. More...
 
MatrixDoublegetTangetAtGaussPts ()
 get tangent vector to edge curve at integration points More...
 
const EdgeElementForcesAndSourcesCoreBasegetEdgeFE ()
 get pointer to this finite element More...
 
FTensor::Tensor1< double, 3 > getTensor1Direction ()
 
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getTensor1Coords ()
 get get coords at gauss points More...
 
FTensor::Tensor1< FTensor::PackPtr< double *, 3 >, 3 > getTensor1TangentAtGaussPts ()
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPLAST, const bool symm=true)
 
 UserDataOperator (const std::string &field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
virtual ~UserDataOperator ()
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
boost::weak_ptr< SideNumber > getSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement ()
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
const std::string & getFEName () const
 Get name of the element. More...
 
Vec getSnesF () const
 
Vec getSnesX () const
 
Mat getSnesA () const
 
Mat getSnesB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTSa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
DEPRECATED MoFEMErrorCode getPorblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true, const bool do_vertices=true, const bool do_edges=true, const bool do_quads=true, const bool do_tris=true, const bool do_tets=true, const bool do_prisms=true)
 
virtual ~DataOperator ()
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, DataForcesAndSourcesCore::EntData &row_data, DataForcesAndSourcesCore::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data, bool symm=true)
 
virtual MoFEMErrorCode opLhs (DataForcesAndSourcesCore &row_data, DataForcesAndSourcesCore &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
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 opRhs (DataForcesAndSourcesCore &data, const bool error_if_no_base=true)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType { OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPLAST = 1 << 3 }
 Controls loop over entities on element. More...
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
bool doVertices
 If false skip vertices. More...
 
bool doEdges
 If false skip edges. More...
 
bool doQuads
 
bool doTris
 
bool doTets
 
bool doPrisms
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

default operator for EDGE element

Examples
forces_and_sources_testing_edge_element.cpp, and hcurl_divergence_operator_2d.cpp.

Definition at line 49 of file EdgeElementForcesAndSourcesCore.hpp.

Member Function Documentation

◆ getConn()

const EntityHandle* MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getConn ( )

get element connectivity

Definition at line 55 of file EdgeElementForcesAndSourcesCore.hpp.

55  {
56  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE)->cOnn;
57  }

◆ getCoords()

VectorDouble& MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getCoords ( )

get edge node coordinates

Definition at line 83 of file EdgeElementForcesAndSourcesCore.hpp.

83  {
84  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE)->cOords;
85  }

◆ getCoordsAtGaussPts()

MatrixDouble& MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getCoordsAtGaussPts ( )

get coordinate at integration point

Definition at line 90 of file EdgeElementForcesAndSourcesCore.hpp.

90  {
91  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE)
92  ->coordsAtGaussPts;
93  }

◆ getDirection()

VectorDouble& MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getDirection ( )

get edge direction

Definition at line 75 of file EdgeElementForcesAndSourcesCore.hpp.

75  {
76  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE)
77  ->dIrection;
78  }

◆ getEdgeFE()

const EdgeElementForcesAndSourcesCoreBase* MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getEdgeFE ( )

get pointer to this finite element

Definition at line 114 of file EdgeElementForcesAndSourcesCore.hpp.

114  {
115  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE);
116  }

◆ getFTensor1CoordsAtGaussPts()

auto MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getFTensor1CoordsAtGaussPts ( )

get coordinates at Gauss pts.

Definition at line 97 of file EdgeElementForcesAndSourcesCore.hpp.

97  {
98  double *ptr = &*getCoordsAtGaussPts().data().begin();
99  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
100  &ptr[2]);
101  }
MatrixDouble & getCoordsAtGaussPts()
get coordinate at integration point

◆ getLength()

double MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getLength ( )

get edge length

Definition at line 62 of file EdgeElementForcesAndSourcesCore.hpp.

62  {
63  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE)->lEngth;
64  }

◆ getMeasure()

double MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getMeasure ( )

get measure of element

Returns
length of face

Definition at line 70 of file EdgeElementForcesAndSourcesCore.hpp.

◆ getTangetAtGaussPts()

MatrixDouble& MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getTangetAtGaussPts ( )

get tangent vector to edge curve at integration points

Definition at line 106 of file EdgeElementForcesAndSourcesCore.hpp.

106  {
107  return static_cast<EdgeElementForcesAndSourcesCoreBase *>(ptrFE)
108  ->tangentAtGaussPts;
109  }

◆ getTensor1Coords()

FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getTensor1Coords ( )

get get coords at gauss points

auto t_center;
auto t_coords = getTensor1Coords();
t_center(i) = 0;
for(int nn = 0;nn!=2;nn++) {
t_center(i) += t_coords(i);
++t_coords;
}
t_center(i) /= 2;

Definition at line 140 of file EdgeElementForcesAndSourcesCore.hpp.

140  {
141  double *ptr = &*getCoords().data().begin();
142  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
143  &ptr[2]);
144  }

◆ getTensor1Direction()

FTensor::Tensor1<double, 3> MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getTensor1Direction ( )

◆ getTensor1TangentAtGaussPts()

FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3> MoFEM::EdgeElementForcesAndSourcesCoreBase::UserDataOperator::getTensor1TangentAtGaussPts ( )

Definition at line 147 of file EdgeElementForcesAndSourcesCore.hpp.

147  {
148  double *ptr = &*getTangetAtGaussPts().data().begin();
149  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
150  &ptr[2]);
151  }
MatrixDouble & getTangetAtGaussPts()
get tangent vector to edge curve at integration points

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