v0.8.23
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 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_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 450 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 457 of file DataOperators.hpp.

461  : cOords_at_GaussPt(coords_at_gausspt),
462  nOrmals_at_GaussPt(normals_at_gausspt),
463  tAngent1_at_GaussPt(tangent1_at_gausspt),
464  tAngent2_at_GaussPt(tangent2_at_gausspt) {}

Member Function Documentation

◆ calculateNormals()

MoFEMErrorCode MoFEM::OpGetCoordsAndNormalsOnFace::calculateNormals ( )

Definition at line 1223 of file DataOperators.cpp.

1223  {
1225 
1226  nOrmals_at_GaussPt.resize(tAngent1_at_GaussPt.size1(), 3, false);
1227 
1228  auto get_ftensor1 = [](MatrixDouble &m) {
1230  &m(0, 0), &m(0, 1), &m(0, 2));
1231  };
1232  auto t_normal = get_ftensor1(nOrmals_at_GaussPt);
1233  auto t_t1 = get_ftensor1(tAngent1_at_GaussPt);
1234  auto t_t2 = get_ftensor1(tAngent2_at_GaussPt);
1235 
1239 
1240  for (unsigned int gg = 0; gg != tAngent1_at_GaussPt.size1(); ++gg) {
1241  t_normal(j) = FTensor::levi_civita(i, j, k) * t_t1(k) * t_t2(i);
1242  ++t_normal;
1243  ++t_t1;
1244  ++t_t2;
1245  }
1246 
1248 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:77
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
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 right hand side.

Reimplemented from MoFEM::DataOperator.

Definition at line 1123 of file DataOperators.cpp.

1124  {
1126 
1127  unsigned int nb_dofs = data.getFieldData().size();
1128  if (nb_dofs == 0)
1130 
1131  int nb_gauss_pts = data.getN().size1();
1132  cOords_at_GaussPt.resize(nb_gauss_pts, 3, false);
1133  tAngent1_at_GaussPt.resize(nb_gauss_pts, 3, false);
1134  tAngent2_at_GaussPt.resize(nb_gauss_pts, 3, false);
1135 
1136  auto get_ftensor1 = [](MatrixDouble &m) {
1138  &m(0, 0), &m(0, 1), &m(0, 2));
1139  };
1140  auto t_coords = get_ftensor1(cOords_at_GaussPt);
1141  auto t_t1 = get_ftensor1(tAngent1_at_GaussPt);
1142  auto t_t2 = get_ftensor1(tAngent2_at_GaussPt);
1144  FTensor::Number<0> N0;
1145  FTensor::Number<1> N1;
1146 
1147  switch (type) {
1148  case MBVERTEX: {
1149  cOords_at_GaussPt.clear();
1150  tAngent1_at_GaussPt.clear();
1151  tAngent2_at_GaussPt.clear();
1152  auto t_base = data.getFTensor0N();
1153  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1154  auto t_data = data.getFTensor1FieldData<3>();
1155  auto t_diff_base = data.getFTensor1DiffN<2>();
1156  for (int nn = 0; nn != 3; nn++) {
1157  t_coords(i) += t_base * t_data(i);
1158  t_t1(i) += t_data(i) * t_diff_base(N0);
1159  t_t2(i) += t_data(i) * t_diff_base(N1);
1160  ++t_data;
1161  ++t_base;
1162  ++t_diff_base;
1163  }
1164  ++t_coords;
1165  ++t_t1;
1166  ++t_t2;
1167  }
1168  } break;
1169  case MBEDGE:
1170  case MBTRI: {
1171  if (2 * data.getN().size2() != data.getDiffN().size2()) {
1172  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1173  }
1174  if (nb_dofs % 3 != 0) {
1175  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "data inconsistency");
1176  }
1177  if (nb_dofs > 3 * data.getN().size2()) {
1178  unsigned int nn = 0;
1179  for (; nn != nb_dofs; nn++) {
1180  if (!data.getFieldDofs()[nn]->getActive())
1181  break;
1182  }
1183  if (nn > 3 * data.getN().size2()) {
1184  SETERRQ1(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1185  "data inconsistency for base %s",
1186  ApproximationBaseNames[data.getBase()]);
1187  } else {
1188  nb_dofs = nn;
1189  if (!nb_dofs)
1191  }
1192  }
1193  const int nb_base_functions = data.getN().size2();
1194  auto t_base = data.getFTensor0N();
1195  auto t_diff_base = data.getFTensor1DiffN<2>();
1196  for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1197  auto t_data = data.getFTensor1FieldData<3>();
1198  int bb = 0;
1199  for (; bb != nb_dofs / 3; ++bb) {
1200  t_coords(i) += t_base * t_data(i);
1201  t_t1(i) += t_data(i) * t_diff_base(N0);
1202  t_t2(i) += t_data(i) * t_diff_base(N1);
1203  ++t_data;
1204  ++t_base;
1205  ++t_diff_base;
1206  }
1207  for (; bb != nb_base_functions; ++bb) {
1208  ++t_base;
1209  ++t_diff_base;
1210  }
1211  ++t_coords;
1212  ++t_t1;
1213  ++t_t2;
1214  }
1215  } break;
1216  default:
1217  SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "not implemented");
1218  }
1219 
1221 }
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:500
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:77
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:507
static const char *const ApproximationBaseNames[]
Definition: definitions.h:157

Member Data Documentation

◆ cOords_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::cOords_at_GaussPt

Definition at line 452 of file DataOperators.hpp.

◆ nOrmals_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::nOrmals_at_GaussPt

Definition at line 453 of file DataOperators.hpp.

◆ sPin

MatrixDouble MoFEM::OpGetCoordsAndNormalsOnFace::sPin

Definition at line 466 of file DataOperators.hpp.

◆ tAngent1_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::tAngent1_at_GaussPt

Definition at line 454 of file DataOperators.hpp.

◆ tAngent2_at_GaussPt

MatrixDouble& MoFEM::OpGetCoordsAndNormalsOnFace::tAngent2_at_GaussPt

Definition at line 455 of file DataOperators.hpp.


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