v0.14.0
Public Member Functions | Protected Member Functions | List of all members
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator Struct Reference

default operator for TET element More...

#include <src/finite_elements/VolumeElementForcesAndSourcesCoreOnSide.hpp>

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

Public Member Functions

VolumeElementForcesAndSourcesCoreOnSidegetVolumeFE () const
 
FaceElementForcesAndSourcesCoregetFaceFE () const
 
DEPRECATED int getFaceSense () const
 get face sense in respect to volume More...
 
int getSkeletonSense () const
 Get the skeleton sense. More...
 
int getFaceSideNumber () const
 get face side number in respect to volume More...
 
bool getEdgeFace (const int ee) const
 
VectorDoublegetNormal ()
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
MatrixDoublegetFaceCoordsAtGaussPts ()
 get face coordinates at Gauss pts. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, 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)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (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 () const
 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...
 
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...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
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...
 

Protected Member Functions

MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0, OPCOL = 1 << 1, OPROWCOL = 1 << 2, OPSPACE = 1 << 3,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > >>
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

default operator for TET element

Examples
continuity_check_on_skeleton_3d.cpp, and hello_world.cpp.

Definition at line 97 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

Member Function Documentation

◆ getEdgeFace()

bool UserDataOperator::getEdgeFace ( const int  ee) const
inline

Definition at line 303 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

304  {
305  constexpr bool edges_on_faces[6][4] = {{true, false, false, true}, // e0
306  {false, true, false, true}, // e1
307  {false, false, true, true}, // e2
308  {true, false, true, false}, // e3
309  {true, true, false, false}, // e4
310  {false, true, true, false}};
311  return edges_on_faces[ee][getFaceSideNumber()];
312 }

◆ getFaceCoordsAtGaussPts()

MatrixDouble & UserDataOperator::getFaceCoordsAtGaussPts ( )
inline

get face coordinates at Gauss pts.

Note
Coordinates should be the same what function getCoordsAtGaussPts on tets is returning. If both coordinates are different it is error, or you do something very unusual.

Definition at line 299 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

299  {
300  return getFaceFE()->coordsAtGaussPts;
301 }

◆ getFaceFE()

FaceElementForcesAndSourcesCore * UserDataOperator::getFaceFE ( ) const
inline

Definition at line 226 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

227  {
228  return static_cast<FaceElementForcesAndSourcesCore *>(
230 }

◆ getFaceSense()

int UserDataOperator::getFaceSense ( ) const
inline

get face sense in respect to volume

Deprecated:
use getSkeletonSense
Returns
error code

Definition at line 254 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

254  {
255  return getSkeletonSense();
256 }

◆ getFaceSideNumber()

int UserDataOperator::getFaceSideNumber ( ) const
inline

get face side number in respect to volume

Returns
error code

Definition at line 289 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

289  {
290  return getVolumeFE()->faceSideNumber;
291 }

◆ getFTensor1Normal()

auto UserDataOperator::getFTensor1Normal ( )
inline

get normal as tensor

Definition at line 264 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

264  {
265  double *ptr = &*getNormal().data().begin();
266  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
267 }

◆ getFTensor1NormalsAtGaussPts()

auto UserDataOperator::getFTensor1NormalsAtGaussPts ( )
inline

get normal at integration points

Example:

double nrm2;
auto t_normal = getFTensor1NormalsAtGaussPts();
for(int gg = gg!=data.getN().size1();gg++) {
nrm2 = sqrt(t_normal(i)*t_normal(i));
++t_normal;
}

Definition at line 282 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

282  {
283  double *ptr = &*getNormalsAtGaussPts().data().begin();
284  return FTensor::Tensor1<FTensor::PackPtr<double *, 3>, 3>(ptr, &ptr[1],
285  &ptr[2]);
286 }

◆ getFTensor1Tangent1()

auto UserDataOperator::getFTensor1Tangent1 ( )
inline

get tangentOne as tensor

Definition at line 270 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

270  {
271  double *ptr = &*getTangent1().data().begin();
272  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
273 }

◆ getFTensor1Tangent2()

auto UserDataOperator::getFTensor1Tangent2 ( )
inline

get tangentTwo as tensor

Definition at line 276 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

276  {
277  double *ptr = &*getTangent2().data().begin();
278  return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
279 }

◆ getNormal()

VectorDouble & UserDataOperator::getNormal ( )
inline

get face normal on side which is this element

Returns
face normal

Definition at line 233 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

233  {
234  return getFaceFE()->nOrmal;
235 }

◆ getNormalsAtGaussPts() [1/2]

MatrixDouble & UserDataOperator::getNormalsAtGaussPts ( )
inline

if higher order geometry return normals at Gauss pts.

Note: returned matrix has size 0 in rows and columns if no HO approximation of geometry is available.

Definition at line 294 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

294  {
295  return getFaceFE()->normalsAtGaussPts;
296 }

◆ getNormalsAtGaussPts() [2/2]

ublas::matrix_row< MatrixDouble > UserDataOperator::getNormalsAtGaussPts ( const int  gg)
inline

if higher order geometry return normals at Gauss pts.

Parameters
gggauss point number

Definition at line 315 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

315  {
316  return ublas::matrix_row<MatrixDouble>(getNormalsAtGaussPts(), gg);
317 }

◆ getSkeletonSense()

int UserDataOperator::getSkeletonSense ( ) const
inline

Get the skeleton sense.

calls getFaceSense()

Returns
int

Definition at line 259 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

259  {
260  return getVolumeFE()->faceSense;
261 }

◆ getTangent1()

VectorDouble & UserDataOperator::getTangent1 ( )
inline

get triangle tangent 1

Definition at line 238 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

238  {
239  return getFaceFE()->tangentOne;
240 }

◆ getTangent2()

VectorDouble & UserDataOperator::getTangent2 ( )
inline

get triangle tangent 2

Definition at line 243 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

243  {
244  return getFaceFE()->tangentTwo;
245 }

◆ getVolumeFE()

VolumeElementForcesAndSourcesCoreOnSide * UserDataOperator::getVolumeFE ( ) const
inline

Definition at line 248 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

249  {
250  return static_cast<VolumeElementForcesAndSourcesCoreOnSide *>(ptrFE);
251 }

◆ setPtrFE()

MoFEMErrorCode UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Reimplemented from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator.

Definition at line 98 of file VolumeElementForcesAndSourcesCoreOnSide.cpp.

99  {
101  if (!(ptrFE =
102  dynamic_cast<VolumeElementForcesAndSourcesCoreOnSide *>(ptr)))
103  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
104  "User operator and finite element do not work together");
106 }

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
MoFEM::FaceElementForcesAndSourcesCore::tangentOne
VectorDouble tangentOne
Definition: FaceElementForcesAndSourcesCore.hpp:79
MoFEM::FaceElementForcesAndSourcesCore::tangentTwo
VectorDouble tangentTwo
Definition: FaceElementForcesAndSourcesCore.hpp:79
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MoFEM::ForcesAndSourcesCore::sidePtrFE
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
Definition: ForcesAndSourcesCore.hpp:507
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getTangent1
VectorDouble & getTangent1()
get triangle tangent 1
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:238
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFTensor1NormalsAtGaussPts
auto getFTensor1NormalsAtGaussPts()
get normal at integration points
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:282
MoFEM::ForcesAndSourcesCore::coordsAtGaussPts
MatrixDouble coordsAtGaussPts
coordinated at gauss points
Definition: ForcesAndSourcesCore.hpp:544
MoFEM::ForcesAndSourcesCore::VolumeElementForcesAndSourcesCoreOnSide
friend class VolumeElementForcesAndSourcesCoreOnSide
Definition: ForcesAndSourcesCore.hpp:534
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getNormalsAtGaussPts
MatrixDouble & getNormalsAtGaussPts()
if higher order geometry return normals at Gauss pts.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:294
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::faceSideNumber
int faceSideNumber
Face side number.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:88
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getVolumeFE
VolumeElementForcesAndSourcesCoreOnSide * getVolumeFE() const
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:248
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceSideNumber
int getFaceSideNumber() const
get face side number in respect to volume
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:289
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceFE
FaceElementForcesAndSourcesCore * getFaceFE() const
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:226
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FaceElementForcesAndSourcesCore
FTensor::Index< 'i', 3 >
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getTangent2
VectorDouble & getTangent2()
get triangle tangent 2
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:243
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getNormal
VectorDouble & getNormal()
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:233
MoFEM::FaceElementForcesAndSourcesCore::nOrmal
VectorDouble nOrmal
Definition: FaceElementForcesAndSourcesCore.hpp:79
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::faceSense
int faceSense
Sense of face, could be 1 or -1.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:87
MoFEM::FaceElementForcesAndSourcesCore::normalsAtGaussPts
MatrixDouble normalsAtGaussPts
Definition: FaceElementForcesAndSourcesCore.hpp:82
MoFEM::VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getSkeletonSense
int getSkeletonSense() const
Get the skeleton sense.
Definition: VolumeElementForcesAndSourcesCoreOnSide.hpp:259
MoFEM::ForcesAndSourcesCore::UserDataOperator::ptrFE
ForcesAndSourcesCore * ptrFE
Definition: ForcesAndSourcesCore.hpp:985