v0.14.0
prism_elements_from_surface.cpp
Go to the documentation of this file.
1 /** \file prism_elements_from_surface.cpp
2  \example prism_elements_from_surface.cpp
3  \brief Adding prims on the surface and checking conformity between quads,
4  triangles and prism
5 
6 */
7 
8 
9 
10 #include <MoFEM.hpp>
11 
12 using namespace MoFEM;
13 
14 static char help[] = "...\n\n";
15 
16 static constexpr int precision_exponent = 5;
17 static constexpr int number_of_prisms_layers = 18;
18 static constexpr double delta =
19  1. / static_cast<double>(number_of_prisms_layers);
20 static constexpr std::array<double, 3> d3 = {0, 0, 0};
21 static constexpr std::array<double, 3> d4 = {0, 0, delta};
22 
24 
26 
27  inline static double getArg(double x) {
28  return std::round(x * pow(10., precision_exponent - 1));
29  };
30 
31  int x, y, z;
33  CoordsAndHandle(const double *coords, EntityHandle v)
34  : x(getArg(coords[0])), y(getArg(coords[1])), z(getArg(coords[2])),
35  node(v) {}
36 };
37 
38 typedef multi_index_container<
40  indexed_by<
41 
42  hashed_unique<composite_key<
43  CoordsAndHandle, member<CoordsAndHandle, int, &CoordsAndHandle::x>,
44  member<CoordsAndHandle, int, &CoordsAndHandle::y>,
45  member<CoordsAndHandle, int, &CoordsAndHandle::z>>>
46 
47  >>
49 
51 
52  PrismOp(moab::Interface &post_proc, MapCoords &map_coords);
53  MoFEMErrorCode doWork(int side, EntityType type,
55 
56 public:
59  std::vector<EntityHandle> nodeHandles;
60 };
61 
63 
64  PrismFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords);
65 
66  int getRuleTrianglesOnly(int order);
67  int getRuleThroughThickness(int order);
68 
69  MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only);
70 
71  MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness);
72 
73 private:
76 };
77 
78 template <typename OP> struct Op : public OP {
79 
80  Op(moab::Interface &post_proc, MapCoords &map_coords, EntityHandle prism);
81  MoFEMErrorCode doWork(int side, EntityType type,
83 
84 public:
88  std::vector<EntityHandle> nodeHandles;
89 };
90 
92 
93  TriFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords,
94  EntityHandle prims);
95 
96  int getRule(int order_row, int order_col, int order_data);
97 
98  MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data);
99 
100 private:
103 };
104 
106 
107  QuadFE(MoFEM::Interface &m_field, std::array<Range, 3> &edges_blocks,
108  EntityHandle prims);
109 
110  int getRule(int order_row, int order_col, int order_data);
111 
112  MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data);
113 
114 private:
115  std::array<Range, 3> &edgeBlocks;
117 };
118 struct EdgeFE : public EdgeEle {
119 
120  EdgeFE(MoFEM::Interface &m_field, std::array<Range, 3> &edges_blocks,
121  EntityHandle prims);
122 
123  int getRule(int order_row, int order_col, int order_data);
124 
125  MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data);
126 
127 private:
128  std::array<Range, 3> &edgeBlocks;
130 };
131 
132 int main(int argc, char *argv[]) {
133 
134  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
135 
136  try {
137 
138  moab::Core mb_instance;
139  moab::Interface &moab = mb_instance;
140  int rank;
141  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
142 
143  // Read parameters from line command
144  PetscBool flg = PETSC_TRUE;
145  char mesh_file_name[255];
146  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
147  255, &flg);
148  if (flg != PETSC_TRUE)
149  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
150  "error -my_file (MESH FILE NEEDED)");
151 
152  // Read mesh to MOAB
153  const char *option;
154  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
155  CHKERR moab.load_file(mesh_file_name, 0, option);
156 
157  Range tris;
158  CHKERR moab.get_entities_by_type(0, MBTRI, tris, false);
159  Range tris_verts;
160  CHKERR moab.get_connectivity(tris, tris_verts);
161  MatrixDouble tri_coords(tris_verts.size(), 3);
162  CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
163 
164  MoFEM::Core core(moab);
165  MoFEM::Interface &m_field = core;
166 
167  std::array<Range, 3> edge_block;
168  for (auto b : {1, 2, 3})
170  b, BLOCKSET, 1, edge_block[b - 1]);
171 
172  PrismsFromSurfaceInterface *prisms_from_surface_interface;
173  CHKERR m_field.getInterface(prisms_from_surface_interface);
174 
175  Range prisms;
176  CHKERR prisms_from_surface_interface->createPrisms(
178  prisms_from_surface_interface->setThickness(prisms, d3.data(), d4.data());
179  Range add_prims_layer;
180  Range extrude_prisms = prisms;
181 
182  for (int ll = 1; ll != number_of_prisms_layers; ++ll) {
183  prisms_from_surface_interface->createdVertices.clear();
184  CHKERR prisms_from_surface_interface->createPrismsFromPrisms(
185  extrude_prisms, false, add_prims_layer);
186  prisms_from_surface_interface->setThickness(add_prims_layer, d3.data(),
187  d4.data());
188  extrude_prisms = add_prims_layer;
189  prisms.merge(add_prims_layer);
190  add_prims_layer.clear();
191  }
192 
193  EntityHandle meshset;
194  CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
195  CHKERR moab.add_entities(meshset, prisms);
196 
197  MapCoords map_coords;
198  Range verts;
199  CHKERR moab.get_connectivity(prisms, verts);
200  MatrixDouble coords(verts.size(), 3);
201  CHKERR moab.get_coords(verts, &coords(0, 0));
202 
203  for (size_t v = 0; v != verts.size(); ++v)
204  map_coords.insert(CoordsAndHandle(&coords(v, 0), verts[v]));
205 
206  EntityHandle one_prism_meshset;
207  CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER,
208  one_prism_meshset);
209 
210  std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
211  0, 0, 1, 1, 0, 1, 0, 1, 1};
212  std::array<EntityHandle, 6> one_prism_nodes;
213  for (int n = 0; n != 6; ++n)
214  CHKERR moab.create_vertex(&one_prism_coords[3 * n], one_prism_nodes[n]);
215  EntityHandle one_prism;
216  CHKERR m_field.get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
217  one_prism);
218  Range one_prism_range;
219  one_prism_range.insert(one_prism);
220  CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
221  Range one_prism_verts;
222  CHKERR moab.get_connectivity(one_prism_range, one_prism_verts);
223  CHKERR moab.add_entities(one_prism_meshset, one_prism_verts);
224  Range one_prism_adj_ents;
225  for (int d = 1; d != 3; ++d)
226  CHKERR moab.get_adjacencies(one_prism_range, d, true, one_prism_adj_ents,
227  moab::Interface::UNION);
228  CHKERR moab.add_entities(one_prism_meshset, one_prism_adj_ents);
229 
230  BitRefLevel bit_level0;
231  bit_level0.set(0);
232  CHKERR m_field.getInterface<BitRefManager>()->setEntitiesBitRefLevel(
233  one_prism_range, bit_level0);
234  CHKERR prisms_from_surface_interface->seedPrismsEntities(one_prism_range,
235  bit_level0);
236 
237  // Fields
238  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
239  CHKERR m_field.add_ents_to_field_by_type(one_prism_meshset, MBPRISM,
240  "FIELD1", VERY_NOISY);
241 
242  CHKERR m_field.set_field_order(one_prism_meshset, MBVERTEX, "FIELD1", 1);
243  CHKERR m_field.set_field_order(one_prism_meshset, MBEDGE, "FIELD1", 5,
244  VERY_NOISY);
245  CHKERR m_field.set_field_order(one_prism_meshset, MBTRI, "FIELD1", 5);
246  CHKERR m_field.set_field_order(one_prism_meshset, MBQUAD, "FIELD1", 5,
247  VERY_NOISY);
248  CHKERR m_field.set_field_order(one_prism_meshset, MBPRISM, "FIELD1", 7,
249  VERY_NOISY);
250  CHKERR m_field.build_fields(VERY_NOISY);
251 
252  // FE
253  CHKERR m_field.add_finite_element("PRISM");
254  CHKERR m_field.add_finite_element("EDGE");
255  CHKERR m_field.add_finite_element("TRI");
256  CHKERR m_field.add_finite_element("QUAD");
257 
258  // Define rows/cols and element data
259  CHKERR m_field.modify_finite_element_add_field_row("PRISM", "FIELD1");
260  CHKERR m_field.modify_finite_element_add_field_col("PRISM", "FIELD1");
261  CHKERR m_field.modify_finite_element_add_field_data("PRISM", "FIELD1");
262  // Define rows/cols and element data
263  CHKERR m_field.modify_finite_element_add_field_row("EDGE", "FIELD1");
264  CHKERR m_field.modify_finite_element_add_field_col("EDGE", "FIELD1");
265  CHKERR m_field.modify_finite_element_add_field_data("EDGE", "FIELD1");
266  // Define rows/cols and element data
267  CHKERR m_field.modify_finite_element_add_field_row("TRI", "FIELD1");
268  CHKERR m_field.modify_finite_element_add_field_col("TRI", "FIELD1");
269  CHKERR m_field.modify_finite_element_add_field_data("TRI", "FIELD1");
270  // Define rows/cols and element data
271  CHKERR m_field.modify_finite_element_add_field_row("QUAD", "FIELD1");
272  CHKERR m_field.modify_finite_element_add_field_col("QUAD", "FIELD1");
273  CHKERR m_field.modify_finite_element_add_field_data("QUAD", "FIELD1");
274 
275  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBEDGE,
276  "EDGE");
277  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBTRI,
278  "TRI");
279  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBQUAD,
280  "QUAD");
281  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset,
282  MBPRISM, "PRISM");
283 
284  // build finite elemnts
285  CHKERR m_field.build_finite_elements();
286  // //build adjacencies
287  CHKERR m_field.build_adjacencies(bit_level0);
288 
289  // Problem
290  CHKERR m_field.add_problem("TEST_PROBLEM");
291 
292  // set finite elements for problem
293  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "PRISM");
294  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE");
295  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TRI");
296  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "QUAD");
297  // set refinement level for problem
298  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
299 
300  // build problem
301  ProblemsManager *prb_mng_ptr;
302  CHKERR m_field.getInterface(prb_mng_ptr);
303  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
304  // partition
305  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
306  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
307  // what are ghost nodes, see Petsc Manual
308  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
309 
310  PrismFE fe_prism(m_field, tri_coords);
311  fe_prism.getOpPtrVector().push_back(new PrismOp(moab, map_coords));
312  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe_prism);
313 
314  EdgeFE fe_edge(m_field, edge_block, one_prism);
315  fe_edge.getOpPtrVector().push_back(
317  moab, map_coords, one_prism));
318  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE", fe_edge);
319 
320  TriFE fe_tri(m_field, tri_coords, one_prism);
321  fe_tri.getOpPtrVector().push_back(
323  moab, map_coords, one_prism));
324  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TRI", fe_tri);
325 
326  QuadFE fe_quad(m_field, edge_block, one_prism);
327  fe_quad.getOpPtrVector().push_back(
329  moab, map_coords, one_prism));
330  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "QUAD", fe_quad);
331 
332  CHKERR moab.write_file("prism_mesh.vtk", "VTK", "", &meshset, 1);
333  CHKERR moab.write_file("one_prism_mesh.vtk", "VTK", "", &one_prism_meshset,
334  1);
335  }
336  CATCH_ERRORS;
337 
339  return 0;
340 }
341 
344  "FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
345  postProc(post_proc), mapCoords(map_coords) {}
346 
347 MoFEMErrorCode PrismOp::doWork(int side, EntityType type,
349  constexpr double def_val[] = {0, 0, 0};
351  switch (type) {
352  case MBVERTEX:
353  case MBEDGE:
354  case MBTRI:
355  case MBQUAD:
356  case MBPRISM:
357  break;
358  default:
360  }
361  if (type == MBTRI && (side != 3 && side != 4))
363  if (type == MBQUAD && (side == 3 || side == 4))
365 
366  const int nb_dofs = data.getIndices().size();
367  for (int dd = 0; dd != nb_dofs; ++dd)
368  data.getFieldDofs()[dd]->getFieldData() = data.getIndices()[dd];
369 
370  if (type == MBVERTEX) {
371  const size_t nb_gauss_pts = getGaussPts().size2();
372  auto &coords_at_pts = getGaussPts();
373  nodeHandles.reserve(nb_gauss_pts);
374  nodeHandles.clear();
375  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
376  auto t = boost::make_tuple(CoordsAndHandle::getArg(coords_at_pts(0, gg)),
377  CoordsAndHandle::getArg(coords_at_pts(1, gg)),
378  CoordsAndHandle::getArg(coords_at_pts(2, gg)));
379 
380  auto it = mapCoords.find(t);
381  if (it != mapCoords.end())
382  nodeHandles.emplace_back(it->node);
383  else
384  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vertex not found");
385  }
386 
387  MatrixDouble node_coords(nb_gauss_pts, 3);
388  CHKERR postProc.get_coords(&nodeHandles[0], nodeHandles.size(),
389  &node_coords(0, 0));
390  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
391  constexpr double eps = 1e-12;
392  for (auto d : {0, 1, 2})
393  if (std::abs(node_coords(gg, d) - getGaussPts()(d, gg)) > eps) {
394  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
395  "Inconsistency between node coords and integration "
396  "points");
397  }
398  for (auto d : {0, 1, 2})
399  if (std::abs(node_coords(gg, d) - getCoordsAtGaussPts()(gg, d)) > eps)
400  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
401  "Inconsistency between node coords and integration "
402  "points");
403  }
404 
405  Tag th;
406  CHKERR postProc.tag_get_handle("Coords", 3, MB_TYPE_DOUBLE, th,
407  MB_TAG_CREAT | MB_TAG_DENSE, def_val);
408  CHKERR postProc.tag_set_data(th, &nodeHandles[0], nodeHandles.size(),
409  &getCoordsAtGaussPts()(0, 0));
410  }
411 
412  auto to_str = [](auto i) { return boost::lexical_cast<std::string>(i); };
413  std::string tag_name_base =
414  "PrismType" + to_str(type) + "Side" + to_str(side);
415  std::cout << tag_name_base << endl;
416 
417  MatrixDouble trans_base = trans(data.getN());
418  if (trans_base.size2() != nodeHandles.size())
419  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "wrong size %d != %d",
420  trans_base.size2(), nodeHandles.size());
421  for (size_t rr = 0; rr != trans_base.size1(); ++rr) {
422  auto tag_name = tag_name_base + "Base" + to_str(rr);
423  Tag th;
424  CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE, th,
425  MB_TAG_CREAT | MB_TAG_DENSE, def_val);
426  CHKERR postProc.tag_set_data(th, &nodeHandles[0], nodeHandles.size(),
427  &trans_base(rr, 0));
428  }
429 
431 }
432 
434  : FatPrismElementForcesAndSourcesCore(m_field), triCoords(tri_coords) {}
435 
436 int PrismFE::getRuleTrianglesOnly(int order) { return -1; };
437 int PrismFE::getRuleThroughThickness(int order) { return -1; };
438 
441  gaussPtsTrianglesOnly.resize(3, triCoords.size1(), false);
442  gaussPtsTrianglesOnly.clear();
443  for (int gg = 0; gg != triCoords.size1(); ++gg)
444  for (int dd = 0; dd != 2; ++dd)
445  gaussPtsTrianglesOnly(dd, gg) = triCoords(gg, dd);
447 }
448 
452  gaussPtsThroughThickness.clear();
453  for (int gg = 0; gg != number_of_prisms_layers + 1; ++gg)
454  gaussPtsThroughThickness(0, gg) = delta * gg;
455 
457 }
458 
460  EntityHandle prism)
461  : FaceElementForcesAndSourcesCore(m_field), triCoords(tri_coords),
462  prism(prism) {}
463 
464 int TriFE::getRule(int order_row, int order_col, int order_data) { return -1; }
465 
466 MoFEMErrorCode TriFE::setGaussPts(int order_row, int order_col,
467  int order_data) {
469 
470  const EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
471  int side, sense, offset;
472  CHKERR mField.get_moab().side_number(prism, ent, side, sense, offset);
473  std::array<int, 2> swap = {0, 1};
474  if (side == 3)
475  swap = std::array<int, 2>{1, 0};
476 
477  gaussPts.resize(3, triCoords.size1(), false);
478  gaussPts.clear();
479  for (int gg = 0; gg != triCoords.size1(); ++gg)
480  for (int dd = 0; dd != 2; ++dd)
481  gaussPts(dd, gg) = triCoords(gg, swap[dd]);
482 
484 }
485 
486 QuadFE::QuadFE(MoFEM::Interface &m_field, std::array<Range, 3> &edge_blocks,
487  EntityHandle prism)
488  : FaceElementForcesAndSourcesCore(m_field), edgeBlocks(edge_blocks),
489  prism(prism) {}
490 
491 int QuadFE::getRule(int order_row, int order_col, int order_data) { return -1; }
492 
493 MoFEMErrorCode QuadFE::setGaussPts(int order_row, int order_col,
494  int order_data) {
496 
497  const EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
498  int side, sense, offset;
499  CHKERR mField.get_moab().side_number(prism, ent, side, sense, offset);
500 
501  Range edge_verts;
502  CHKERR mField.get_moab().get_connectivity(edgeBlocks[side], edge_verts);
503  MatrixDouble edge_coords(edge_verts.size(), 3);
504  CHKERR mField.get_moab().get_coords(edge_verts, &edge_coords(0, 0));
505 
506  constexpr double normal[3][2] = {{-1, 0}, {-0.5, 0.5}, {0, 1}};
507  constexpr double origin[3][2] = {{1, 0}, {1, 0}, {0, 0}};
508  constexpr int swap[3][2] = {{0, 1}, {0, 1}, {1, 0}};
509  gaussPts.resize(3, edge_verts.size() * (number_of_prisms_layers + 1), false);
510  int gg = 0;
511  for (size_t rr = 0; rr != edge_verts.size(); ++rr) {
512  const double x = edge_coords(rr, 0) - origin[side][0];
513  const double y = edge_coords(rr, 1) - origin[side][1];
514  const double edge_dist = x * normal[side][0] + y * normal[side][1];
515  for (size_t cc = 0; cc != number_of_prisms_layers + 1; ++cc, ++gg) {
516  gaussPts(swap[side][0], gg) = edge_dist;
517  gaussPts(swap[side][1], gg) = delta * cc;
518  }
519  }
520 
522 }
523 
524 EdgeFE::EdgeFE(MoFEM::Interface &m_field, std::array<Range, 3> &edge_blocks,
525  EntityHandle prism)
526  : EdgeEle(m_field), edgeBlocks(edge_blocks),
527  prism(prism) {}
528 
529 int EdgeFE::getRule(int order_row, int order_col, int order_data) { return -1; }
530 
531 MoFEMErrorCode EdgeFE::setGaussPts(int order_row, int order_col,
532  int order_data) {
534 
535  const EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
536  int side, sense, offset;
537  CHKERR mField.get_moab().side_number(prism, ent, side, sense, offset);
538 
539  if (side >= 3 && side <= 5) {
540 
541  gaussPts.resize(2, number_of_prisms_layers + 1, false);
542  for (size_t gg = 0; gg != number_of_prisms_layers + 1; ++gg)
543  gaussPts(0, gg) = delta * gg;
544 
545  // cerr << gaussPts << endl;
546 
547  } else {
548 
549  constexpr double normal[3][2] = {{1, 0}, {-0.5, 0.5}, {0, 1}};
550  constexpr double origin[3][2] = {{0, 1}, {1, 0}, {0, 0}};
551  constexpr int side_map[9] = {0, 1, 2, -1, -1, -1, 0, 1, 2};
552 
553  int num_nodes;
554  const EntityHandle *conn;
555  CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
556  MatrixDouble coords(num_nodes, 3);
557  CHKERR mField.get_moab().get_coords(conn, num_nodes,
558  &*coords.data().begin());
559  side = side_map[side];
560 
561  Range edge_verts;
562  CHKERR mField.get_moab().get_connectivity(edgeBlocks[side], edge_verts);
563  MatrixDouble edge_coords(edge_verts.size(), 3);
564  CHKERR mField.get_moab().get_coords(edge_verts, &edge_coords(0, 0));
565  gaussPts.resize(2, edge_verts.size(), false);
566 
567  for (size_t gg = 0; gg != edge_verts.size(); ++gg) {
568  const double x = edge_coords(gg, 0) - origin[side][0];
569  const double y = edge_coords(gg, 1) - origin[side][1];
570  const double edge_dist = x * normal[side][0] + y * normal[side][1];
571  gaussPts(0, gg) = edge_dist;
572  }
573  }
574 
576 }
577 
578 template <typename OP>
579 Op<OP>::Op(moab::Interface &post_proc, MapCoords &map_coords,
580  EntityHandle prism)
581  : OP("FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
582  postProc(post_proc), mapCoords(map_coords), prism(prism) {}
583 
584 template <typename OP>
585 MoFEMErrorCode Op<OP>::doWork(int side, EntityType type,
587  constexpr double def_val[] = {0, 0, 0};
589  switch (type) {
590  case MBVERTEX:
591  case MBEDGE:
592  case MBTRI:
593  case MBQUAD:
594  break;
595  default:
597  }
598 
599  const int nb_dofs = data.getIndices().size();
600  for (int dd = 0; dd != nb_dofs; ++dd)
601  if (data.getFieldData()[dd] != data.getIndices()[dd]) {
602  std::cerr << "Indices: " << data.getIndices() << std::endl;
603  std::cerr << "Local indices: " << data.getLocalIndices() << std::endl;
604  std::cerr << "Data: " << data.getFieldData() << std::endl;
605  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
606  "Indicices/data inconsistency %3.1f != %d",
607  data.getFieldData()[dd], data.getIndices()[dd]);
608  }
609 
610  const EntityHandle fe_ent = OP::getFEEntityHandle();
611  const EntityHandle ent = OP::getSideEntity(side, type);
612  int side_prism, sense, offset;
613  if (type == MBVERTEX) {
614  CHKERR postProc.side_number(prism, fe_ent, side_prism, sense, offset);
615  } else
616  CHKERR postProc.side_number(prism, ent, side_prism, sense, offset);
617 
618  if (type == MBVERTEX) {
619  auto &coords_at_pts = OP::getCoordsAtGaussPts();
620  const size_t nb_gauss_pts = coords_at_pts.size1();
621 
622  nodeHandles.reserve(nb_gauss_pts);
623  nodeHandles.clear();
624  for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
625  auto t = boost::make_tuple(CoordsAndHandle::getArg(coords_at_pts(gg, 0)),
626  CoordsAndHandle::getArg(coords_at_pts(gg, 1)),
627  CoordsAndHandle::getArg(coords_at_pts(gg, 2)));
628 
629  auto it = mapCoords.find(t);
630  if (it != mapCoords.end())
631  nodeHandles.emplace_back(it->node);
632  else
633  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vertex not found");
634  }
635  }
636 
637  auto to_str = [](auto i) { return boost::lexical_cast<std::string>(i); };
638  std::string tag_name_base =
639  "FEType" + to_str(OP::getNumeredEntFiniteElementPtr()->getEntType()) +
640  "Type" + to_str(type) + "Side" + to_str(side_prism);
641 
642  std::string tag_prism_name_base =
643  "PrismType" + to_str(type) + "Side" + to_str(side_prism);
644 
645  MatrixDouble trans_base = trans(data.getN());
646  MatrixDouble prism_base(trans_base.size1(), trans_base.size2());
647  if (trans_base.size2() != nodeHandles.size())
648  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "wrong size %d != %d",
649  trans_base.size2(), nodeHandles.size());
650 
651  for (size_t rr = 0; rr != trans_base.size1(); ++rr) {
652 
653  std::string tag_name = tag_name_base + "Base";
654  if (type == MBVERTEX) {
655  EntityHandle node = data.getFieldDofs()[rr]->getEnt();
656  int prism_mode_side;
657  CHKERR postProc.side_number(prism, node, prism_mode_side, sense, offset);
658  tag_name += to_str(prism_mode_side);
659  } else {
660  tag_name += to_str(rr);
661  }
662 
663  std::cout << tag_name << endl;
664 
665  Tag th;
666  CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE, th,
667  MB_TAG_CREAT | MB_TAG_DENSE, def_val);
668  CHKERR postProc.tag_set_data(th, &nodeHandles[0], nodeHandles.size(),
669  &trans_base(rr, 0));
670 
671  if (type != MBVERTEX) {
672  auto tag_prism_name = tag_prism_name_base + "Base" + to_str(rr);
673  Tag th_prism;
674  CHKERR postProc.tag_get_handle(tag_prism_name.c_str(), th_prism);
675  CHKERR postProc.tag_get_data(th_prism, &nodeHandles[0],
676  nodeHandles.size(), &prism_base(rr, 0));
677  }
678  }
679 
680  auto sum_matrix = [](MatrixDouble &m) {
681  double s = 0;
682  for (unsigned int ii = 0; ii < m.size1(); ii++) {
683  for (unsigned int jj = 0; jj < m.size2(); jj++) {
684  s += std::abs(m(ii, jj));
685  }
686  }
687  return s;
688  };
689 
690  if (type != MBVERTEX) {
691  prism_base -= trans_base;
692  double sum = sum_matrix(prism_base);
693  constexpr double eps = 1e-6;
694 
695  if (std::abs(sum) > eps)
696  cout << "Inconsistent base " << tag_prism_name_base << " "
697  << tag_name_base << " sum " << sum << endl;
698 
699  if (std::abs(sum) > eps)
700  SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
701  "Inconsistent base %s sum %6.4e", tag_prism_name_base.c_str(),
702  sum);
703  }
704 
706 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
main
int main(int argc, char *argv[])
Definition: prism_elements_from_surface.cpp:132
MoFEM::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
EdgeFE::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:129
QuadFE::edgeBlocks
std::array< Range, 3 > & edgeBlocks
Definition: prism_elements_from_surface.cpp:115
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::PrismsFromSurfaceInterface
merge node from two bit levels
Definition: PrismsFromSurfaceInterface.hpp:15
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
TriFE::setGaussPts
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: prism_elements_from_surface.cpp:466
H1
@ H1
continuous field
Definition: definitions.h:85
PrismOp::PrismOp
PrismOp(moab::Interface &post_proc, MapCoords &map_coords)
Definition: prism_elements_from_surface.cpp:342
MoFEM::PetscData::x
Vec x
Definition: LoopMethods.hpp:51
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator
default operator for Flat Prism element
Definition: FatPrismElementForcesAndSourcesCore.hpp:54
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
PrismFE::setGaussPtsThroughThickness
MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
Definition: prism_elements_from_surface.cpp:449
PrismFE::setGaussPtsTrianglesOnly
MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
Definition: prism_elements_from_surface.cpp:439
d3
static constexpr std::array< double, 3 > d3
Definition: prism_elements_from_surface.cpp:20
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
MoFEM::EntitiesFieldData::EntData::getLocalIndices
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
Definition: EntitiesFieldData.hpp:1229
PrismOp::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: prism_elements_from_surface.cpp:347
PrismFE::getRuleTrianglesOnly
int getRuleTrianglesOnly(int order)
Definition: prism_elements_from_surface.cpp:436
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
PrismOp
Definition: prism_elements_from_surface.cpp:50
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
VERY_NOISY
@ VERY_NOISY
Definition: definitions.h:225
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1254
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
MatrixDouble gaussPtsTrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:189
EdgeFE::EdgeFE
EdgeFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
Definition: prism_elements_from_surface.cpp:524
EdgeFE::getRule
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: prism_elements_from_surface.cpp:529
MoFEM.hpp
EdgeFE
Definition: prism_elements_from_surface.cpp:118
CoordsAndHandle::node
EntityHandle node
Definition: prism_elements_from_surface.cpp:32
TriFE
Definition: prism_elements_from_surface.cpp:91
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
PrismOp::postProc
moab::Interface & postProc
Definition: prism_elements_from_surface.cpp:57
TriFE::TriFE
TriFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords, EntityHandle prims)
Definition: prism_elements_from_surface.cpp:459
QuadFE::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:116
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
CoordsAndHandle
Definition: prism_elements_from_surface.cpp:25
QuadFE::getRule
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: prism_elements_from_surface.cpp:491
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
Op::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Definition: prism_elements_from_surface.cpp:585
TriFE::getRule
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: prism_elements_from_surface.cpp:464
Op::Op
Op(moab::Interface &post_proc, MapCoords &map_coords, EntityHandle prism)
Definition: prism_elements_from_surface.cpp:579
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
TriFE::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:102
CoordsAndHandle::getArg
static double getArg(double x)
Definition: prism_elements_from_surface.cpp:27
number_of_prisms_layers
static constexpr int number_of_prisms_layers
Definition: prism_elements_from_surface.cpp:17
PrismFE::triCoords
MatrixDouble & triCoords
Definition: prism_elements_from_surface.cpp:74
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
QuadFE::QuadFE
QuadFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
Definition: prism_elements_from_surface.cpp:486
PrismFE::getRuleThroughThickness
int getRuleThroughThickness(int order)
Definition: prism_elements_from_surface.cpp:437
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
delta
static constexpr double delta
Definition: prism_elements_from_surface.cpp:18
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
MoFEM::MeshsetsManager::getEntitiesByDimension
MoFEMErrorCode getEntitiesByDimension(const int ms_id, const unsigned int cubit_bc_type, const int dimension, Range &entities, const bool recursive=true) const
get entities from CUBIT/meshset of a particular entity dimension
Definition: MeshsetsManager.cpp:669
MoFEM::PrismsFromSurfaceInterface::seedPrismsEntities
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
Definition: PrismsFromSurfaceInterface.cpp:110
PrismFE::prims
EntityHandle prims
Definition: prism_elements_from_surface.cpp:75
convert.type
type
Definition: convert.py:64
MoFEM::PrismsFromSurfaceInterface::createPrismsFromPrisms
MoFEMErrorCode createPrismsFromPrisms(const Range &prisms, bool from_down, Range &out_prisms, int verb=-1)
Make prisms by extruding top or bottom prisms.
Definition: PrismsFromSurfaceInterface.cpp:143
EdgeFE::setGaussPts
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: prism_elements_from_surface.cpp:531
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1214
help
static char help[]
Definition: prism_elements_from_surface.cpp:14
precision_exponent
static constexpr int precision_exponent
Definition: prism_elements_from_surface.cpp:16
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
PrismFE
Definition: prism_elements_from_surface.cpp:62
MoFEM::EntitiesFieldData::EntData::getFieldDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
Definition: EntitiesFieldData.hpp:1269
MoFEM::FatPrismElementForcesAndSourcesCore
FatPrism finite element.
Definition: FatPrismElementForcesAndSourcesCore.hpp:31
t
constexpr double t
plate stiffness
Definition: plate.cpp:58
Op::postProc
moab::Interface & postProc
Definition: prism_elements_from_surface.cpp:85
TriFE::triCoords
MatrixDouble & triCoords
Definition: prism_elements_from_surface.cpp:101
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::ProblemsManager::partitionFiniteElements
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
Definition: ProblemsManager.cpp:2167
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts
MatrixDouble & getCoordsAtGaussPts()
get coordinates at Gauss pts.
Definition: FatPrismElementForcesAndSourcesCore.hpp:260
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
Op::nodeHandles
std::vector< EntityHandle > nodeHandles
Definition: prism_elements_from_surface.cpp:88
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
convert.n
n
Definition: convert.py:82
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
MatrixDouble gaussPtsThroughThickness
Definition: FatPrismElementForcesAndSourcesCore.hpp:191
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
edge_coords
static const double edge_coords[6][6]
Definition: forces_and_sources_h1_continuity_check.cpp:18
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
CoordsAndHandle::CoordsAndHandle
CoordsAndHandle(const double *coords, EntityHandle v)
Definition: prism_elements_from_surface.cpp:33
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
FTensor::dd
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
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
PrismOp::nodeHandles
std::vector< EntityHandle > nodeHandles
Definition: prism_elements_from_surface.cpp:59
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1318
PrismFE::PrismFE
PrismFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords)
Definition: prism_elements_from_surface.cpp:433
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
Op
Definition: prism_elements_from_surface.cpp:78
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
MoFEM::PrismsFromSurfaceInterface::NO_SWAP
@ NO_SWAP
Definition: PrismsFromSurfaceInterface.hpp:38
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
QuadFE
Definition: prism_elements_from_surface.cpp:105
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
sdf_hertz_2d_axisymm_plane.d
float d
Definition: sdf_hertz_2d_axisymm_plane.py:4
CoordsAndHandle::z
int z
Definition: prism_elements_from_surface.cpp:31
MoFEM::PrismsFromSurfaceInterface::createdVertices
std::map< EntityHandle, EntityHandle > createdVertices
Definition: PrismsFromSurfaceInterface.hpp:24
d4
static constexpr std::array< double, 3 > d4
Definition: prism_elements_from_surface.cpp:21
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
QuadFE::setGaussPts
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: prism_elements_from_surface.cpp:493
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEM::PrismsFromSurfaceInterface::createPrisms
MoFEMErrorCode createPrisms(const Range &ents, const SwapType swap_type, Range &prisms, int verb=-1)
Make prisms from triangles.
Definition: PrismsFromSurfaceInterface.cpp:22
MoFEM::CoreInterface::set_field_order
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
sum_matrix
static double sum_matrix(MatrixDouble &m)
Definition: edge_and_bubble_shape_functions_on_quad.cpp:12
PrismOp::mapCoords
MapCoords & mapCoords
Definition: prism_elements_from_surface.cpp:58
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
OP
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
Op::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:87
EdgeFE::edgeBlocks
std::array< Range, 3 > & edgeBlocks
Definition: prism_elements_from_surface.cpp:128
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
get Gauss pts. in the prism
Definition: FatPrismElementForcesAndSourcesCore.hpp:243
MoFEM::CoreInterface::add_field
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
MoFEM::PrismsFromSurfaceInterface::setThickness
MoFEMErrorCode setThickness(const Range &prisms, const double director3[], const double director4[])
Definition: PrismsFromSurfaceInterface.cpp:161
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:24
MapCoords
multi_index_container< CoordsAndHandle, indexed_by< hashed_unique< composite_key< CoordsAndHandle, member< CoordsAndHandle, int, &CoordsAndHandle::x >, member< CoordsAndHandle, int, &CoordsAndHandle::y >, member< CoordsAndHandle, int, &CoordsAndHandle::z > > > > > MapCoords
Definition: prism_elements_from_surface.cpp:48
Op::mapCoords
MapCoords & mapCoords
Definition: prism_elements_from_surface.cpp:86