v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | List of all members
MoFEM::FlatPrismElementForcesAndSourcesCore::UserDataOperator Struct Reference

default operator for Flat Prism element More...

#include <src/finite_elements/FlatPrismElementForcesAndSourcesCore.hpp>

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

Public Member Functions

double getArea (const int dd)
 get face aRea More...
 
double getAreaF3 ()
 
double getAreaF4 ()
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorAdaptor getNormalF3 ()
 
VectorAdaptor getNormalF4 ()
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
MatrixDoublegetCoordsAtGaussPts ()
 get coordinates at Gauss pts. More...
 
MatrixDoublegetHOCoordsAtGaussPtsF3 ()
 coordinate at Gauss points on face 3 (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetHOCoordsAtGaussPtsF4 ()
 coordinate at Gauss points on face 4 (if hierarchical approximation of element geometry) More...
 
MatrixDoublegetNormalsAtGaussPtsF3 ()
 if higher order geometry return normals at face F3 at Gauss pts. More...
 
MatrixDoublegetNormalsAtGaussPtsF4 ()
 if higher order geometry return normals at face F4 at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPtsF3 (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPtsF4 (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPtF3 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPtF3 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPtF4 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPtF4 ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
const FlatPrismElementForcesAndSourcesCoregetFlatPrismElementForcesAndSourcesCore ()
 return pointer to triangle 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)
 
virtual 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 Flat Prism element

Examples
forces_and_sources_testing_flat_prism_element.cpp.

Definition at line 34 of file FlatPrismElementForcesAndSourcesCore.hpp.

Member Function Documentation

◆ getArea()

double FlatPrismElementForcesAndSourcesCore::UserDataOperator::getArea ( const int  dd)
inline

get face aRea

Parameters
ddif dd == 0 it is for face F3 if dd == 1 is for face F4

Definition at line 175 of file FlatPrismElementForcesAndSourcesCore.hpp.

175 {
176 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
177}

◆ getAreaF3()

double FlatPrismElementForcesAndSourcesCore::UserDataOperator::getAreaF3 ( )
inline

Definition at line 180 of file FlatPrismElementForcesAndSourcesCore.hpp.

180 {
181 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
182}

◆ getAreaF4()

double FlatPrismElementForcesAndSourcesCore::UserDataOperator::getAreaF4 ( )
inline

Definition at line 184 of file FlatPrismElementForcesAndSourcesCore.hpp.

184 {
185 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[1];
186}

◆ getCoords()

VectorDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getCoords ( )
inline

get triangle coordinates

Vector has 6 elements, i.e. coordinates on face F3 and F4

Definition at line 208 of file FlatPrismElementForcesAndSourcesCore.hpp.

208 {
209 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->coords;
210}

◆ getCoordsAtGaussPts()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts ( )
inline

get coordinates at Gauss pts.

Matrix has size (nb integration points)x(coordinates on F3 and F4 = 6), i.e. coordinates on face F3 and F4

Definition at line 213 of file FlatPrismElementForcesAndSourcesCore.hpp.

213 {
214 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
216}

◆ getFlatPrismElementForcesAndSourcesCore()

const FlatPrismElementForcesAndSourcesCore * FlatPrismElementForcesAndSourcesCore::UserDataOperator::getFlatPrismElementForcesAndSourcesCore ( )
inline

return pointer to triangle finite element object

Definition at line 285 of file FlatPrismElementForcesAndSourcesCore.hpp.

286 {
287 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE);
288}

◆ getHOCoordsAtGaussPtsF3()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getHOCoordsAtGaussPtsF3 ( )
inline

coordinate at Gauss points on face 3 (if hierarchical approximation of element geometry)

Definition at line 218 of file FlatPrismElementForcesAndSourcesCore.hpp.

219 {
220 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
222}

◆ getHOCoordsAtGaussPtsF4()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getHOCoordsAtGaussPtsF4 ( )
inline

coordinate at Gauss points on face 4 (if hierarchical approximation of element geometry)

Definition at line 224 of file FlatPrismElementForcesAndSourcesCore.hpp.

225 {
226 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
228}

◆ getNormal()

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

get triangle normal

Normal has 6 elements, first 3 are for face F3 another three for face F4

Definition at line 189 of file FlatPrismElementForcesAndSourcesCore.hpp.

189 {
190 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->normal;
191}

◆ getNormalF3()

VectorAdaptor FlatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalF3 ( )
inline

Definition at line 194 of file FlatPrismElementForcesAndSourcesCore.hpp.

194 {
195 double *data =
196 &(static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[0]);
197 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
198}
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115

◆ getNormalF4()

VectorAdaptor FlatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalF4 ( )
inline

Definition at line 201 of file FlatPrismElementForcesAndSourcesCore.hpp.

201 {
202 double *data =
203 &(static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[3]);
204 return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
205}

◆ getNormalsAtGaussPtsF3() [1/2]

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF3 ( )
inline

if higher order geometry return normals at face F3 at Gauss pts.

Face 3 is top face in canonical triangle numeration, see [48]

Definition at line 230 of file FlatPrismElementForcesAndSourcesCore.hpp.

231 {
232 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
234}

◆ getNormalsAtGaussPtsF3() [2/2]

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

if higher order geometry return normals at Gauss pts.

Face 3 is top face in canonical triangle numeration, see [48]

Parameters
gggauss point number

Definition at line 243 of file FlatPrismElementForcesAndSourcesCore.hpp.

244 {
245 return ublas::matrix_row<MatrixDouble>(
246 static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
248 gg);
249}

◆ getNormalsAtGaussPtsF4() [1/2]

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getNormalsAtGaussPtsF4 ( )
inline

if higher order geometry return normals at face F4 at Gauss pts.

Face 4 is top face in canonical triangle numeration, see [48]

Definition at line 236 of file FlatPrismElementForcesAndSourcesCore.hpp.

237 {
238 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
240}

◆ getNormalsAtGaussPtsF4() [2/2]

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

if higher order geometry return normals at Gauss pts.

Face 3 is top face in canonical triangle numeration, see [48]

Parameters
gggauss point number

Definition at line 252 of file FlatPrismElementForcesAndSourcesCore.hpp.

253 {
254 return ublas::matrix_row<MatrixDouble>(
255 static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
257 gg);
258}

◆ getTangent1AtGaussPtF3()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent1AtGaussPtF3 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 260 of file FlatPrismElementForcesAndSourcesCore.hpp.

261 {
262 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
264}

◆ getTangent1AtGaussPtF4()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent1AtGaussPtF4 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 272 of file FlatPrismElementForcesAndSourcesCore.hpp.

273 {
274 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
276}

◆ getTangent2AtGaussPtF3()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent2AtGaussPtF3 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 266 of file FlatPrismElementForcesAndSourcesCore.hpp.

267 {
268 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
270}

◆ getTangent2AtGaussPtF4()

MatrixDouble & FlatPrismElementForcesAndSourcesCore::UserDataOperator::getTangent2AtGaussPtF4 ( )
inline

if higher order geometry return tangent vector to triangle at Gauss pts.

Definition at line 278 of file FlatPrismElementForcesAndSourcesCore.hpp.

279 {
280 return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
282}

◆ setPtrFE()

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

Reimplemented from MoFEM::ForcesAndSourcesCore::UserDataOperator.

Definition at line 237 of file FlatPrismElementForcesAndSourcesCore.cpp.

238 {
240 if (!(ptrFE = dynamic_cast<FlatPrismElementForcesAndSourcesCore *>(ptr)))
241 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
242 "User operator and finite element do not work together");
244}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
@ 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 ...
Definition: definitions.h:440

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