v0.8.4
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 left 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 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

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 651 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 662 of file DataOperators.hpp.

670  :
671  cOords_at_GaussPtF3(coords_at_gaussptf3),
672  nOrmals_at_GaussPtF3(normals_at_gaussptf3),
673  tAngent1_at_GaussPtF3(tangent1_at_gaussptf3),
674  tAngent2_at_GaussPtF3(tangent2_at_gaussptf3),
675  cOords_at_GaussPtF4(coords_at_gaussptf4),
676  nOrmals_at_GaussPtF4(normals_at_gaussptf4),
677  tAngent1_at_GaussPtF4(tangent1_at_gaussptf4),
678  tAngent2_at_GaussPtF4(tangent2_at_gaussptf4) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnPrism::calculateNormals ( )

Definition at line 1363 of file DataOperators.cpp.

1363  {
1365 
1366  sPin.resize(3, 3);
1367  sPin.clear();
1368  nOrmals_at_GaussPtF3.resize(tAngent1_at_GaussPtF3.size1(), 3, false);
1369  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF3.size1(); ++gg) {
1370  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF3(gg, 0));
1371  CHKERRG(ierr);
1372  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1373  &tAngent2_at_GaussPtF3(gg, 0), 1, 0.,
1374  &nOrmals_at_GaussPtF3(gg, 0), 1);
1375  }
1376  sPin.clear();
1377  nOrmals_at_GaussPtF4.resize(tAngent1_at_GaussPtF4.size1(), 3, false);
1378  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF4.size1(); ++gg) {
1379  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF4(gg, 0));
1380  CHKERRG(ierr);
1381  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1382  &tAngent2_at_GaussPtF4(gg, 0), 1, 0.,
1383  &nOrmals_at_GaussPtF4(gg, 0), 1);
1384  }
1386 }
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:498
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:443

◆ doWork()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnPrism::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 1279 of file DataOperators.cpp.

1280  {
1282 
1283  if (data.getFieldData().size() == 0)
1285  const int valid_edges3[] = {1, 1, 1, 0, 0, 0, 0, 0, 0};
1286  const int valid_faces3[] = {0, 0, 0, 1, 0, 0, 0, 0, 0};
1287  const int valid_edges4[] = {0, 0, 0, 0, 0, 0, 1, 1, 1};
1288  const int valid_faces4[] = {0, 0, 0, 0, 1, 0, 0, 0, 0};
1289 
1290  if (type == MBEDGE) {
1291  if (!valid_edges3[side] || valid_edges4[side])
1293  } else if (type == MBTRI) {
1294  if (!valid_faces3[side] || valid_faces4[side])
1296  }
1297 
1298  switch (type) {
1299  case MBVERTEX: {
1300  for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
1301  for (int dd = 0; dd < 3; dd++) {
1302  cOords_at_GaussPtF3(gg, dd) =
1303  cblas_ddot(3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1305  3, &data.getDiffN()(gg, 0), 2, &data.getFieldData()[dd], 3);
1307  3, &data.getDiffN()(gg, 1), 2, &data.getFieldData()[dd], 3);
1309  3, &data.getN(gg)[0], 1, &data.getFieldData()[9 + dd], 3);
1311  3, &data.getDiffN()(gg, 6 + 0), 2, &data.getFieldData()[9 + dd], 3);
1313  3, &data.getDiffN()(gg, 6 + 1), 2, &data.getFieldData()[9 + dd], 3);
1314  }
1315  }
1316  } break;
1317  case MBEDGE:
1318  case MBTRI: {
1319  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1320  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1321  }
1322  unsigned int nb_dofs = data.getFieldData().size();
1323  if (nb_dofs % 3 != 0) {
1324  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1325  }
1326  if (nb_dofs > 3 * data.getN().size2()) {
1327  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1328  "data inconsistency, side %d type %d", side, type);
1329  }
1330  for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
1331  for (int dd = 0; dd < 3; dd++) {
1332  if ((type == MBTRI && valid_faces3[side]) ||
1333  (type == MBEDGE && valid_edges3[side])) {
1335  nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1336  tAngent1_at_GaussPtF3(gg, dd) +=
1337  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1338  &data.getFieldData()[dd], 3);
1339  tAngent2_at_GaussPtF3(gg, dd) +=
1340  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1341  &data.getFieldData()[dd], 3);
1342  } else if ((type == MBTRI && valid_faces4[side]) ||
1343  (type == MBEDGE && valid_edges4[side])) {
1345  nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1346  tAngent1_at_GaussPtF4(gg, dd) +=
1347  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1348  &data.getFieldData()[dd], 3);
1349  tAngent2_at_GaussPtF4(gg, dd) +=
1350  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1351  &data.getFieldData()[dd], 3);
1352  }
1353  }
1354  }
1355  } break;
1356  default:
1357  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1358  }
1359 
1361 }
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:498
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
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:443

Member Data Documentation

◆ cOords_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF3

Definition at line 653 of file DataOperators.hpp.

◆ cOords_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF4

Definition at line 657 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF3

Definition at line 654 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF4

Definition at line 658 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnPrism::sPin

Definition at line 680 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF3

Definition at line 655 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF4

Definition at line 659 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF3

Definition at line 656 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF4

Definition at line 660 of file DataOperators.hpp.


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