v0.15.0
Loading...
Searching...
No Matches
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
 Access the owning volume finite element.
 
FaceElementForcesAndSourcesCoregetFaceFE () const
 Access the face finite element Face element, execute side element, that is volume, on which the volume side operators are executed.
 
DEPRECATED int getFaceSense () const
 Get face sense with respect to the volume.
 
int getSkeletonSense () const
 Get the skeleton sense.
 
int getFaceSideNumber () const
 Get face side number with respect to the volume.
 
bool getEdgeFace (const int ee) const
 Query if an edge belongs to the current face.
 
VectorDoublegetNormal ()
 Get the face normal for this side element.
 
VectorDoublegetTangent1 ()
 Get first triangle tangent.
 
VectorDoublegetTangent2 ()
 Get second triangle tangent.
 
auto getFTensor1Normal ()
 Get normal as tensor.
 
auto getFTensor1Tangent1 ()
 Get tangent one as tensor.
 
auto getFTensor1Tangent2 ()
 Get tangent two as tensor.
 
MatrixDoublegetNormalsAtGaussPts ()
 If higher-order geometry is available, return normals at Gauss points.
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 If higher-order geometry is available, return normals at Gauss points.
 
auto getFTensor1NormalsAtGaussPts ()
 Get normal at integration points.
 
MatrixDoublegetFaceCoordsAtGaussPts ()
 Get face coordinates at Gauss points.
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes
 
const EntityHandlegetConn ()
 get element connectivity
 
double getVolume () const
 element volume (linear geometry)
 
doublegetVolume ()
 element volume (linear geometry)
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian
 
VectorDoublegetCoords ()
 nodal coordinates
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 Constructor for operators working on finite element spaces.
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 Constructor for operators working on a single field.
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 Constructor for operators working on two fields (bilinear forms)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
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
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, 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.
 
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.
 
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.
 
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.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- 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.
 
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.
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Protected Member Functions

MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 Assign base class pointer for the operator.
 

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.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 
- 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
mofem/atom_tests/continuity_check_on_skeleton_3d.cpp, mofem/tutorials/fun-0/hello_world.cpp, and mofem/users_modules/eshelbian_plasticity/src/impl/EshelbianPlasticity.cpp.

Definition at line 146 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

Member Function Documentation

◆ getEdgeFace()

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

Query if an edge belongs to the current face.

Parameters
eeLocal edge index (0-5 for tetrahedra)
Returns
True when the edge lies on the active face

Definition at line 370 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

371 {
372 constexpr bool edges_on_faces[6][4] = {{true, false, false, true}, // e0
373 {false, true, false, true}, // e1
374 {false, false, true, true}, // e2
375 {true, false, true, false}, // e3
376 {true, true, false, false}, // e4
377 {false, true, true, false}};
378 return edges_on_faces[ee][getFaceSideNumber()];
379}
int getFaceSideNumber() const
Get face side number with respect to the volume.

◆ getFaceCoordsAtGaussPts()

MatrixDouble & VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceCoordsAtGaussPts ( )
inline

Get face coordinates at Gauss points.

Note
Coordinates should match the values returned by getCoordsAtGaussPts on tetrahedra. If they differ, it indicates an error or an unusual setup.

Definition at line 365 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

366 {
367 return getFaceFE()->coordsAtGaussPts;
368}
MatrixDouble coordsAtGaussPts
coordinated at gauss points
FaceElementForcesAndSourcesCore * getFaceFE() const
Access the face finite element Face element, execute side element, that is volume,...

◆ getFaceFE()

FaceElementForcesAndSourcesCore * VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceFE ( ) const
inline

Access the face finite element Face element, execute side element, that is volume, on which the volume side operators are executed.

Returns
Pointer to the face finite element attached to this side

Definition at line 295 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

