v0.9.0
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 476 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 487 of file DataOperators.hpp.

492  : cOords_at_GaussPtF3(coords_at_gaussptf3),
493  nOrmals_at_GaussPtF3(normals_at_gaussptf3),
494  tAngent1_at_GaussPtF3(tangent1_at_gaussptf3),
495  tAngent2_at_GaussPtF3(tangent2_at_gaussptf3),
496  cOords_at_GaussPtF4(coords_at_gaussptf4),
497  nOrmals_at_GaussPtF4(normals_at_gaussptf4),
498  tAngent1_at_GaussPtF4(tangent1_at_gaussptf4),
499  tAngent2_at_GaussPtF4(tangent2_at_gaussptf4) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnPrism::calculateNormals ( )

Definition at line 1327 of file DataOperators.cpp.

1327  {
1329 
1330  sPin.resize(3, 3);
1331  sPin.clear();
1332  nOrmals_at_GaussPtF3.resize(tAngent1_at_GaussPtF3.size1(), 3, false);
1333  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF3.size1(); ++gg) {
1334  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF3(gg, 0));
1335  CHKERRG(ierr);
1336  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1337  &tAngent2_at_GaussPtF3(gg, 0), 1, 0.,
1338  &nOrmals_at_GaussPtF3(gg, 0), 1);
1339  }
1340  sPin.clear();
1341  nOrmals_at_GaussPtF4.resize(tAngent1_at_GaussPtF4.size1(), 3, false);
1342  for (unsigned int gg = 0; gg < tAngent1_at_GaussPtF4.size1(); ++gg) {
1343  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPtF4(gg, 0));
1344  CHKERRG(ierr);
1345  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1346  &tAngent2_at_GaussPtF4(gg, 0), 1, 0.,
1347  &nOrmals_at_GaussPtF4(gg, 0), 1);
1348  }
1350 }
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:558
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:477
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:544
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:407

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

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

Member Data Documentation

◆ cOords_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF3

Definition at line 478 of file DataOperators.hpp.

◆ cOords_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::cOords_at_GaussPtF4

Definition at line 482 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF3

Definition at line 479 of file DataOperators.hpp.

◆ nOrmals_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::nOrmals_at_GaussPtF4

Definition at line 483 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnPrism::sPin

Definition at line 501 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF3

Definition at line 480 of file DataOperators.hpp.

◆ tAngent1_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent1_at_GaussPtF4

Definition at line 484 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF3

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF3

Definition at line 481 of file DataOperators.hpp.

◆ tAngent2_at_GaussPtF4

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnPrism::tAngent2_at_GaussPtF4

Definition at line 485 of file DataOperators.hpp.


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