v0.9.1
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...
 
MatrixDoublegetNormalsAtGaussPtF3 ()
 if higher order geometry return normals at face F3 at Gauss pts. More...
 
MatrixDoublegetNormalsAtGaussPtF4 ()
 if higher order geometry return normals at face F4 at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPtF3 (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPtF4 (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=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)
 
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...
 
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...
 
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
 
DEPRECATED Vec getSnesF () const
 
DEPRECATED Vec getSnesX () const
 
DEPRECATED Mat getSnesA () const
 
DEPRECATED 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 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
 
- 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, 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)
 
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 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)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 

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...
 
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...
 
- 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 46 of file FlatPrismElementForcesAndSourcesCore.hpp.

Member Function Documentation

◆ getArea()

double UserDataOperator::getArea ( const int  dd)

get face aRea

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

Definition at line 187 of file FlatPrismElementForcesAndSourcesCore.hpp.

187  {
188  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
189 }

◆ getAreaF3()

double UserDataOperator::getAreaF3 ( )

Definition at line 192 of file FlatPrismElementForcesAndSourcesCore.hpp.

192  {
193  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[0];
194 }

◆ getAreaF4()

double UserDataOperator::getAreaF4 ( )

Definition at line 196 of file FlatPrismElementForcesAndSourcesCore.hpp.

196  {
197  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->aRea[1];
198 }

◆ getCoords()

VectorDouble & UserDataOperator::getCoords ( )

get triangle coordinates

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

Definition at line 220 of file FlatPrismElementForcesAndSourcesCore.hpp.

220  {
221  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->coords;
222 }

◆ getCoordsAtGaussPts()

MatrixDouble & UserDataOperator::getCoordsAtGaussPts ( )

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 225 of file FlatPrismElementForcesAndSourcesCore.hpp.

225  {
226  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
227  ->coordsAtGaussPts;
228 }

◆ getFlatPrismElementForcesAndSourcesCore()

const FlatPrismElementForcesAndSourcesCore * UserDataOperator::getFlatPrismElementForcesAndSourcesCore ( )

return pointer to triangle finite element object

Definition at line 298 of file FlatPrismElementForcesAndSourcesCore.hpp.

298  {
299  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE);
300 }

◆ getHoCoordsAtGaussPtsF3()

MatrixDouble & UserDataOperator::getHoCoordsAtGaussPtsF3 ( )

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

Definition at line 231 of file FlatPrismElementForcesAndSourcesCore.hpp.

231  {
232  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
233  ->hoCoordsAtGaussPtsF3;
234 }

◆ getHoCoordsAtGaussPtsF4()

MatrixDouble & UserDataOperator::getHoCoordsAtGaussPtsF4 ( )

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

Definition at line 237 of file FlatPrismElementForcesAndSourcesCore.hpp.

237  {
238  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
239  ->hoCoordsAtGaussPtsF4;
240 }

◆ getNormal()

VectorDouble & UserDataOperator::getNormal ( )

get triangle normal

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

Definition at line 201 of file FlatPrismElementForcesAndSourcesCore.hpp.

201  {
202  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->normal;
203 }

◆ getNormalF3()

VectorAdaptor UserDataOperator::getNormalF3 ( )

Definition at line 206 of file FlatPrismElementForcesAndSourcesCore.hpp.

206  {
207  double *data =
208  &(static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[0]);
209  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
210 }
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:104

◆ getNormalF4()

VectorAdaptor UserDataOperator::getNormalF4 ( )

Definition at line 213 of file FlatPrismElementForcesAndSourcesCore.hpp.

213  {
214  double *data =
215  &(static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)->normal[3]);
216  return VectorAdaptor(3, ublas::shallow_array_adaptor<double>(3, data));
217 }
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:104

◆ getNormalsAtGaussPtF3() [1/2]

MatrixDouble & UserDataOperator::getNormalsAtGaussPtF3 ( )

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

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

Definition at line 243 of file FlatPrismElementForcesAndSourcesCore.hpp.

243  {
244  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
245  ->nOrmals_at_GaussPtF3;
246 }

◆ getNormalsAtGaussPtF3() [2/2]

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

if higher order geometry return normals at Gauss pts.

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

Parameters
gggauss point number

Definition at line 255 of file FlatPrismElementForcesAndSourcesCore.hpp.

256  {
257  return ublas::matrix_row<MatrixDouble>(
258  static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
260  gg);
261 }

◆ getNormalsAtGaussPtF4() [1/2]

MatrixDouble & UserDataOperator::getNormalsAtGaussPtF4 ( )

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

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

Definition at line 249 of file FlatPrismElementForcesAndSourcesCore.hpp.

249  {
250  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
251  ->nOrmals_at_GaussPtF4;
252 }

◆ getNormalsAtGaussPtF4() [2/2]

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

if higher order geometry return normals at Gauss pts.

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

Parameters
gggauss point number

Definition at line 264 of file FlatPrismElementForcesAndSourcesCore.hpp.

265  {
266  return ublas::matrix_row<MatrixDouble>(
267  static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
269  gg);
270 }

◆ getTangent1AtGaussPtF3()

MatrixDouble & UserDataOperator::getTangent1AtGaussPtF3 ( )

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

Definition at line 273 of file FlatPrismElementForcesAndSourcesCore.hpp.

273  {
274  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
275  ->tAngent1_at_GaussPtF3;
276 }

◆ getTangent1AtGaussPtF4()

MatrixDouble & UserDataOperator::getTangent1AtGaussPtF4 ( )

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

Definition at line 285 of file FlatPrismElementForcesAndSourcesCore.hpp.

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

◆ getTangent2AtGaussPtF3()

MatrixDouble & UserDataOperator::getTangent2AtGaussPtF3 ( )

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

Definition at line 279 of file FlatPrismElementForcesAndSourcesCore.hpp.

279  {
280  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
281  ->tAngent2_at_GaussPtF3;
282 }

◆ getTangent2AtGaussPtF4()

MatrixDouble & UserDataOperator::getTangent2AtGaussPtF4 ( )

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

Definition at line 291 of file FlatPrismElementForcesAndSourcesCore.hpp.

291  {
292  return static_cast<FlatPrismElementForcesAndSourcesCore *>(ptrFE)
293  ->tAngent2_at_GaussPtF4;
294 }

◆ setPtrFE()

MoFEMErrorCode UserDataOperator::setPtrFE ( ForcesAndSourcesCore ptr)
protectedvirtual

Reimplemented from MoFEM::ForcesAndSourcesCore::UserDataOperator.

Definition at line 230 of file FlatPrismElementForcesAndSourcesCore.cpp.

231  {
233  if (!(ptrFE = dynamic_cast<FlatPrismElementForcesAndSourcesCore *>(ptr)))
234  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
235  "User operator and finite element do not work together");
237 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:506
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:513

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