v0.8.20
Public Member Functions | Public Attributes | List of all members
MoFEM::OpGetCoordsAndNormalsOnPrism Struct Reference

calculate normals at Gauss points of triangle element More...

#include <src/finite_elements/DataOperators.hpp>

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

Public Member Functions

 OpGetCoordsAndNormalsOnPrism (MatrixDouble &coords_at_gaussptf3, MatrixDouble &normals_at_gaussptf3, MatrixDouble &tangent1_at_gaussptf3, MatrixDouble &tangent2_at_gaussptf3, MatrixDouble &coords_at_gaussptf4, MatrixDouble &normals_at_gaussptf4, MatrixDouble &tangent1_at_gaussptf4, MatrixDouble &tangent2_at_gaussptf4)
 
MoFEMErrorCode doWork (int side, EntityType type, DataForcesAndSourcesCore::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
MoFEMErrorCode calculateNormals ()
 
- 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 left 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

MatrixDoublecOords_at_GaussPtF3
 
MatrixDoublenOrmals_at_GaussPtF3
 
MatrixDoubletAngent1_at_GaussPtF3
 
MatrixDoubletAngent2_at_GaussPtF3
 
MatrixDoublecOords_at_GaussPtF4
 
MatrixDoublenOrmals_at_GaussPtF4
 
MatrixDoubletAngent1_at_GaussPtF4
 
MatrixDoubletAngent2_at_GaussPtF4
 
MatrixDouble sPin
 
- 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
 

Detailed Description

calculate normals at Gauss points of triangle element

Definition at line 628 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ OpGetCoordsAndNormalsOnPrism()

MoFEM::OpGetCoordsAndNormalsOnPrism::OpGetCoordsAndNormalsOnPrism ( MatrixDouble coords_at_gaussptf3,
MatrixDouble normals_at_gaussptf3,
MatrixDouble tangent1_at_gaussptf3,
MatrixDouble tangent2_at_gaussptf3,
MatrixDouble coords_at_gaussptf4,
MatrixDouble normals_at_gaussptf4,
MatrixDouble tangent1_at_gaussptf4,
MatrixDouble tangent2_at_gaussptf4 
)

Definition at line 639 of file DataOperators.hpp.

644  : cOords_at_GaussPtF3(coords_at_gaussptf3),
645  nOrmals_at_GaussPtF3(normals_at_gaussptf3),
646  tAngent1_at_GaussPtF3(tangent1_at_gaussptf3),
647  tAngent2_at_GaussPtF3(tangent2_at_gaussptf3),
648  cOords_at_GaussPtF4(coords_at_gaussptf4),
649  nOrmals_at_GaussPtF4(normals_at_gaussptf4),
650  tAngent1_at_GaussPtF4(tangent1_at_gaussptf4),
651  tAngent2_at_GaussPtF4(tangent2_at_gaussptf4) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnPrism::calculateNormals ( )

Definition at line 1335 of file DataOperators.cpp.

1335  {
1337 
1338  sPin.resize(3, 3);
1339  sPin.clear();
1340  nOrmals_at_GaussPtF3.resize(tAngent1_at_GaussPtF3.size1(), 3, false);
1341  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF3.size1(); ++gg) {
1342  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF3(gg, 0));
1343  CHKERRG(ierr);
1344  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1345  &tAngent2_at_GaussPtF3(gg, 0), 1, 0.,
1346  &nOrmals_at_GaussPtF3(gg, 0), 1);
1347  }
1348  sPin.clear();
1349  nOrmals_at_GaussPtF4.resize(tAngent1_at_GaussPtF4.size1(), 3, false);
1350  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF4.size1(); ++gg) {
1351  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF4(gg, 0));
1352  CHKERRG(ierr);
1353  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1354  &tAngent2_at_GaussPtF4(gg, 0), 1, 0.,
1355  &nOrmals_at_GaussPtF4(gg, 0), 1);
1356  }
1358 }
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
PetscErrorCode Spin(double *spinOmega, double *vecOmega)
calculate spin matrix from vector
Definition: fem_tools.c:480
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:542
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