295 {
296 return static_cast<FaceElementForcesAndSourcesCore *>(
298}
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
VolumeElementForcesAndSourcesCoreOnSide * getVolumeFE() const
Access the owning volume finite element.

◆ getFaceSense()

int VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceSense ( ) const
inline

Get face sense with respect to the volume.

Deprecated:
Use getSkeletonSense()
Returns
Face sense (+1 or -1)

Definition at line 320 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

◆ getFaceSideNumber()

int VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFaceSideNumber ( ) const
inline

Get face side number with respect to the volume.

Returns
Face side number

Definition at line 355 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

◆ getFTensor1Normal()

auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFTensor1Normal ( )
inline

Get normal as tensor.

Definition at line 330 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

331 {
332 double *ptr = &*getNormal().data().begin();
333 return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
334}

◆ getFTensor1NormalsAtGaussPts()

auto VolumeElementForcesAndSourcesCoreOnSide::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;
}
FTensor::Index< 'i', SPACE_DIM > i

Definition at line 348 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

349 {
350 double *ptr = &*getNormalsAtGaussPts().data().begin();
352 &ptr[2]);
353}
MatrixDouble & getNormalsAtGaussPts()
If higher-order geometry is available, return normals at Gauss points.

◆ getFTensor1Tangent1()

auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFTensor1Tangent1 ( )
inline

Get tangent one as tensor.

Definition at line 336 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

337 {
338 double *ptr = &*getTangent1().data().begin();
339 return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
340}

◆ getFTensor1Tangent2()

auto VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getFTensor1Tangent2 ( )
inline

Get tangent two as tensor.

Definition at line 342 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

343 {
344 double *ptr = &*getTangent2().data().begin();
345 return FTensor::Tensor1<double *, 3>(ptr, &ptr[1], &ptr[2]);
346}

◆ getNormal()

VectorDouble & VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getNormal ( )
inline

Get the face normal for this side element.

Returns
Face normal vector

Definition at line 301 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

◆ getNormalsAtGaussPts() [1/2]

MatrixDouble & VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getNormalsAtGaussPts ( )
inline

If higher-order geometry is available, return normals at Gauss points.

Note
The returned matrix has zero rows and columns when no higher-order approximation of geometry is available.

Definition at line 360 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

◆ getNormalsAtGaussPts() [2/2]

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

If higher-order geometry is available, return normals at Gauss points.

Parameters
ggGauss point number

Definition at line 382 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

383 {
384 return ublas::matrix_row<MatrixDouble>(getNormalsAtGaussPts(), gg);
385}

◆ getSkeletonSense()

int VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getSkeletonSense ( ) const
inline

Get the skeleton sense.

Calls getFaceSense().

Returns
Skeleton sense (+1 or -1)

Definition at line 325 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

326 {
327 return getVolumeFE()->faceSense;
328}

◆ getTangent1()

VectorDouble & VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getTangent1 ( )
inline

Get first triangle tangent.

Definition at line 306 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

◆ getTangent2()

VectorDouble & VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getTangent2 ( )
inline

Get second triangle tangent.

Definition at line 311 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

◆ getVolumeFE()

VolumeElementForcesAndSourcesCoreOnSide * VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::getVolumeFE ( ) const
inline

Access the owning volume finite element.

Volume element is side element, to face element.

Returns
Pointer to the parent volume finite element

Definition at line 316 of file VolumeElementForcesAndSourcesCoreOnSide.hpp.

316 {
317 return static_cast<VolumeElementForcesAndSourcesCoreOnSide *>(ptrFE);
318}

◆ setPtrFE()

MoFEMErrorCode VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Assign base class pointer for the operator.

Parameters
ptrPointer to the forces and sources core instance
Returns
Error code (0 on success)

Reimplemented from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator.

Definition at line 70 of file VolumeElementForcesAndSourcesCoreOnSide.cpp.

71 {
73 if (!(ptrFE = dynamic_cast<VolumeElementForcesAndSourcesCoreOnSide *>(ptr)))
74 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
75 "User operator and finite element do not work together");
77}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

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