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

1287  {
1289 
1290  nOrmals_at_GaussPt.resize(tAngent1_at_GaussPt.size1(), 3, false);
1291 
1292  auto get_ftensor1 = [](MatrixDouble &m) {
1294  &m(0, 0), &m(0, 1), &m(0, 2));
1295  };
1296  auto t_normal = get_ftensor1(nOrmals_at_GaussPt);
1297  auto t_t1 = get_ftensor1(tAngent1_at_GaussPt);
1298  auto t_t2 = get_ftensor1(tAngent2_at_GaussPt);
1299 
1303 
1304  for (unsigned int gg = 0; gg != tAngent1_at_GaussPt.size1(); ++gg) {
1305  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
1306  ++t_normal;
1307  ++t_t1;
1308  ++t_t2;
1309  }
1310 
1312 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:519
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Common.hpp:212
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use

◆ 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  tAngent1_at_GaussPt.resize(nb_gauss_pts, 3, false);
1198  tAngent2_at_GaussPt.resize(nb_gauss_pts, 3, false);
1199 
1200  auto get_ftensor1 = [](MatrixDouble &m) {
1202  &m(0, 0), &m(0, 1), &m(0, 2));
1203  };
1204  auto t_coords = get_ftensor1(cOords_at_GaussPt);
1205  auto t_t1 = get_ftensor1(tAngent1_at_GaussPt);
1206  auto t_t2 = get_ftensor1(tAngent2_at_GaussPt);
1208  FTensor::Number<0> N0;
1209  FTensor::Number<1> N1;
1210 
1211  switch (type) {
1212  case MBVERTEX: {
1213  cOords_at_GaussPt.clear();
1214  tAngent1_at_GaussPt.clear();
1215  tAngent2_at_GaussPt.clear();
1216  auto t_base = data.getFTensor0N();
1217  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1218  auto t_data = data.getFTensor1FieldData<3>();
1219  auto t_diff_base = data.getFTensor1DiffN<2>();
1220  for (int nn = 0; nn != 3; nn++) {
1221  t_coords(i) += t_base * t_data(i);
1222  t_t1(i) += t_data(i) * t_diff_base(N0);
1223  t_t2(i) += t_data(i) * t_diff_base(N1);
1224  ++t_data;
1225  ++t_base;
1226  ++t_diff_base;
1227  }
1228  ++t_coords;
1229  ++t_t1;
1230  ++t_t2;
1231  }
1232  } break;
1233  case MBEDGE:
1234  case MBTRI: {
1235  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1236  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1237  }
1238  if (nb_dofs % 3 != 0) {
1239  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1240  }
1241  if (nb_dofs > 3 * data.getN().size2()) {
1242  unsigned int nn = 0;
1243  for (; nn != nb_dofs; nn++) {
1244  if (!data.getFieldDofs()[nn]->getActive())
1245  break;
1246  }
1247  if (nn > 3 * data.getN().size2()) {
1248  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1249  "data inconsistency for base %s",
1250  ApproximationBaseNames[data.getBase()]);
1251  } else {
1252  nb_dofs = nn;
1253  if (!nb_dofs)
1255  }
1256  }
1257  const int nb_base_functions = data.getN().size2();
1258  auto t_base = data.getFTensor0N();
1259  auto t_diff_base = data.getFTensor1DiffN<2>();
1260  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1261  auto t_data = data.getFTensor1FieldData<3>();
1262  int bb = 0;
1263  for (; bb != nb_dofs / 3; ++bb) {
1264  t_coords(i) += t_base * t_data(i);
1265  t_t1(i) += t_data(i) * t_diff_base(N0);
1266  t_t2(i) += t_data(i) * t_diff_base(N1);
1267  ++t_data;
1268  ++t_base;
1269  ++t_diff_base;
1270  }
1271  for (; bb != nb_base_functions;++bb) {
1272  ++t_base;
1273  ++t_diff_base;
1274  }
1275  ++t_coords;
1276  ++t_t1;
1277  ++t_t2;
1278  }
1279  } break;
1280  default:
1281  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1282  }
1283 
1285 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:519
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:526
static const char *const ApproximationBaseNames[]
Definition: definitions.h:155
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Common.hpp:212

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: