v0.8.13
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 1399 of file DataOperators.cpp.

1399  {
1401 
1402  sPin.resize(3, 3);
1403  sPin.clear();
1404  nOrmals_at_GaussPtF3.resize(tAngent1_at_GaussPtF3.size1(), 3, false);
1405  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF3.size1(); ++gg) {
1406  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF3(gg, 0));
1407  CHKERRG(ierr);
1408  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1409  &tAngent2_at_GaussPtF3(gg, 0), 1, 0.,
1410  &nOrmals_at_GaussPtF3(gg, 0), 1);
1411  }
1412  sPin.clear();
1413  nOrmals_at_GaussPtF4.resize(tAngent1_at_GaussPtF4.size1(), 3, false);
1414  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF4.size1(); ++gg) {
1415  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF4(gg, 0));
1416  CHKERRG(ierr);
1417  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1418  &tAngent2_at_GaussPtF4(gg, 0), 1, 0.,
1419  &nOrmals_at_GaussPtF4(gg, 0), 1);
1420  }
1422 }
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:495
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:562
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:439

◆ 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 1315 of file DataOperators.cpp.

1316  {
1318 
1319  if (data.getFieldData().size() == 0)
1321  const int valid_edges3[] = {1, 1, 1, 0, 0, 0, 0, 0, 0};
1322  const int valid_faces3[] = {0, 0, 0, 1, 0, 0, 0, 0, 0};
1323  const int valid_edges4[] = {0, 0, 0, 0, 0, 0, 1, 1, 1};
1324  const int valid_faces4[] = {0, 0, 0, 0, 1, 0, 0, 0, 0};
1325 
1326  if (type == MBEDGE) {
1327  if (!valid_edges3[side] || valid_edges4[side])
1329  } else if (type == MBTRI) {
1330  if (!valid_faces3[side] || valid_faces4[side])
1332  }
1333 
1334  switch (type) {
1335  case MBVERTEX: {
1336  for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
1337  for (int dd = 0; dd < 3; dd++) {
1338  cOords_at_GaussPtF3(gg, dd) =
1339  cblas_ddot(3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1341  3, &data.getDiffN()(gg, 0), 2, &data.getFieldData()[dd], 3);
1343  3, &data.getDiffN()(gg, 1), 2, &data.getFieldData()[dd], 3);
1345  3, &data.getN(gg)[0], 1, &data.getFieldData()[9 + dd], 3);
1347  3, &data.getDiffN()(gg, 6 + 0), 2, &data.getFieldData()[9 + dd], 3);
1349  3, &data.getDiffN()(gg, 6 + 1), 2, &data.getFieldData()[9 + dd], 3);
1350  }
1351  }
1352  } break;
1353  case MBEDGE:
1354  case MBTRI: {
1355  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1356  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1357  }
1358  unsigned int nb_dofs = data.getFieldData().size();
1359  if (nb_dofs % 3 != 0) {
1360  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1361  }
1362  if (nb_dofs > 3 * data.getN().size2()) {
1363  SETERRQ2(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1364  "data inconsistency, side %d type %d", side, type);
1365  }
1366  for (unsigned int gg = 0; gg < data.getN().size1(); ++gg) {
1367  for (int dd = 0; dd < 3; dd++) {
1368  if ((type == MBTRI && valid_faces3[side]) ||
1369  (type == MBEDGE && valid_edges3[side])) {
1371  nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1372  tAngent1_at_GaussPtF3(gg, dd) +=
1373  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1374  &data.getFieldData()[dd], 3);
1375  tAngent2_at_GaussPtF3(gg, dd) +=
1376  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1377  &data.getFieldData()[dd], 3);
1378  } else if ((type == MBTRI && valid_faces4[side]) ||
1379  (type == MBEDGE && valid_edges4[side])) {
1381  nb_dofs / 3, &data.getN(gg)[0], 1, &data.getFieldData()[dd], 3);
1382  tAngent1_at_GaussPtF4(gg, dd) +=
1383  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1384  &data.getFieldData()[dd], 3);
1385  tAngent2_at_GaussPtF4(gg, dd) +=
1386  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1387  &data.getFieldData()[dd], 3);
1388  }
1389  }
1390  }
1391  } break;
1392  default:
1393  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1394  }
1395 
1397 }
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:495
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
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:439

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: