v0.8.19
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 602 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 609 of file DataOperators.hpp.

613  : cOords_at_GaussPt(coords_at_gausspt),
614  nOrmals_at_GaussPt(normals_at_gausspt),
615  tAngent1_at_GaussPt(tangent1_at_gausspt),
616  tAngent2_at_GaussPt(tangent2_at_gausspt) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnFace::calculateNormals ( )

Definition at line 1258 of file DataOperators.cpp.

1258  {
1260 
1261  nOrmals_at_GaussPt.resize(tAngent1_at_GaussPt.size1(), 3, false);
1262 
1263  auto get_ftensor1 = [](MatrixDouble &m) {
1265  &m(0, 0), &m(0, 1), &m(0, 2));
1266  };
1267  auto t_normal = get_ftensor1(nOrmals_at_GaussPt);
1268  auto t_t1 = get_ftensor1(tAngent1_at_GaussPt);
1269  auto t_t2 = get_ftensor1(tAngent2_at_GaussPt);
1270 
1274 
1275  for (unsigned int gg = 0; gg != tAngent1_at_GaussPt.size1(); ++gg) {
1276  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
1277  ++t_normal;
1278  ++t_t1;
1279  ++t_t2;
1280  }
1281 
1283 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
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 1158 of file DataOperators.cpp.

1159  {
1161 
1162  unsigned int nb_dofs = data.getFieldData().size();
1163  if (nb_dofs == 0)
1165 
1166  int nb_gauss_pts = data.getN().size1();
1167  cOords_at_GaussPt.resize(nb_gauss_pts, 3, false);
1168  tAngent1_at_GaussPt.resize(nb_gauss_pts, 3, false);
1169  tAngent2_at_GaussPt.resize(nb_gauss_pts, 3, false);
1170 
1171  auto get_ftensor1 = [](MatrixDouble &m) {
1173  &m(0, 0), &m(0, 1), &m(0, 2));
1174  };
1175  auto t_coords = get_ftensor1(cOords_at_GaussPt);
1176  auto t_t1 = get_ftensor1(tAngent1_at_GaussPt);
1177  auto t_t2 = get_ftensor1(tAngent2_at_GaussPt);
1179  FTensor::Number<0> N0;
1180  FTensor::Number<1> N1;
1181 
1182  switch (type) {
1183  case MBVERTEX: {
1184  cOords_at_GaussPt.clear();
1185  tAngent1_at_GaussPt.clear();
1186  tAngent2_at_GaussPt.clear();
1187  auto t_base = data.getFTensor0N();
1188  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1189  auto t_data = data.getFTensor1FieldData<3>();
1190  auto t_diff_base = data.getFTensor1DiffN<2>();
1191  for (int nn = 0; nn != 3; nn++) {
1192  t_coords(i) += t_base * t_data(i);
1193  t_t1(i) += t_data(i) * t_diff_base(N0);
1194  t_t2(i) += t_data(i) * t_diff_base(N1);
1195  ++t_data;
1196  ++t_base;
1197  ++t_diff_base;
1198  }
1199  ++t_coords;
1200  ++t_t1;
1201  ++t_t2;
1202  }
1203  } break;
1204  case MBEDGE:
1205  case MBTRI: {
1206  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1207  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1208  }
1209  if (nb_dofs % 3 != 0) {
1210  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1211  }
1212  if (nb_dofs > 3 * data.getN().size2()) {
1213  unsigned int nn = 0;
1214  for (; nn != nb_dofs; nn++) {
1215  if (!data.getFieldDofs()[nn]->getActive())
1216  break;
1217  }
1218  if (nn > 3 * data.getN().size2()) {
1219  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1220  "data inconsistency for base %s",
1221  ApproximationBaseNames[data.getBase()]);
1222  } else {
1223  nb_dofs = nn;
1224  if (!nb_dofs)
1226  }
1227  }
1228  const int nb_base_functions = data.getN().size2();
1229  auto t_base = data.getFTensor0N();
1230  auto t_diff_base = data.getFTensor1DiffN<2>();
1231  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1232  auto t_data = data.getFTensor1FieldData<3>();
1233  int bb = 0;
1234  for (; bb != nb_dofs / 3; ++bb) {
1235  t_coords(i) += t_base * t_data(i);
1236  t_t1(i) += t_data(i) * t_diff_base(N0);
1237  t_t2(i) += t_data(i) * t_diff_base(N1);
1238  ++t_data;
1239  ++t_base;
1240  ++t_diff_base;
1241  }
1242  for (; bb != nb_base_functions; ++bb) {
1243  ++t_base;
1244  ++t_diff_base;
1245  }
1246  ++t_coords;
1247  ++t_t1;
1248  ++t_t2;
1249  }
1250  } break;
1251  default:
1252  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1253  }
1254 
1256 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:499
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return() ...
Definition: definitions.h:506
static const char *const ApproximationBaseNames[]
Definition: definitions.h:156
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 604 of file DataOperators.hpp.

◆ nOrmals_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::nOrmals_at_GaussPt

Definition at line 605 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnFace::sPin

Definition at line 618 of file DataOperators.hpp.

◆ tAngent1_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::tAngent1_at_GaussPt

Definition at line 606 of file DataOperators.hpp.

◆ tAngent2_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::tAngent2_at_GaussPt

Definition at line 607 of file DataOperators.hpp.


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