◆ doWork()

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

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

Reimplemented from MoFEM::DataOperator.

Definition at line 1251 of file DataOperators.cpp.

1252  {
1254 
1255  if (data.getFieldData().size() == 0)
1257  const int valid_edges3[] = {1, 1, 1, 0, 0, 0, 0, 0, 0};
1258  const int valid_faces3[] = {0, 0, 0, 1, 0, 0, 0, 0, 0};
1259  const int valid_edges4[] = {0, 0, 0, 0, 0, 0, 1, 1, 1};
1260  const int valid_faces4[] = {0, 0, 0, 0, 1, 0, 0, 0, 0};
1261 
1262  if (type == MBEDGE) {
1263  if (!valid_edges3[side] || valid_edges4[side])
1265  } else if (type == MBTRI) {
1266  if (!valid_faces3[side] || valid_faces4[side])
1268  }
1269 
1270  switch (type) {
1271  case MBVERTEX: {
1272  for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
1273  for (int dd = 0; dd < 3; dd++) {
1274  cOords_at_GaussPtF3(gg, dd) =
1275  cblas_ddot(3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1277  3, &data.getDiffN()(gg, 0), 2, &data.getFieldData()[dd], 3);
1279  3, &data.getDiffN()(gg, 1), 2, &data.getFieldData()[dd], 3);
1281  3, &data.getN(gg)[0], 1, &data.getFieldData()[9 + dd], 3);
1283  3, &data.getDiffN()(gg, 6 + 0), 2, &data.getFieldData()[9 + dd], 3);
1285  3, &data.getDiffN()(gg, 6 + 1), 2, &data.getFieldData()[9 + dd], 3);
1286  }
1287  }
1288  } break;
1289  case MBEDGE:
1290  case MBTRI: {
1291  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1292  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1293  }
1294  unsigned int nb_dofs = data.getFieldData().size();
1295  if (nb_dofs % 3 != 0) {
1296  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1297  }
1298  if (nb_dofs > 3 * data.getN().size2()) {
1299  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1300  "data inconsistency, side %d type %d", side, type);
1301  }
1302  for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
1303  for (int dd = 0; dd < 3; dd++) {
1304  if ((type == MBTRI && valid_faces3[side]) ||
1305  (type == MBEDGE && valid_edges3[side])) {
1307  nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1308  tAngent1_at_GaussPtF3(gg, dd) +=
1309  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1310  &data.getFieldData()[dd], 3);
1311  tAngent2_at_GaussPtF3(gg, dd) +=
1312  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1313  &data.getFieldData()[dd], 3);
1314  } else if ((type == MBTRI && valid_faces4[side]) ||
1315  (type == MBEDGE && valid_edges4[side])) {
1317  nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1318  tAngent1_at_GaussPtF4(gg, dd) +=
1319  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1320  &data.getFieldData()[dd], 3);
1321  tAngent2_at_GaussPtF4(gg, dd) +=
1322  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1323  &data.getFieldData()[dd], 3);
1324  }
1325  }
1326  }
1327  } break;
1328  default:
1329  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1330  }
1331 
1333 }
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 MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:475
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
double cblas_ddot(const int N, const double *X, const int incX, const double *Y, const int incY)
Definition: cblas_ddot.c:12
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:405

Member Data Documentation

◆ cOords_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF3

Definition at line 630 of file DataOperators.hpp.

◆ cOords_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF4

Definition at line 634 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF3

Definition at line 631 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF4

Definition at line 635 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnPrism::sPin

Definition at line 653 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF3

Definition at line 632 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF4

Definition at line 636 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF3

Definition at line 633 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF4

Definition at line 637 of file DataOperators.hpp.


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