v0.8.13
Public Member Functions | Public Attributes | List of all members
MoFEM::OpSetInvJacH1ForFace Struct Reference

Transform local reference derivatives of shape functions to global derivatives. More...

#include <src/finite_elements/FaceElementForcesAndSourcesCore.hpp>

Inheritance diagram for MoFEM::OpSetInvJacH1ForFace:
[legend]
Collaboration diagram for MoFEM::OpSetInvJacH1ForFace:
[legend]

Public Member Functions

 OpSetInvJacH1ForFace (MatrixDouble &inv_jac)
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on left hand side. More...
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space)
 
 UserDataOperator (const std::string &field_name, const char type)
 
 UserDataOperator (const std::string &row_field_name, const std::string &col_field_name, const char type, const bool symm=true)
 
double getArea ()
 get area of face More...
 
double getMeasure ()
 get measure of element More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
DEPRECATED auto getTensor1Normal ()
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
DEPRECATED auto getTensor1Tangent1 ()
 
auto getFTensor2Tangent1 ()
 get tangentTwo as tensor More...
 
DEPRECATED auto getTensor2Tangent1 ()
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
DEPRECATED auto getTensor1Coords ()
 
MatrixDoublegetGaussPts ()
 get matrix of integration (Gauss) points on Face Element where columns 0,1 are x,y coordinates respectively and column 2 is a value of weight for example getGaussPts()(1,13) returns y coordinate of 13th Gauss point on particular face 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 Gauss pts. More...
 
DEPRECATED auto getTensor1CoordsAtGaussPts ()
 
MatrixDoublegetHoCoordsAtGaussPts ()
 coordinate at Gauss points (if hierarchical approximation of element geometry) More...
 
auto getFTensor1HoCoordsAtGaussPts ()
 get coordinates at Gauss pts (takes in account ho approx. of geometry) More...
 
DEPRECATED auto getTensor1HoCoordsAtGaussPts ()
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
DEPRECATED MatrixDoublegetNormalsAtGaussPt ()
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
DEPRECATED auto getNormalsAtGaussPt (const int gg)
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
DEPRECATED MatrixDoublegetTangent1AtGaussPt ()
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
DEPRECATED MatrixDoublegetTangent2AtGaussPt ()
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
DEPRECATED auto getTensor1NormalsAtGaussPt ()
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
DEPRECATED auto getTensor1Tangent1AtGaussPt ()
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
DEPRECATED auto getTensor1Tangent2AtGaussPt ()
 
DEPRECATED const FaceElementForcesAndSourcesCoregetFaceElementForcesAndSourcesCore ()
 
DEPRECATED const FaceElementForcesAndSourcesCoregetTriFE ()
 
const FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string &fe_name, VolumeElementForcesAndSourcesCoreOnSide &method)
 
- 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...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
DEPRECATED MoFEMErrorCode getPorblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
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)
 
DEPRECATED MoFEMErrorCode getPorblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 
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
 
- 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 right 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 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...
 

Public Attributes

MatrixDoubleinvJac
 
MatrixDouble diffNinvJac
 
- 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
 

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...
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Transform local reference derivatives of shape functions to global derivatives.

Examples:
cell_forces.cpp, and minimal_surface_area.cpp.

Definition at line 530 of file FaceElementForcesAndSourcesCore.hpp.

Constructor & Destructor Documentation

◆ OpSetInvJacH1ForFace()

MoFEM::OpSetInvJacH1ForFace::OpSetInvJacH1ForFace ( MatrixDouble inv_jac)
Deprecated:
Field name do not needed to construct class, change v0.5.17.

Definition at line 541 of file FaceElementForcesAndSourcesCore.hpp.

543  }
ForcesAndSourcesCore::UserDataOperator UserDataOperator
continuous field
Definition: definitions.h:179

Member Function Documentation

◆ doWork()

MoFEMErrorCode MoFEM::OpSetInvJacH1ForFace::doWork ( int  side,
EntityType  type,
DataForcesAndSourcesCore::EntData data 
)
virtual

Operator for linear form, usually to calculate values on left hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 688 of file FaceElementForcesAndSourcesCore.cpp.

689  {
691  //
692 
693  if (getNumeredEntFiniteElementPtr()->getEntType() != MBTRI &&
694  getNumeredEntFiniteElementPtr()->getEntType() != MBQUAD) {
695  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
696  "This operator can be used only with element which is triangle");
697  }
698 
699  for (int b = AINSWORTH_LEGENDRE_BASE; b != USER_BASE; b++) {
700 
701  FieldApproximationBase base = ApproximationBaseArray[b];
702 
703  try {
704 
705  unsigned int nb_dofs = data.getN(base).size2();
706  if (nb_dofs == 0)
708  unsigned int nb_gauss_pts = data.getN(base).size1();
709  diffNinvJac.resize(nb_gauss_pts, 2 * nb_dofs, false);
710 
711  if (type != MBVERTEX) {
712  if (nb_dofs != data.getDiffN(base).size2() / 2) {
713  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
714  "data inconsistency nb_dofs != data.diffN.size2()/2 ( %u != "
715  "%u/2 )",
716  nb_dofs, data.getDiffN(base).size2());
717  }
718  }
719 
720  // std::cerr << type << std::endl;
721  // std::cerr << data.getDiffN() << std::endl;
722  // std::cerr << std::endl;
723  switch (type) {
724  case MBVERTEX:
725  case MBEDGE:
726  case MBTRI: {
727  for (unsigned int gg = 0; gg < nb_gauss_pts; gg++) {
728  for (unsigned int dd = 0; dd < nb_dofs; dd++) {
730  &*invJac.data().begin(), 2,
731  &data.getDiffN(base)(gg, 2 * dd), 1, 0,
732  &diffNinvJac(gg, 2 * dd), 1);
733  }
734  }
735  data.getDiffN(base).data().swap(diffNinvJac.data());
736  } break;
737  default:
738  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
739  }
740 
741  } catch (std::exception &ex) {
742  std::ostringstream ss;
743  ss << "Error in OpSetInvJacH1ForFace "
744  << "throw in method: " << ex.what() << " at line " << __LINE__
745  << " in file " << __FILE__;
746  SETERRQ(PETSC_COMM_SELF, MOFEM_STD_EXCEPTION_THROW, ss.str().c_str());
747  }
748  }
749 
751 }
void cblas_dgemv(const enum CBLAS_ORDER order, const enum CBLAS_TRANSPOSE TransA, const int M, const int N, const double alpha, const double *A, const int lda, const double *X, const int incX, const double beta, double *Y, const int incY)
Definition: cblas_dgemv.c:11
user implemented approximation base
Definition: definitions.h:151
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T *> &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:519
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
FieldApproximationBase
approximation base
Definition: definitions.h:140
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:142

Member Data Documentation

◆ diffNinvJac

MatrixDouble MoFEM::OpSetInvJacH1ForFace::diffNinvJac

Definition at line 545 of file FaceElementForcesAndSourcesCore.hpp.

◆ invJac

MatrixDouble& MoFEM::OpSetInvJacH1ForFace::invJac

Definition at line 532 of file FaceElementForcesAndSourcesCore.hpp.


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