v0.8.4
Public Member Functions | Public Attributes | List of all members
MoFEM::OpGetCoordsAndNormalsOnFace Struct Reference

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

#include <src/finite_elements/DataOperators.hpp>

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

Public Member Functions

 OpGetCoordsAndNormalsOnFace (MatrixDouble &coords_at_gausspt, MatrixDouble &normals_at_gausspt, MatrixDouble &tangent1_at_gausspt, MatrixDouble &tangent2_at_gausspt)
 
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_GaussPt
 
MatrixDoublenOrmals_at_GaussPt
 
MatrixDoubletAngent1_at_GaussPt
 
MatrixDoubletAngent2_at_GaussPt
 
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 620 of file DataOperators.hpp.

Constructor & Destructor Documentation

◆ OpGetCoordsAndNormalsOnFace()

MoFEM::OpGetCoordsAndNormalsOnFace::OpGetCoordsAndNormalsOnFace ( MatrixDouble coords_at_gausspt,
MatrixDouble normals_at_gausspt,
MatrixDouble tangent1_at_gausspt,
MatrixDouble tangent2_at_gausspt 
)

Definition at line 627 of file DataOperators.hpp.

631  :
632  cOords_at_GaussPt(coords_at_gausspt),
633  nOrmals_at_GaussPt(normals_at_gausspt),
634  tAngent1_at_GaussPt(tangent1_at_gausspt),
635  tAngent2_at_GaussPt(tangent2_at_gausspt) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnFace::calculateNormals ( )

Definition at line 1261 of file DataOperators.cpp.

1261  {
1263 
1264  sPin.resize(3, 3);
1265  sPin.clear();
1266  nOrmals_at_GaussPt.resize(tAngent1_at_GaussPt.size1(), 3, false);
1267  for (unsigned int gg = 0; gg < tAngent1_at_GaussPt.size1(); ++gg) {
1268  ierr = Spin(&*sPin.data().begin(), &tAngent1_at_GaussPt(gg, 0));
1269  CHKERRG(ierr);
1270  cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., &*sPin.data().begin(), 3,
1271  &tAngent2_at_GaussPt(gg, 0), 1, 0., &nOrmals_at_GaussPt(gg, 0),
1272  1);
1273  }
1274 
1276 }
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 MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:522
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:565
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:528
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Common.hpp:80

◆ doWork()

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

1188  {
1190 
1191  unsigned int nb_dofs = data.getFieldData().size();
1192  if (nb_dofs == 0)
1194 
1195  int nb_gauss_pts = data.getN().size1();
1196  cOords_at_GaussPt.resize(nb_gauss_pts, 3, false);
1197  nOrmals_at_GaussPt.resize(nb_gauss_pts, 3, false);
1198  tAngent1_at_GaussPt.resize(nb_gauss_pts, 3, false);
1199  tAngent2_at_GaussPt.resize(nb_gauss_pts, 3, false);
1200 
1201  // std::cerr << type << " " << side << " " <<
1202  // ApproximationBaseNames[data.getBase()] << std::endl;
1203 
1204  switch (type) {
1205  case MBVERTEX: {
1206  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1207  for (int nn = 0; nn < 3; nn++) {
1208  cOords_at_GaussPt(gg, nn) =
1209  cblas_ddot(3, &data.getN(gg)[0], 1, &data.getFieldData()[nn], 3);
1210  tAngent1_at_GaussPt(gg, nn) = cblas_ddot(3, &data.getDiffN()(0, 0), 2,
1211  &data.getFieldData()[nn], 3);
1212  tAngent2_at_GaussPt(gg, nn) = cblas_ddot(3, &data.getDiffN()(0, 1), 2,
1213  &data.getFieldData()[nn], 3);
1214  }
1215  }
1216  } break;
1217  case MBEDGE:
1218  case MBTRI: {
1219  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1220  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1221  }
1222  if (nb_dofs % 3 != 0) {
1223  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1224  }
1225  if (nb_dofs > 3 * data.getN().size2()) {
1226  unsigned int nn = 0;
1227  for (; nn != nb_dofs; nn++) {
1228  if (!data.getFieldDofs()[nn]->getActive())
1229  break;
1230  }
1231  if (nn > 3 * data.getN().size2()) {
1232  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1233  "data inconsistency for base %s",
1234  ApproximationBaseNames[data.getBase()]);
1235  } else {
1236  nb_dofs = nn;
1237  if (!nb_dofs)
1239  }
1240  }
1241  for (int gg = 0; gg < nb_gauss_pts; ++gg) {
1242  for (int dd = 0; dd < 3; dd++) {
1243  cOords_at_GaussPt(gg, dd) += cblas_ddot(nb_dofs / 3, &data.getN(gg)[0],
1244  1, &data.getFieldData()[dd], 3);
1245  tAngent1_at_GaussPt(gg, dd) +=
1246  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 0), 2,
1247  &data.getFieldData()[dd], 3);
1248  tAngent2_at_GaussPt(gg, dd) +=
1249  cblas_ddot(nb_dofs / 3, &data.getDiffN()(gg, 1), 2,
1250  &data.getFieldData()[dd], 3);
1251  }
1252  }
1253  } break;
1254  default:
1255  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1256  }
1257 
1259 }
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:522
#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
static const char *const ApproximationBaseNames[]
Definition: definitions.h:153

Member Data Documentation

◆ cOords_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::cOords_at_GaussPt

Definition at line 622 of file DataOperators.hpp.

◆ nOrmals_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::nOrmals_at_GaussPt

Definition at line 623 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnFace::sPin

Definition at line 637 of file DataOperators.hpp.

◆ tAngent1_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::tAngent1_at_GaussPt

Definition at line 624 of file DataOperators.hpp.

◆ tAngent2_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::tAngent2_at_GaussPt

Definition at line 625 of file DataOperators.hpp.


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