v0.15.0
Loading...
Searching...
No Matches
CutMeshInterface.cpp
Go to the documentation of this file.
1/** \file CutMeshInterface.cpp
2 * \brief Cut mesh by surface
3 */
4
5
6
7#define CutMeshFunctionBegin \
8 MoFEMFunctionBegin; \
9 MOFEM_LOG_CHANNEL("WORLD"); \
10 MOFEM_LOG_FUNCTION(); \
11 MOFEM_LOG_TAG("WORLD", "CutMesh");
12
13namespace MoFEM {
14
16CutMeshInterface::query_interface(boost::typeindex::type_index type_index,
17 UnknownInterface **iface) const {
19 *iface = const_cast<CutMeshInterface *>(this);
21}
22
24 : cOre(const_cast<Core &>(core)) {
25 lineSearchSteps = 10;
28}
29
35
41
47
49 double *shift, double *origin,
50 double *transform,
51 const std::string save_mesh) {
52 CoreInterface &m_field = cOre;
53 moab::Interface &moab = m_field.get_moab();
55 sUrface.clear();
56 std::map<EntityHandle, EntityHandle> verts_map;
57 for (Range::const_iterator tit = surface.begin(); tit != surface.end();
58 tit++) {
59 int num_nodes;
60 const EntityHandle *conn;
61 CHKERR moab.get_connectivity(*tit, conn, num_nodes, true);
62 MatrixDouble coords(num_nodes, 3);
63 if (th) {
64 CHKERR moab.tag_get_data(th, conn, num_nodes, &coords(0, 0));
65 } else {
66 CHKERR moab.get_coords(conn, num_nodes, &coords(0, 0));
67 }
68 EntityHandle new_verts[num_nodes];
69 for (int nn = 0; nn != num_nodes; nn++) {
70 if (verts_map.find(conn[nn]) != verts_map.end()) {
71 new_verts[nn] = verts_map[conn[nn]];
72 } else {
73 if (transform) {
74 ublas::matrix_row<MatrixDouble> mr(coords, nn);
75 if (origin) {
76 VectorAdaptor vec_origin(
77 3, ublas::shallow_array_adaptor<double>(3, origin));
78 mr = mr - vec_origin;
79 }
80 MatrixAdaptor mat_transform = MatrixAdaptor(
81 3, 3, ublas::shallow_array_adaptor<double>(9, transform));
82 mr = prod(mat_transform, mr);
83 if (origin) {
84 VectorAdaptor vec_origin(
85 3, ublas::shallow_array_adaptor<double>(3, origin));
86 mr = mr + vec_origin;
87 }
88 }
89 if (shift) {
90 ublas::matrix_row<MatrixDouble> mr(coords, nn);
91 VectorAdaptor vec_shift(
92 3, ublas::shallow_array_adaptor<double>(3, shift));
93 mr = mr + vec_shift;
94 }
95 CHKERR moab.create_vertex(&coords(nn, 0), new_verts[nn]);
96 verts_map[conn[nn]] = new_verts[nn];
97 }
98 }
99 EntityHandle ele;
100 CHKERR moab.create_element(MBTRI, new_verts, num_nodes, ele);
101 sUrface.insert(ele);
102 }
103 if (!save_mesh.empty())
104 CHKERR SaveData(m_field.get_moab())(save_mesh, sUrface);
106}
107
113
119
125
131
133 const double rel_tol,
134 const double abs_tol,
135 Tag th,
136 const bool debug) {
137 CoreInterface &m_field = cOre;
138 moab::Interface &moab = m_field.get_moab();
140
141 // Get cutting surface skin
142 Skinner skin(&moab);
143 Range surface_skin;
144 CHKERR skin.find_skin(0, sUrface, false, surface_skin);
145
146 CHKERR snapSurfaceToEdges(surface_skin, fixed_edges, rel_tol, abs_tol, th,
147 debug);
148
150}
151
153 const Range fixed_edges,
154 const double rel_tol,
155 const double abs_tol,
156 Tag th, const bool debug) {
157 CoreInterface &m_field = cOre;
158 moab::Interface &moab = m_field.get_moab();
159 FTensor::Index<'i', 3> i;
161
162 map<EntityHandle, double> map_verts_length;
163
164 for (auto f : surface_edges) {
165 int num_nodes;
166 const EntityHandle *conn_skin;
167 CHKERR moab.get_connectivity(f, conn_skin, num_nodes, true);
168 VectorDouble6 coords_skin(6);
169 if (th)
170 CHKERR moab.tag_get_data(th, conn_skin, num_nodes, &coords_skin[0]);
171 else
172 CHKERR moab.get_coords(conn_skin, num_nodes, &coords_skin[0]);
174 &coords_skin[0], &coords_skin[1], &coords_skin[2]);
176 &coords_skin[3], &coords_skin[4], &coords_skin[5]);
177 FTensor::Tensor1<double, 3> t_skin_delta;
178 t_skin_delta(i) = t_n1(i) - t_n0(i);
179 const double skin_edge_length = sqrt(t_skin_delta(i) * t_skin_delta(i));
180 for (int nn = 0; nn != 2; ++nn) {
181 auto it = map_verts_length.find(conn_skin[nn]);
182 if (it != map_verts_length.end())
183 it->second = std::min(it->second, skin_edge_length);
184 else
185 map_verts_length[conn_skin[nn]] = skin_edge_length;
186 }
187 }
188
189 for (auto m : map_verts_length) {
190
192 if (th)
193 CHKERR moab.tag_get_data(th, &m.first, 1, &t_n(0));
194 else
195 CHKERR moab.get_coords(&m.first, 1, &t_n(0));
196
197 double min_dist = rel_tol * m.second;
198 FTensor::Tensor1<double, 3> t_min_coords;
199 CHKERR cOre.getInterface<Tools>()->findMinDistanceFromTheEdges(
200 &t_n(0), 1, fixed_edges, &min_dist, &t_min_coords(0));
201
202 if (min_dist < rel_tol * m.second || min_dist < abs_tol) {
203 if (debug)
204 MOFEM_LOG("WORLD", Sev::noisy) << "Snap " << min_dist;
205 if (th)
206 CHKERR moab.tag_set_data(th, &m.first, 1, &t_min_coords(0));
207 else
208 CHKERR moab.set_coords(&m.first, 1, &t_min_coords(0));
209 }
210 }
211
213}
214
216 CoreInterface &m_field = cOre;
217 moab::Interface &moab = m_field.get_moab();
219 treeSurfPtr = boost::shared_ptr<OrientedBoxTreeTool>(
220 new OrientedBoxTreeTool(&moab, "ROOTSETSURF", true));
223}
224
227 const double tol_geometry, const double tol_cut_close,
228 Range *fixed_edges, Range *corner_nodes,
229 const bool update_meshsets, const bool debug) {
230 CoreInterface &m_field = cOre;
231 moab::Interface &moab = m_field.get_moab();
233
234 // cut mesh
236 CHKERR projectZeroDistanceEnts(fixed_edges, corner_nodes, tol_geometry,
237 tol_cut_close, QUIET, debug);
240
241 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(constrainSurface,
243 if (fixed_edges)
244 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(*fixed_edges,
245 *fixed_edges);
246 if (corner_nodes)
247 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(*corner_nodes,
248 *corner_nodes);
249 if (update_meshsets)
251 ->updateAllMeshsetsByEntitiesChildren(cut_bit);
253
254 if (debug) {
256 if (fixed_edges)
257 CHKERR SaveData(moab)("out_cut_new_fixed_edges.vtk", *fixed_edges);
258 }
259
261}
262
264 const double tol_trim_close,
265 Range *fixed_edges,
266 Range *corner_nodes,
267 const bool update_meshsets,
268 const bool debug) {
269 CoreInterface &m_field = cOre;
270 moab::Interface &moab = m_field.get_moab();
272
273 // trim mesh
274 CHKERR findEdgesToTrim(fixed_edges, corner_nodes, th, tol_trim_close, debug);
276
277 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(constrainSurface,
279 if (fixed_edges)
280 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(*fixed_edges,
281 *fixed_edges);
282
283 if (corner_nodes)
284 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(*corner_nodes,
285 *corner_nodes);
286
287 if (update_meshsets)
289 ->updateAllMeshsetsByEntitiesChildren(trim_bit);
290
291 // move nodes
293
294 // remove faces
295 CHKERR trimSurface(fixed_edges, corner_nodes, tol_trim_close, th, debug);
296
297 if (debug) {
299 Range bit2_edges;
300 CHKERR cOre.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
301 trim_bit, BitRefLevel().set(), MBEDGE, bit2_edges);
302 CHKERR SaveData(moab)("trim_fixed_edges.vtk",
303 intersect(*fixed_edges, bit2_edges));
304 }
305
307}
308
310 int &first_bit, Tag th, const double tol_geometry,
311 const double tol_cut_close, const double tol_trim_close, Range *fixed_edges,
312 Range *corner_nodes, const bool update_meshsets, const bool debug) {
313 CoreInterface &m_field = cOre;
315
316 std::vector<BitRefLevel> bit_levels;
317
318 auto get_back_bit_levels = [&]() {
319 bit_levels.push_back(BitRefLevel());
320 bit_levels.back().set(first_bit + bit_levels.size() - 1);
321 return bit_levels.back();
322 };
323
324 BitRefLevel cut_bit = get_back_bit_levels();
325
327 tol_geometry, tol_cut_close, fixed_edges, corner_nodes,
328 update_meshsets, debug);
329
330 auto get_min_quality = [&m_field](const BitRefLevel bit, Tag th) {
331 Range tets_level; // test at level
332 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
333 bit, BitRefLevel().set(), MBTET, tets_level);
334 double min_q = 1;
335 CHKERR m_field.getInterface<Tools>()->minTetsQuality(tets_level, min_q, th);
336 return min_q;
337 };
338
339 MOFEM_LOG_C("WORLD", Sev::inform, "Min quality cut %6.4g",
340 get_min_quality(cut_bit, th));
341
342 Range starting_volume = cutNewVolumes;
343 Range starting_volume_nodes;
344 CHKERR m_field.get_moab().get_connectivity(starting_volume,
345 starting_volume_nodes, true);
346 std::vector<double> staring_volume_coords(3 * starting_volume_nodes.size());
347 CHKERR m_field.get_moab().get_coords(starting_volume_nodes,
348 &*staring_volume_coords.begin());
349
350 if (th) {
351 std::vector<double> staring_volume_th_coords(3 *
352 starting_volume_nodes.size());
353 CHKERR m_field.get_moab().tag_get_data(th, starting_volume_nodes,
354 &*staring_volume_th_coords.begin());
355 CHKERR m_field.get_moab().set_coords(starting_volume_nodes,
356 &*staring_volume_th_coords.begin());
357 }
358
359 if (th)
361
362 BitRefLevel trim_bit = get_back_bit_levels();
363
364 CHKERR trimOnly(trim_bit, th, tol_trim_close, fixed_edges, corner_nodes,
365 update_meshsets, debug);
366
367 MOFEM_LOG_C("WORLD", Sev::inform, "Min quality trim %3.2g",
368 get_min_quality(trim_bit, th));
369
370 first_bit += bit_levels.size() - 1;
371
372 if (th)
373 CHKERR m_field.get_moab().set_coords(starting_volume_nodes,
374 &*staring_volume_coords.begin());
375
377}
378
380 int &first_bit, const int fraction_level, Tag th, const double tol_geometry,
381 const double tol_cut_close, const double tol_trim_close, Range &fixed_edges,
382 Range &corner_nodes, const bool update_meshsets, const bool debug) {
383 CoreInterface &m_field = cOre;
385
386 std::vector<BitRefLevel> bit_levels;
387
388 auto get_back_bit_levels = [&]() {
389 bit_levels.push_back(BitRefLevel());
390 bit_levels.back().set(first_bit + bit_levels.size() - 1);
391 return bit_levels.back();
392 };
393
394 if (debug) {
395 CHKERR cOre.getInterface<BitRefManager>()->writeEntitiesNotInDatabase(
396 "ents_not_in_database.vtk", "VTK", "");
397 }
398
399 CHKERR cutAndTrim(first_bit, th, tol_geometry, tol_cut_close, tol_trim_close,
400 &fixed_edges, &corner_nodes, update_meshsets, debug);
401 if (debug)
402 CHKERR cOre.getInterface<BitRefManager>()->writeEntitiesNotInDatabase(
403 "cut_trim_ents_not_in_database.vtk", "VTK", "");
404
405 BitRefLevel bit_level1 = BitRefLevel().set(first_bit - 1);
406 BitRefLevel bit_level2 = get_back_bit_levels();
407 BitRefLevel bit_level3 = get_back_bit_levels();
408
409 CHKERR mergeBadEdges(fraction_level, bit_level1, bit_level2, bit_level3,
410 getNewTrimSurfaces(), fixed_edges, corner_nodes, th,
411 update_meshsets, debug);
412
413 auto get_min_quality = [&m_field, debug](const BitRefLevel bit, Tag th) {
414 Range tets_level; // test at level
415 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
416 bit, BitRefLevel().set(), MBTET, tets_level);
417 double min_q = 1;
418 CHKERR m_field.getInterface<Tools>()->minTetsQuality(tets_level, min_q, th);
419 if (min_q < 0 && debug) {
420 CHKERR m_field.getInterface<Tools>()->writeTetsWithQuality(
421 "negative_tets.vtk", "VTK", "", tets_level, th);
422 }
423 return min_q;
424 };
425
426 MOFEM_LOG_C("WORLD", Sev::inform, "Min quality node merge %6.4g",
427 get_min_quality(bit_level3, th));
428
429 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(constrainSurface,
431 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(fixed_edges,
432 fixed_edges);
433 CHKERR cOre.getInterface<BitRefManager>()->updateRangeByChildren(corner_nodes,
434 corner_nodes);
435
436 first_bit += bit_levels.size() - 1;
437
438 if (debug) {
439 CHKERR cOre.getInterface<BitRefManager>()->writeBitLevelByType(
440 bit_level3, BitRefLevel().set(), MBTET, "out_tets_merged.vtk", "VTK",
441 "");
442 CHKERR cOre.getInterface<BitRefManager>()->writeEntitiesNotInDatabase(
443 "cut_trim_merge_ents_not_in_database.vtk", "VTK", "");
444 CHKERR SaveData(m_field.get_moab())("merged_surface.vtk", mergedSurfaces);
445 }
446
448}
449
451 CoreInterface &m_field = cOre;
452 moab::Interface &moab = m_field.get_moab();
454 Skinner skin(&moab);
455 Range tets_skin;
456 CHKERR skin.find_skin(0, vOlume, false, tets_skin);
457 Range tets_skin_edges;
458 ErrorCode tmp_result;
459 tmp_result = moab.get_adjacencies(tets_skin, 1, false, tets_skin_edges,
460 moab::Interface::UNION);
461
462 if (MB_SUCCESS != tmp_result)
463 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
464 "Duplicated edges: most likely the source of error is comming from "
465 "adding the vertices of the cracking "
466 "volume to a BLOCKSET rather than NODESET (corresponding to the "
467 "input parameter-vertex_block_set)");
468
469 Range surface_skin;
470 CHKERR skin.find_skin(0, sUrface, false, surface_skin);
471 fRont = subtract(surface_skin, tets_skin_edges);
472 if (fRont.empty()) {
473 MOFEM_LOG("WORLD", Sev::warning)
474 << "The front not found. Using edges of the cut surface";
475 CHKERR moab.get_adjacencies(sUrface, 1, true, fRont,
476 moab::Interface::UNION);
477 }
478 if (debug)
479 CHKERR SaveData(moab)("front.vtk", fRont);
481}
482
484 const bool debug) {
485 CoreInterface &m_field = cOre;
486 moab::Interface &moab = m_field.get_moab();
488 auto tools_interface = m_field.getInterface<Tools>();
489
490 auto create_tag = [&](const std::string name, const int dim) {
491 Tag th;
492 rval = moab.tag_get_handle(name.c_str(), th);
493 if (rval == MB_SUCCESS)
494 return th;
495 std::vector<double> def_val(dim, 0);
496 CHKERR moab.tag_get_handle(name.c_str(), dim, MB_TYPE_DOUBLE, th,
497 MB_TAG_CREAT | MB_TAG_SPARSE, &*def_val.begin());
498
499 return th;
500 };
501
502 auto set_vol = [&](const Range &vol_verts, std::vector<double> &coords,
503 std::vector<double> &dist_surface_vec,
504 std::vector<double> &dist_surface_normal_vec,
505 std::vector<double> &dist_surface_signed_dist_vec) {
507
508 coords.resize(3 * vol_verts.size());
509 dist_surface_vec.resize(3 * vol_verts.size());
510 dist_surface_normal_vec.resize(3 * vol_verts.size());
511 dist_surface_signed_dist_vec.resize(vol_verts.size());
512
513 CHKERR moab.get_coords(vol_verts, &*coords.begin());
514 std::srand(0);
515
516 for (auto v : vol_verts) {
517
518 const int index = vol_verts.index(v);
519 auto point_in = getVectorAdaptor(&coords[3 * index], 3);
520 VectorDouble3 point_out(3);
521 EntityHandle facets_out;
522 CHKERR treeSurfPtr->closest_to_location(&point_in[0], rootSetSurf,
523 &point_out[0], facets_out);
524
525 VectorDouble3 delta = point_out - point_in;
526 VectorDouble3 n_first(3);
527 CHKERR tools_interface->getTriNormal(facets_out, &*n_first.begin());
528
529 // Check of only triangle in proximity of point out is one triangle, if is
530 // more than one, use normal of lager triangle to set orientation.
531 auto check_triangle_orientation = [&](auto n) {
532 int num_nodes;
533 const EntityHandle *conn;
534 CHKERR moab.get_connectivity(facets_out, conn, num_nodes, true);
535 MatrixDouble3by3 coords(3, 3);
536 CHKERR moab.get_coords(conn, 3, &coords(0, 0));
537 VectorDouble3 center(3);
538 center.clear();
539 for (auto ii : {0, 1, 2})
540 for (auto jj : {0, 1, 2})
541 center(jj) += coords(ii, jj) / 3;
542
543 std::vector<EntityHandle> triangles_out;
544 std::vector<double> distance_out;
545 auto a_max = norm_2(n);
546 VectorDouble3 ray = -n / a_max;
547 VectorDouble3 pt = center - ray * sqrt(a_max);
548 const double ray_length = 2 * sqrt(a_max);
549
550 CHKERR treeSurfPtr->ray_intersect_triangles(
551 distance_out, triangles_out, rootSetSurf,
552 std::numeric_limits<float>::epsilon(), &pt[0], &ray[0],
553 &ray_length);
554
555 if (triangles_out.size() > 1) {
556 VectorDouble3 nt(3);
557 for (auto t : triangles_out) {
558 CHKERR tools_interface->getTriNormal(t, &*nt.begin());
559 double at = norm_2(nt);
560 if (at > a_max) {
561 n = nt;
562 a_max = at;
563 }
564 }
565 }
566
567 return n;
568 };
569
570 auto n = check_triangle_orientation(n_first);
571 n /= norm_2(n);
572 const double dot = inner_prod(delta, n);
573
574
575 if (std::abs(dot) < std::numeric_limits<double>::epsilon()) {
576 if (std::rand() % 2 == 0)
577 delta += n * std::numeric_limits<double>::epsilon();
578 else
579 delta -= n * std::numeric_limits<double>::epsilon();
580 }
581
582 auto dist_vec = getVectorAdaptor(&dist_surface_vec[3 * index], 3);
583 noalias(dist_vec) = delta;
584
585 auto dist_normal_vec =
586 getVectorAdaptor(&dist_surface_normal_vec[3 * index], 3);
587
588 dist_surface_signed_dist_vec[index] = dot;
589 noalias(dist_normal_vec) = dot * n;
590 }
591
593 };
594
595 std::vector<double> coords;
596 std::vector<double> dist_surface_vec;
597 std::vector<double> dist_surface_normal_vec;
598 std::vector<double> dist_surface_signed_dist_vec;
599 Range vol_verts;
600 CHKERR moab.get_connectivity(vOlume, vol_verts, true);
601
602 CHKERR set_vol(vol_verts, coords, dist_surface_vec, dist_surface_normal_vec,
603 dist_surface_signed_dist_vec);
604
605 CHKERR moab.tag_set_data(create_tag("DIST_SURFACE_VECTOR", 3), vol_verts,
606 &*dist_surface_vec.begin());
607 CHKERR moab.tag_set_data(create_tag("DIST_SURFACE_NORMAL_VECTOR", 3),
608 vol_verts, &*dist_surface_normal_vec.begin());
609 CHKERR moab.tag_set_data(create_tag("DIST_SURFACE_NORMAL_SIGNED", 1),
610 vol_verts, &*dist_surface_signed_dist_vec.begin());
611
613}
614
616 int verb,
617 const bool debug) {
618 CoreInterface &m_field = cOre;
619 moab::Interface &moab = m_field.get_moab();
621
622 auto create_tag = [&](const std::string name, const int dim) {
623 Tag th;
624 rval = moab.tag_get_handle(name.c_str(), th);
625 if (rval == MB_SUCCESS)
626 return th;
627 std::vector<double> def_val(dim, 0);
628 CHKERR moab.tag_get_handle(name.c_str(), dim, MB_TYPE_DOUBLE, th,
629 MB_TAG_CREAT | MB_TAG_SPARSE, &*def_val.begin());
630 return th;
631 };
632
633 Range vol_vertices;
634 CHKERR moab.get_connectivity(vol, vol_vertices, true);
635
636 std::vector<double> min_distances_from_front(vol_vertices.size(), -1);
637 std::vector<double> points_on_edges(3 * vol_vertices.size(), 0);
638 std::vector<EntityHandle> closest_edges(vol_vertices.size(), 0);
639
640 std::vector<double> coords(3 * vol_vertices.size());
641 if (th)
642 CHKERR moab.tag_get_data(th, vol_vertices, &*coords.begin());
643 else
644 CHKERR moab.get_coords(vol_vertices, &*coords.begin());
645
646 CHKERR cOre.getInterface<Tools>()->findMinDistanceFromTheEdges(
647 &*coords.begin(), vol_vertices.size(), fRont,
648 &*min_distances_from_front.begin(), &*points_on_edges.begin(),
649 &*closest_edges.begin());
650
651 if (!points_on_edges.empty()) {
652 for (int i = 0; i != min_distances_from_front.size(); ++i) {
653 Range faces;
654 CHKERR moab.get_adjacencies(&closest_edges[0], 1, 2, false, faces);
655 auto point_in = getVectorAdaptor(&coords[3 * i], 3);
656 auto point_out = getVectorAdaptor(&points_on_edges[3 * i], 3);
657 point_out -= point_in;
658 }
659 }
660
661 auto th_dist_front_vec = create_tag("DIST_FRONT_VECTOR", 3);
662 CHKERR moab.tag_set_data(th_dist_front_vec, vol_vertices,
663 &*points_on_edges.begin());
664
666}
667
669 Tag th, Range &vol_edges, const bool remove_adj_prims_edges, int verb,
670 const bool debug, const std::string edges_file_name) {
671 CoreInterface &m_field = cOre;
672 moab::Interface &moab = m_field.get_moab();
674
675 auto get_tag_dim = [&](auto th) {
676 int dim;
677 CHKERR moab.tag_get_length(th, dim);
678 return dim;
679 };
680 auto dim = get_tag_dim(th);
681
682 auto get_tag_data = [&](auto th, auto conn) {
683 const void *ptr;
684 CHKERR moab.tag_get_by_ptr(th, &conn, 1, &ptr);
685 return getVectorAdaptor(
686 const_cast<double *>(static_cast<const double *>(ptr)), 3);
687 };
688
689 auto get_edge_ray = [&](auto conn) {
690 VectorDouble6 coords(6);
691 CHKERR moab.get_coords(conn, 2, &coords[0]);
692 VectorAdaptor n0 = getVectorAdaptor(&coords[0], 3);
693 VectorAdaptor n1 = getVectorAdaptor(&coords[3], 3);
694 VectorDouble3 ray = n1 - n0;
695 return ray;
696 };
697
698 Range edges;
699 CHKERR moab.get_adjacencies(vOlume, 1, true, edges, moab::Interface::UNION);
700
701 auto remove_prisms_edges = [&](Range &edges) {
703 Range prisms;
704 CHKERR moab.get_adjacencies(edges, 3, false, prisms,
705 moab::Interface::UNION);
706 prisms = prisms.subset_by_type(MBPRISM);
707 Range prisms_verts;
708 CHKERR moab.get_connectivity(prisms, prisms_verts, true);
709 Range prism_edges;
710 CHKERR moab.get_adjacencies(prisms_verts, 1, false, prism_edges,
711 moab::Interface::UNION);
712 edges = subtract(edges, prism_edges);
714 };
715 if (remove_adj_prims_edges)
716 CHKERR remove_prisms_edges(edges);
717
718 std::vector<EntityHandle> cut_edges;
719 cut_edges.reserve(edges.size());
720
721 auto get_cut_edges_vec = [&](auto th, auto &cut_edges, auto e, auto &&conn) {
723
724 auto ray = get_edge_ray(conn);
725 const double length = norm_2(ray);
726 ray /= length;
727 const auto dist0 = get_tag_data(th, conn[0]);
728 const auto dist1 = get_tag_data(th, conn[1]);
729 const double max_dist = std::max(norm_2(dist0), norm_2(dist1));
730 if (max_dist < length) {
731 cut_edges.push_back(e);
732 }
733
735 };
736
737 auto get_cut_edges_signed_dist = [&](auto th, auto &cut_edges, auto e,
738 auto &&conn) {
740 auto get_tag_signed_dist = [&](auto conn) {
741 double dist;
742 CHKERR moab.tag_get_data(th, &conn, 1, &dist);
743 return dist;
744 };
745 const auto dist0 = get_tag_signed_dist(conn[0]);
746 const auto dist1 = get_tag_signed_dist(conn[1]);
747 const double opposite = dist0 * dist1;
748 if (opposite <= 0)
749 cut_edges.push_back(e);
751 };
752
753 auto get_conn = [&](auto e) {
754 int num_nodes;
755 const EntityHandle *conn;
756 CHKERR moab.get_connectivity(e, conn, num_nodes, true);
757 return conn;
758 };
759
760 if (dim == 3)
761 for (auto e : edges)
762 CHKERR get_cut_edges_vec(th, cut_edges, e, get_conn(e));
763 else
764 for (auto e : edges)
765 CHKERR get_cut_edges_signed_dist(th, cut_edges, e, get_conn(e));
766
767 CHKERR moab.get_adjacencies(&*cut_edges.begin(), cut_edges.size(), 3, false,
768 vol_edges, moab::Interface::UNION);
769 vol_edges = intersect(vol_edges, vOlume);
770
771 auto convert_to_ranger = [](auto &v) {
772 Range r;
773 r.insert_list(v.begin(), v.end());
774 return r;
775 };
776
777 if (debug && !edges_file_name.empty())
778 CHKERR SaveData(m_field.get_moab())(edges_file_name,
779 convert_to_ranger(cut_edges));
780
782}
783
785 const bool debug) {
786 CoreInterface &m_field = cOre;
787 moab::Interface &moab = m_field.get_moab();
789
790 CHKERR createFrontLevelSets(vOlume, nullptr, verb, debug);
791 Tag th_dist_front_vec;
792 CHKERR moab.tag_get_handle("DIST_FRONT_VECTOR", th_dist_front_vec);
793 CHKERR findLevelSetVolumes(th_dist_front_vec, cutFrontVolumes, true, verb,
794 debug, "cutFrontEdges.vtk");
795
797
798 Tag th_dist_surface_vec;
799 CHKERR moab.tag_get_handle("DIST_SURFACE_VECTOR", th_dist_surface_vec);
800 cutSurfaceVolumes.clear();
801 CHKERR findLevelSetVolumes(th_dist_surface_vec, cutSurfaceVolumes, true, verb,
802 debug, "cutSurfaceEdges.vtk");
803
804 if (debug)
805 CHKERR SaveData(m_field.get_moab())("level_sets.vtk", vOlume);
806 if (debug)
807 CHKERR SaveData(m_field.get_moab())("cutSurfaceVolumes.vtk",
809 if (debug)
810 CHKERR SaveData(m_field.get_moab())("cutFrontVolumes.vtk", cutFrontVolumes);
811
813}
814
816 const int surf_levels,
817 const int front_levels,
818 Range *fixed_edges, int verb,
819 const bool debug) {
820 CoreInterface &m_field = cOre;
821 moab::Interface &moab = m_field.get_moab();
822 MeshRefinement *refiner;
823 BitRefManager *bit_ref_manager;
825 CHKERR m_field.getInterface(refiner);
826 CHKERR m_field.getInterface(bit_ref_manager);
827
828 auto add_bit = [&](const int bit) {
830 CHKERR bit_ref_manager->addBitRefLevel(vOlume, BitRefLevel().set(bit),
831 verb);
832 Range adj_ents;
833 for (auto d : {2, 1, 0}) {
834 CHKERR moab.get_adjacencies(vOlume, d, true, adj_ents,
835 moab::Interface::UNION);
836 CHKERR bit_ref_manager->addBitRefLevel(adj_ents, BitRefLevel().set(bit),
837 verb);
838 }
840 };
841 CHKERR add_bit(init_bit_level);
842
843 auto update_range = [&](Range *r_ptr) {
845 if (r_ptr) {
846 Range childs;
847 CHKERR bit_ref_manager->updateRangeByChildren(*r_ptr, childs);
848 r_ptr->merge(childs);
849 }
851 };
852
853 auto refine = [&](const BitRefLevel bit, const Range tets) {
855 Range verts;
856 CHKERR moab.get_connectivity(tets, verts, true);
857 Range ref_edges;
858 CHKERR moab.get_adjacencies(verts, 1, true, ref_edges,
859 moab::Interface::UNION);
860
861 CHKERR refiner->addVerticesInTheMiddleOfEdges(ref_edges, bit);
862 CHKERR refiner->refineTets(vOlume, bit, verb);
863
864 CHKERR update_range(fixed_edges);
865 CHKERR update_range(&vOlume);
867 ->updateAllMeshsetsByEntitiesChildren(bit);
868
869 Range bit_tets;
870 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
871 bit, BitRefLevel().set(), MBTET, bit_tets);
872 vOlume = intersect(vOlume, bit_tets);
873
874 Range bit_edges;
875 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
876 bit, BitRefLevel().set(), MBEDGE, bit_edges);
877 if (fixed_edges)
878 *fixed_edges = intersect(*fixed_edges, bit_edges);
879
881 };
882
883 for (int ll = init_bit_level; ll != init_bit_level + surf_levels; ++ll) {
885 CHKERR refine(BitRefLevel().set(ll + 1),
887 }
888
889 for (int ll = init_bit_level + surf_levels;
890 ll != init_bit_level + surf_levels + front_levels; ++ll) {
892 CHKERR refine(BitRefLevel().set(ll + 1), cutFrontVolumes);
893 }
894
896
897 if (debug)
898 CHKERR SaveData(m_field.get_moab())("refinedVolume.vtk", vOlume);
899
901}
902
904 const bool debug) {
905 CoreInterface &m_field = cOre;
906 moab::Interface &moab = m_field.get_moab();
908
909 edgesToCut.clear();
910 cutEdges.clear();
911
912 Tag th_signed_dist;
913 CHKERR moab.tag_get_handle("DIST_SURFACE_NORMAL_SIGNED", th_signed_dist);
914
915 auto get_tag_edge_dist = [&](auto th, auto conn) {
916 std::array<double, 2> r;
917 CHKERR moab.tag_get_data(th, conn, 2, r.data());
918 return r;
919 };
920
921 auto get_conn = [&](const auto e) {
922 int num_nodes;
923 const EntityHandle *conn;
924 CHKERR moab.get_connectivity(e, conn, num_nodes, true);
925 return conn;
926 };
927
928 auto get_adj = [&moab](const Range r, const int dim) {
929 Range a;
930 if (dim)
931 CHKERR moab.get_adjacencies(r, dim, false, a, moab::Interface::UNION);
932 else
933 CHKERR moab.get_connectivity(r, a, true);
934 return a;
935 };
936
937 auto vol_edges = get_adj(vol, 1);
938
939 aveLength = 0;
940 maxLength = 0;
941 int nb_ave_length = 0;
942 for (auto e : vol_edges) {
943
944 auto conn = get_conn(e);
945 auto dist = get_tag_edge_dist(th_signed_dist, conn);
946 const auto dist_max = std::max(dist[0], dist[1]);
947 const auto dist_min = std::min(dist[0], dist[1]);
948 const auto dot = dist[0] * dist[1];
949
950 VectorDouble6 coords(6);
951 CHKERR moab.get_coords(conn, 2, &coords[0]);
952 VectorAdaptor n0 = getVectorAdaptor(&coords[0], 3);
953 VectorAdaptor n1 = getVectorAdaptor(&coords[3], 3);
954 VectorDouble3 ray = n1 - n0;
955 const double ray_length = norm_2(ray);
956 ray /= ray_length;
957
958 if (
959
960 (dot < 0 && dist_max > std::numeric_limits<double>::epsilon()) ||
961 (std::abs(dist_min) < std::numeric_limits<double>::epsilon() &&
962 dist_max > std::numeric_limits<double>::epsilon())
963
964 ) {
965
966 // Edges is on two sides of the surface
967 const double s = dist[0] / (dist[0] - dist[1]);
968 const double dist = s * ray_length;
969
970 auto add_edge = [&](auto dist) {
971 edgesToCut[e].dIst = dist;
972 edgesToCut[e].lEngth = ray_length;
973 edgesToCut[e].unitRayDir = ray;
974 edgesToCut[e].rayPoint = n0;
975 cutEdges.insert(e);
976
977 aveLength += norm_2(ray);
978 maxLength = fmax(maxLength, norm_2(ray));
979 ++nb_ave_length;
980 };
981
982 add_edge(dist);
983 }
984 }
985 aveLength /= nb_ave_length;
986
987 if (debug)
988 CHKERR SaveData(m_field.get_moab())("cut_edges.vtk", cutEdges);
989
990 if (debug)
991 CHKERR SaveData(m_field.get_moab())("cut_edges_vol.vtk",
992 get_adj(cutEdges, 3));
993
995}
996
998 Range *fixed_edges, Range *corner_nodes, const double geometry_tol,
999 const double close_tol, const int verb, const bool debug) {
1000 CoreInterface &m_field = cOre;
1001 moab::Interface &moab = m_field.get_moab();
1002 Skinner skin(&moab);
1004
1005 auto get_ent_adj = [&moab](const EntityHandle v, const int dim) {
1006 Range a;
1007 if (dim)
1008 CHKERR moab.get_adjacencies(&v, 1, dim, true, a);
1009 return a;
1010 };
1011
1012 auto get_adj = [&moab](const Range r, const int dim) {
1013 Range a;
1014 if (dim)
1015 CHKERR moab.get_adjacencies(r, dim, false, a, moab::Interface::UNION);
1016 else
1017 CHKERR moab.get_connectivity(r, a, true);
1018 return a;
1019 };
1020
1021 auto get_skin = [&skin](const auto r) {
1022 Range s;
1023 CHKERR skin.find_skin(0, r, false, s);
1024 return s;
1025 };
1026
1027
1028 auto get_range = [](std::vector<EntityHandle> v) {
1029 Range r;
1030 r.insert_list(v.begin(), v.end());
1031 return r;
1032 };
1033
1034 auto get_coords = [&](const auto v) {
1035 VectorDouble3 coords(3);
1036 CHKERR moab.get_coords(&v, 1, &coords[0]);
1037 return coords;
1038 };
1039
1040 Tag th_normal_dist;
1041 CHKERR moab.tag_get_handle("DIST_SURFACE_NORMAL_VECTOR", th_normal_dist);
1042
1043 auto get_normal_dits = [&](const EntityHandle v) {
1044 VectorDouble3 dist_vec(3);
1045 CHKERR moab.tag_get_data(th_normal_dist, &v, 1, &*dist_vec.begin());
1046 return dist_vec;
1047 };
1048
1049 auto get_prj_point = [&](const EntityHandle v, const Range edges,
1050 const double geometry_tol) {
1051 auto get_tuple = [](const EntityHandle e, const double dist,
1052 const double l) { return std::make_tuple(e, dist, l); };
1053
1054 std::tuple<EntityHandle, double, double> min_tuple{0, 0, 0};
1055
1056 for (auto e : edges) {
1057 auto eit = edgesToCut.find(e);
1058 if (eit != edgesToCut.end()) {
1059
1060 auto get_dist = [&](auto eit) {
1061 int num_nodes;
1062 const EntityHandle *conn;
1063 CHKERR moab.get_connectivity(eit->first, conn, num_nodes, true);
1064 if (conn[0] == v)
1065 return eit->second.dIst;
1066 else if (conn[1] == v)
1067 return (eit->second.lEngth - eit->second.dIst);
1068 else
1069 THROW_MESSAGE("Data inconsistency");
1070 };
1071
1072 const auto d = get_dist(eit);
1073 if (std::get<0>(min_tuple) == 0) {
1074 min_tuple = get_tuple(e, d, eit->second.lEngth);
1075 } else if (std::get<1>(min_tuple) > d) {
1076 min_tuple = get_tuple(e, d, eit->second.lEngth);
1077 }
1078 }
1079 }
1080
1081 const auto geom_dist_vec = get_normal_dits(v);
1082 const auto geom_tol = norm_2(geom_dist_vec);
1083 const auto l = std::get<2>(min_tuple);
1084
1085 if (geom_tol < l * geometry_tol) {
1086 return std::make_pair(get_coords(v), l);
1087
1088 } else {
1089 const auto &d = edgesToCut.at(std::get<0>(min_tuple));
1090 return std::make_pair(VectorDouble3(d.rayPoint + d.dIst * d.unitRayDir),
1091 l);
1092 }
1093 };
1094
1095 auto get_in_range = [](auto v, auto &r) { return (r.find(v) != r.end()); };
1096
1097 auto project_nodes = [&](auto nodes_to_check) {
1098 if (debug)
1099 CHKERR SaveData(moab)("nodes_to_check.vtk", nodes_to_check);
1100
1101 auto get_fix_e = [](auto fixed_edges) {
1102 if (fixed_edges)
1103 return *fixed_edges;
1104 else
1105 return Range();
1106 };
1107
1108 const auto fix_e = get_fix_e(fixed_edges);
1109 const auto skin_e = get_adj(unite(get_skin(vOlume), constrainSurface), 1);
1110 const auto cut_fix = intersect(cutEdges, fix_e);
1111 const auto cut_skin = intersect(cutEdges, skin_e);
1112
1113 Range corner_n;
1114 if (corner_nodes)
1115 corner_n = intersect(*corner_nodes, nodes_to_check);
1116 auto fixe_n = intersect(get_adj(fix_e, 0), nodes_to_check);
1117 auto skin_n = intersect(get_adj(skin_e, 0), nodes_to_check);
1118
1119 std::vector<std::pair<EntityHandle, TreeData>> vertices_on_cut_edges;
1120 vertices_on_cut_edges.reserve(nodes_to_check.size());
1121
1122 auto add_zero_vertex = [&](const EntityHandle v, auto &&new_pos) {
1123 auto coords = get_coords(v);
1124 auto ray = new_pos - coords;
1125 auto l = norm_2(ray);
1126 TreeData data;
1127 data.dIst = l;
1128 if (l > std::numeric_limits<double>::epsilon()) {
1129 data.unitRayDir = ray / l;
1130 } else {
1131 data.dIst = 0;
1132 data.unitRayDir = ray;
1133 }
1134 data.rayPoint = coords;
1135 return std::make_pair(v, data);
1136 };
1137
1138 auto intersect_v = [&](const auto v, const auto r) {
1139 return intersect(r, get_ent_adj(v, 1));
1140 };
1141
1142 for (auto v : nodes_to_check) {
1143
1144 const auto e = intersect_v(v, cutEdges);
1145 if (!e.empty()) {
1146
1147 if (get_in_range(v, corner_n)) {
1148 auto p = get_prj_point(v, e, 0);
1149 if (norm_2(get_coords(v) - p.first) < close_tol * p.second) {
1150 vertices_on_cut_edges.push_back(add_zero_vertex(v, get_coords(v)));
1151 continue;
1152 }
1153
1154 } else if (get_in_range(v, fixe_n)) {
1155 const auto b = intersect_v(v, cut_fix);
1156 if (!b.empty()) {
1157 auto p = get_prj_point(v, b, geometry_tol);
1158 if (norm_2(get_coords(v) - p.first) < close_tol * p.second) {
1159 vertices_on_cut_edges.push_back(add_zero_vertex(v, p.first));
1160 continue;
1161 }
1162
1163 } else if (get_in_range(v, skin_n)) {
1164 const auto b = intersect_v(v, cut_skin);
1165 if (!b.empty()) {
1166 auto p = get_prj_point(v, b, geometry_tol);
1167 if (norm_2(get_coords(v) - p.first) < close_tol * p.second) {
1168 vertices_on_cut_edges.push_back(add_zero_vertex(v, p.first));
1169 continue;
1170 }
1171 }
1172
1173 }
1174 }
1175
1176 auto p = get_prj_point(v, e, geometry_tol);
1177 if (norm_2(get_coords(v) - p.first) < close_tol * p.second) {
1178 if (get_in_range(v, fixe_n) || get_in_range(v, skin_n))
1179 vertices_on_cut_edges.push_back(add_zero_vertex(v, get_coords(v)));
1180 else
1181 vertices_on_cut_edges.push_back(add_zero_vertex(v, p.first));
1182 }
1183 }
1184
1185 }
1186
1187 auto get_distances = [&](auto &data) {
1188 std::map<EntityHandle, double> dist_map;
1189 if (!data.empty()) {
1190 Tag th;
1191 CHKERR moab.tag_get_handle("DIST_SURFACE_NORMAL_SIGNED", th);
1192
1193 std::vector<EntityHandle> verts;
1194 verts.reserve(data.size());
1195 for (auto d : data)
1196 verts.emplace_back(d.first);
1197 std::vector<double> distances(verts.size());
1198 CHKERR moab.tag_get_data(th, &*verts.begin(), verts.size(),
1199 &*distances.begin());
1200 for (size_t i = 0; i != verts.size(); ++i)
1201 dist_map[verts[i]] = distances[i];
1202 }
1203 return dist_map;
1204 };
1205
1206 auto dist_map = get_distances(vertices_on_cut_edges);
1207 if (!dist_map.empty()) {
1208 auto cmp = [&dist_map](const auto &a, const auto &b) {
1209 return dist_map[a.first] < dist_map[b.first];
1210 };
1211 std::sort(vertices_on_cut_edges.begin(), vertices_on_cut_edges.end(),
1212 cmp);
1213 }
1214
1215 return vertices_on_cut_edges;
1216 };
1217
1218 auto get_min_quality =
1219 [&](const Range &adj_tets,
1220 const map<EntityHandle, TreeData> &vertices_on_cut_edges,
1221 const std::pair<EntityHandle, TreeData> *test_data_ptr) {
1222 double min_q = 1;
1223 for (auto t : adj_tets) {
1224 int num_nodes;
1225 const EntityHandle *conn;
1226 CHKERR m_field.get_moab().get_connectivity(t, conn, num_nodes, true);
1227 VectorDouble12 coords(12);
1228 CHKERR moab.get_coords(conn, num_nodes, &coords[0]);
1229
1230 auto set_new_coord = [](auto d) {
1231 return d.rayPoint + d.dIst * d.unitRayDir;
1232 };
1233
1234 for (auto n : {0, 1, 2, 3}) {
1235 auto n_coords = getVectorAdaptor(&coords[3 * n], 3);
1236 if (test_data_ptr) {
1237 if (test_data_ptr->first == conn[n]) {
1238 noalias(n_coords) = set_new_coord(test_data_ptr->second);
1239 }
1240 }
1241 auto mit = vertices_on_cut_edges.find(conn[n]);
1242 if (mit != vertices_on_cut_edges.end()) {
1243 noalias(n_coords) = set_new_coord(mit->second);
1244 }
1245 }
1246 min_q = std::min(min_q, Tools::volumeLengthQuality(&coords[0]));
1247 }
1248 return min_q;
1249 };
1250
1251 auto get_zero_distance_verts = [&](const auto &&vertices_on_cut_edges) {
1252 std::vector<EntityHandle> zero_dist_vec;
1253 zero_dist_vec.reserve(vertices_on_cut_edges.size());
1254 for (auto t : vertices_on_cut_edges) {
1255
1256 auto adj_tet = intersect(vOlume, get_ent_adj(t.first, 3));
1257 const double q0 = get_min_quality(adj_tet, verticesOnCutEdges, nullptr);
1258 const double q = get_min_quality(adj_tet, verticesOnCutEdges, &t);
1259 if (std::isnormal(q) && (q / q0) > projectEntitiesQualityTrashold) {
1260 verticesOnCutEdges[t.first] = t.second;
1261 zero_dist_vec.push_back(t.first);
1262 }
1263 }
1264 return zero_dist_vec;
1265 };
1266
1267 auto vol_cut_ents = intersect(vOlume, get_adj(cutEdges, 3));
1268
1269 auto get_zero_distant_ents = [&](auto bridge_ents, const int dim,
1270 const int bridge_dim) {
1271 auto ents =
1272 intersect(get_adj(bridge_ents, dim), get_adj(vol_cut_ents, dim));
1273 auto ents_to_remove =
1274 get_adj(subtract(get_adj(ents, bridge_dim), bridge_ents), dim);
1275 return subtract(ents, ents_to_remove);
1276 };
1277
1278 verticesOnCutEdges.clear();
1279 zeroDistanceVerts.clear();
1280 zeroDistanceVerts.merge(
1281 get_range(get_zero_distance_verts(project_nodes(get_adj(cutEdges, 0)))));
1282 zeroDistanceEnts = subtract(get_zero_distant_ents(zeroDistanceVerts, 2, 0),
1283 get_skin(vOlume));
1284
1285 if (debug)
1286 CHKERR SaveData(moab)("zero_distance_verts.vtk", zeroDistanceVerts);
1287 if (debug)
1288 CHKERR SaveData(moab)("zero_distance_ents.vtk", zeroDistanceEnts);
1289
1290 for (auto f : zeroDistanceVerts) {
1291 auto adj_edges = get_ent_adj(f, 1);
1292 for (auto e : adj_edges) {
1293 cutEdges.erase(e);
1294 edgesToCut.erase(e);
1295 }
1296 }
1297
1298 if (debug)
1299 CHKERR SaveData(moab)("cut_edges_to_cut.vtk", cutEdges);
1300
1302}
1303
1305 Range &cut_vols,
1306 Range &cut_surf,
1307 Range &cut_verts,
1308 const bool debug) {
1309 CoreInterface &m_field = cOre;
1310 moab::Interface &moab = m_field.get_moab();
1311 MeshRefinement *refiner;
1312 auto ref_ents_ptr = m_field.get_ref_ents();
1314
1315 if (cutEdges.size() != edgesToCut.size())
1316 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Data inconsistency");
1317
1318 auto refine_mesh = [&]() {
1320 CHKERR m_field.getInterface(refiner);
1322 CHKERR refiner->refineTets(vOlume, bit, QUIET, debug);
1324 };
1325
1326 CHKERR refine_mesh();
1327
1328 if (debug)
1329 if (cutEdges.size() != edgesToCut.size())
1330 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Data inconsistency");
1331
1332 cut_vols.clear();
1333 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
1334 bit, BitRefLevel().set(), MBTET, cut_vols);
1335 cut_surf.clear();
1336 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
1337 bit, bit, MBTRI, cut_surf);
1338
1339 if (debug)
1340 CHKERR SaveData(moab)("cut_surf_from_bit.vtk", cut_surf);
1341
1342 // Find new vertices on cut edges
1343 cut_verts.clear();
1344 cut_verts.merge(zeroDistanceVerts);
1345
1346 for (auto &m : edgesToCut) {
1347 auto vit = ref_ents_ptr->get<Composite_ParentEnt_And_EntType_mi_tag>().find(
1348 boost::make_tuple(MBVERTEX, m.first));
1349 if (vit ==
1350 ref_ents_ptr->get<Composite_ParentEnt_And_EntType_mi_tag>().end()) {
1351 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1352 "No vertex on cut edges, that make no sense");
1353 }
1354 const boost::shared_ptr<RefEntity> &ref_ent = *vit;
1355 if ((ref_ent->getBitRefLevel() & bit).any()) {
1356 EntityHandle vert = ref_ent->getEnt();
1357 cut_verts.insert(vert);
1358 verticesOnCutEdges[vert] = m.second;
1359 } else {
1360 SETERRQ(
1361 PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1362 "Vertex has wrong bit ref level %s",
1363 boost::lexical_cast<std::string>(ref_ent->getBitRefLevel()).c_str());
1364 }
1365 }
1366
1367 // Add zero distance entities faces
1368 Range tets_skin;
1369 Skinner skin(&moab);
1370 CHKERR skin.find_skin(0, cut_vols, false, tets_skin);
1371 tets_skin.merge(constrainSurface);
1372
1373 // At that point cut_surf has all newly created faces, now take all
1374 // nodes on those faces and subtract nodes on cut edges. Faces adjacent to
1375 // nodes which left are not part of surface.
1376 cut_surf.merge(zeroDistanceEnts.subset_by_type(MBTRI));
1377 Range diff_verts;
1378 CHKERR moab.get_connectivity(cut_surf, diff_verts, true);
1379 diff_verts = subtract(diff_verts, cut_verts);
1380
1381 Range subtract_faces;
1382 CHKERR moab.get_adjacencies(diff_verts, 2, false, subtract_faces,
1383 moab::Interface::UNION);
1384 cut_surf = subtract(cut_surf, unite(subtract_faces, tets_skin));
1385
1386 cut_verts.clear();
1387 CHKERR moab.get_connectivity(cut_surf, cut_verts, true);
1388
1389 // Check non-manifolds
1390 auto check_for_non_manifold = [&]() {
1392 Range surf_edges;
1393 CHKERR moab.get_adjacencies(cut_surf, 1, false, surf_edges,
1394 moab::Interface::UNION);
1395 for (auto e : surf_edges) {
1396
1397 Range faces;
1398 CHKERR moab.get_adjacencies(&e, 1, 2, false, faces);
1399 faces = intersect(faces, cut_surf);
1400 if (faces.size() > 2) {
1401
1402 bool resolved = false;
1403
1404 // Check for haning node
1405 Range nodes;
1406 CHKERR moab.get_connectivity(faces, nodes, true);
1407 for (auto n : nodes) {
1408 Range adj_faces;
1409 CHKERR moab.get_adjacencies(&n, 1, 2, false, adj_faces);
1410 adj_faces = intersect(adj_faces, cut_surf);
1411 if (adj_faces.size() == 1) {
1412 cut_surf.erase(adj_faces[0]);
1413 resolved = true;
1414 }
1415 }
1416
1417 // Check for two edges minfold
1418 Range adj_edges;
1419 CHKERR moab.get_adjacencies(faces, 1, false, adj_edges,
1420 moab::Interface::UNION);
1421 adj_edges = intersect(adj_edges, surf_edges);
1422 adj_edges.erase(e);
1423 for (auto other_e : adj_edges) {
1424 Range other_faces;
1425 CHKERR moab.get_adjacencies(&other_e, 1, 2, false, other_faces);
1426 other_faces = intersect(other_faces, cut_surf);
1427 if (other_faces.size() > 2) {
1428 other_faces = intersect(other_faces, faces);
1429 cut_surf = subtract(cut_surf, other_faces);
1430 resolved = true;
1431 }
1432 }
1433
1434 if (!resolved && !debug)
1435 MOFEM_LOG("WORLD", Sev::warning) << "Non-mainfold surface after cut.";
1436
1437 cut_verts.clear();
1438 CHKERR moab.get_connectivity(cut_surf, cut_verts, true);
1439 }
1440 }
1442 };
1443
1444 CHKERR check_for_non_manifold();
1445
1446 if (debug)
1447 CHKERR SaveData(moab)("cut_surf.vtk", cut_surf);
1448
1450}
1451
1453 CoreInterface &m_field = cOre;
1454 moab::Interface &moab = m_field.get_moab();
1456
1457 // Range out_side_vertices;
1458 for (auto m : verticesOnCutEdges) {
1459 double dist = m.second.dIst;
1460 VectorDouble3 new_coors = m.second.rayPoint + dist * m.second.unitRayDir;
1461 if (th)
1462 CHKERR moab.tag_set_data(th, &m.first, 1, &new_coors[0]);
1463 else
1464 CHKERR moab.set_coords(&m.first, 1, &new_coors[0]);
1465 }
1466
1468}
1469
1471 CoreInterface &m_field = cOre;
1472 moab::Interface &moab = m_field.get_moab();
1474 for (auto &v : verticesOnTrimEdges) {
1475 double dist = v.second.dIst;
1476 VectorDouble3 new_coors = v.second.rayPoint + dist * v.second.unitRayDir;
1477 if (th)
1478 CHKERR moab.tag_set_data(th, &v.first, 1, &new_coors[0]);
1479 else
1480 CHKERR moab.set_coords(&v.first, 1, &new_coors[0]);
1481 }
1483}
1484
1486 Range *corner_nodes, Tag th,
1487 const double tol_trim_close,
1488 const bool debug) {
1489 CoreInterface &m_field = cOre;
1490 moab::Interface &moab = m_field.get_moab();
1492
1493 // takes body skin
1494 Skinner skin(&moab);
1495 Range tets_skin;
1496 CHKERR skin.find_skin(0, cutNewVolumes, false, tets_skin);
1497 tets_skin.merge(constrainSurface);
1498
1499 // vertices on the skin
1500 Range tets_skin_verts;
1501 CHKERR moab.get_connectivity(tets_skin, tets_skin_verts, true);
1502 // edges on the skin
1503 Range tets_skin_edges;
1504 CHKERR moab.get_adjacencies(tets_skin, 1, false, tets_skin_edges,
1505 moab::Interface::UNION);
1506
1507 // get edges on new surface
1508 Range cut_surface_edges;
1509 CHKERR moab.get_adjacencies(cutNewSurfaces, 1, false, cut_surface_edges,
1510 moab::Interface::UNION);
1511 Range cut_surface_verts;
1512 CHKERR moab.get_connectivity(cutNewSurfaces, cut_surface_verts, true);
1513
1514 Range corners;
1515 if (corner_nodes)
1516 corners = *corner_nodes;
1517
1518 Range fix_edges;
1519 if (fixed_edges)
1520 fix_edges = *fixed_edges;
1521
1522 Range fixed_edges_verts;
1523 if (fixed_edges)
1524 CHKERR moab.get_connectivity(*fixed_edges, fixed_edges_verts, true);
1525
1526 Range surface_skin;
1527 if (fRont.empty())
1528 CHKERR skin.find_skin(0, sUrface, false, surface_skin);
1529 else
1530 surface_skin = fRont;
1531
1532 auto get_point_coords = [&](EntityHandle v) {
1533 VectorDouble3 point(3);
1534 if (th)
1535 CHKERR moab.tag_get_data(th, &v, 1, &point[0]);
1536 else
1537 CHKERR moab.get_coords(&v, 1, &point[0]);
1538 return point;
1539 };
1540
1541 struct VertMap {
1542 const double d;
1543 const EntityHandle v;
1544 const EntityHandle e;
1545 VertMap(const double d, const EntityHandle v, const EntityHandle e)
1546 : d(d), v(v), e(e) {}
1547 };
1548
1549 typedef multi_index_container<
1550 VertMap,
1551 indexed_by<
1552 ordered_non_unique<member<VertMap, const double, &VertMap::d>>,
1553 ordered_non_unique<member<VertMap, const EntityHandle, &VertMap::v>>,
1554 ordered_non_unique<member<VertMap, const EntityHandle, &VertMap::e>>
1555
1556 >>
1557 VertMapMultiIndex;
1558
1559 VertMapMultiIndex verts_map;
1560
1561 auto add_vert = [&](const EntityHandle v, const EntityHandle e,
1562 const double dist) {
1563 verts_map.insert(VertMap(dist, v, e));
1564 };
1565
1566 // clear data ranges
1567 trimEdges.clear();
1568 edgesToTrim.clear();
1569 verticesOnTrimEdges.clear();
1570 trimNewVertices.clear();
1571
1572 auto cut_this_edge = [&](const EntityHandle e, const double length, auto &ray,
1573 auto &ray_point) {
1574 trimEdges.insert(e);
1575 edgesToTrim[e].dIst = 1;
1576 edgesToTrim[e].lEngth = length;
1577 edgesToTrim[e].unitRayDir.resize(3, false);
1578 edgesToTrim[e].rayPoint.resize(3, false);
1579 for (int ii = 0; ii != 3; ++ii)
1580 edgesToTrim[e].unitRayDir[ii] = ray(ii);
1581 for (int ii = 0; ii != 3; ++ii)
1582 edgesToTrim[e].rayPoint[ii] = ray_point(ii);
1583 };
1584
1585 FTensor::Index<'i', 3> i;
1586 FTensor::Index<'j', 3> j;
1587 FTensor::Index<'k', 3> k;
1588
1589 int num_nodes;
1590
1591 auto get_tag_data = [&](auto th, auto conn) {
1593 CHKERR moab.tag_get_data(th, &conn, 1, &t(0));
1594 return t;
1595 };
1596
1597 double max_edge_length = 0;
1598
1599 /// Project front entities on on the cut surface plane
1600 if (!fRont.empty()) {
1601 // Calculate distances
1602 treeSurfPtr = boost::shared_ptr<OrientedBoxTreeTool>(
1603 new OrientedBoxTreeTool(&moab, "ROOTSETSURF", true));
1605
1606 for (auto s : surface_skin) {
1607
1608 auto project_on_surface = [&](auto &t) {
1610
1611 EntityHandle facets_out;
1612 CHKERR treeSurfPtr->closest_to_location(&t(0), rootSetSurf, &t_p(0),
1613 facets_out);
1614
1616 CHKERR m_field.getInterface<Tools>()->getTriNormal(facets_out,
1617 &t_normal(0));
1618 t_normal(i) /= sqrt(t_normal(i) * t_normal(i));
1619 const double dot = t_normal(i) * (t_p(i) - t(i));
1620 t_p(i) = t(i) + dot * t_normal(i);
1621
1622 return t_p;
1623 };
1624
1625 const EntityHandle *conn;
1626 CHKERR moab.get_connectivity(s, conn, num_nodes, true);
1627
1628 VectorDouble6 coords_front(6);
1629 CHKERR moab.get_coords(conn, num_nodes, &coords_front[0]);
1630
1631 FTensor::Tensor1<double *, 3> t_f0(&coords_front[0], &coords_front[1],
1632 &coords_front[2]);
1633 FTensor::Tensor1<double *, 3> t_f1(&coords_front[3], &coords_front[4],
1634 &coords_front[5]);
1635
1636 auto t_p_f0 = project_on_surface(t_f0);
1637 auto t_p_f1 = project_on_surface(t_f1);
1638
1639 CHKERR moab.set_coords(&conn[0], 1, &t_p_f0(0));
1640 CHKERR moab.set_coords(&conn[1], 1, &t_p_f1(0));
1641 }
1642 }
1643
1644 if (debug)
1645 CHKERR SaveData(moab)("surface_skin_to_trim.vtk", surface_skin);
1646
1648 Tag th_dist_front_vec;
1649 CHKERR moab.tag_get_handle("DIST_FRONT_VECTOR", th_dist_front_vec);
1650
1651 // iterate over edges on cut surface
1652 for (auto e : cut_surface_edges) {
1653
1654 auto get_conn_edge = [&moab](auto e) {
1655 // Get edge connectivity and coords
1656 const EntityHandle *conn_edge;
1657 int num_nodes;
1658 CHKERR moab.get_connectivity(e, conn_edge, num_nodes, true);
1659 return conn_edge;
1660 };
1661
1662 auto get_coords_edge = [&moab](auto conn_edge) {
1663 std::array<double, 6> coords_edge;
1664 CHKERR moab.get_coords(conn_edge, 2, coords_edge.data());
1665 return coords_edge;
1666 };
1667
1668 auto get_ftensor_coords = [](const double *ptr) {
1669 return FTensor::Tensor1<double, 3>{ptr[0], ptr[1], ptr[2]};
1670 };
1671
1672 auto conn_edge = get_conn_edge(e);
1673 auto coords_edge = get_coords_edge(conn_edge);
1674
1675 auto t_dist0 = get_tag_data(th_dist_front_vec, conn_edge[0]);
1676 auto t_dist1 = get_tag_data(th_dist_front_vec, conn_edge[1]);
1677 auto t_e0 = get_ftensor_coords(&coords_edge[0]);
1678 auto t_e1 = get_ftensor_coords(&coords_edge[3]);
1679
1680 FTensor::Tensor1<double, 3> t_edge_delta;
1681 t_edge_delta(i) = t_e1(i) - t_e0(i);
1682 const double edge_length2 = t_edge_delta(i) * t_edge_delta(i);
1683 const double edge_length = sqrt(edge_length2);
1684 if (edge_length == 0)
1685 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Zero edge length");
1686
1687 auto add_edge = [&](auto t_cut) {
1689 if (t_cut < 0.5) {
1690 t_ray(i) = t_cut * t_edge_delta(i);
1691 add_vert(conn_edge[0], e, fabs(t_cut));
1692 add_vert(conn_edge[1], e, fabs(t_cut - 1));
1693 cut_this_edge(e, edge_length, t_ray, t_e0);
1694 } else {
1695 FTensor::Tensor1<double, 3> t_edge_point;
1696 t_edge_point(i) = t_e0(i) + t_cut * t_edge_delta(i);
1697 t_ray(i) = t_edge_point(i) - t_e1(i);
1698 add_vert(conn_edge[0], e, fabs(t_cut));
1699 add_vert(conn_edge[1], e, fabs(t_cut - 1));
1700 cut_this_edge(e, edge_length, t_ray, t_e1);
1701 }
1702 max_edge_length = std::max(max_edge_length, edge_length);
1703 };
1704
1705 auto trim_cross_edges = [&]() {
1706 if (std::min(t_dist0(i) * t_dist0(i), t_dist1(i) * t_dist1(i)) <
1707 edge_length2) {
1708
1709 auto make_pair = [&](const double d, const double te) {
1710 return std::make_pair(d, te);
1711 };
1712
1713 auto min_pair = make_pair(-1, -1);
1714
1715 for (auto f : surface_skin) {
1716
1717 auto conn_front = get_conn_edge(f);
1718 auto coords_front = get_coords_edge(conn_front);
1719 auto t_f0 = get_ftensor_coords(&coords_front[0]);
1720 auto t_f1 = get_ftensor_coords(&coords_front[3]);
1721
1722 double te;
1723 double tf;
1724
1726 &t_e0(0), &t_e1(0), &t_f0(0), &t_f1(0), &te, &tf);
1727
1728 if (res != Tools::NO_SOLUTION) {
1729
1730 auto check = [](auto v) {
1731 return (v > -std::numeric_limits<double>::epsilon() &&
1732 (v - 1) < std::numeric_limits<double>::epsilon());
1733 };
1734
1735 if (check(te) && check(tf)) {
1736
1737 FTensor::Tensor1<double, 3> t_delta, t_cross;
1738 t_delta(i) = (t_e0(i) + te * t_edge_delta(i)) -
1739 (t_f0(0) + tf * (t_f1(i) - t_f0(i)));
1740
1741 t_delta(i) /= edge_length;
1742
1743 const double dot = t_delta(i) * t_delta(i);
1744
1745 if (min_pair.first < 0 || min_pair.first > dot) {
1746
1747 if (dot < tol_trim_close * tol_trim_close)
1748 min_pair = make_pair(dot, te);
1749 }
1750 }
1751 }
1752 }
1753
1754 if (min_pair.first > -std::numeric_limits<double>::epsilon()) {
1755 add_edge(min_pair.second);
1756 return true;
1757 }
1758 }
1759
1760 return false;
1761 };
1762
1763 if (!trim_cross_edges()) {
1764
1765 const double dot = t_dist0(i) * t_dist1(i);
1766 const double dot_direction0 = t_dist0(i) * t_edge_delta(i);
1767 const double dot_direction1 = t_dist1(i) * t_edge_delta(i);
1768
1769 if (dot < std::numeric_limits<double>::epsilon() &&
1770 dot_direction0 > -std::numeric_limits<double>::epsilon() &&
1771 dot_direction1 < std::numeric_limits<double>::epsilon()) {
1772
1773 const double s0 = t_dist0(i) * t_edge_delta(i) / edge_length;
1774 const double s1 = t_dist1(i) * t_edge_delta(i) / edge_length;
1775 const double t_cut = s0 / (s0 - s1);
1776
1777 add_edge(t_cut);
1778 }
1779 }
1780 }
1781
1782 if (debug)
1783 CHKERR SaveData(moab)("edges_potentially_to_trim.vtk", cut_surface_edges);
1784
1785 if (debug)
1786 CHKERR SaveData(moab)("edges_selected_to_trim.vtk", trimEdges);
1787
1788 auto get_quality_change = [&](const Range &adj_tets, const EntityHandle &v,
1789 const VectorDouble3 &new_pos) {
1790 double q0 = 1;
1791 CHKERR m_field.getInterface<Tools>()->minTetsQuality(adj_tets, q0, th);
1792
1793 double min_q = 1;
1794 for (auto t : adj_tets) {
1795 int num_nodes;
1796 const EntityHandle *conn;
1797 CHKERR m_field.get_moab().get_connectivity(t, conn, num_nodes, true);
1798 VectorDouble12 coords(12);
1799 if (th)
1800 CHKERR moab.tag_get_data(th, conn, num_nodes, &coords[0]);
1801 else
1802 CHKERR moab.get_coords(conn, num_nodes, &coords[0]);
1803
1804 for (int n = 0; n != 4; ++n) {
1805 auto n_coords = getVectorAdaptor(&coords[3 * n], 3);
1806 if (conn[n] == v) {
1807 noalias(n_coords) = new_pos;
1808 } else {
1809 auto m = verticesOnTrimEdges.find(conn[n]);
1810 if (m != verticesOnTrimEdges.end())
1811 noalias(n_coords) =
1812 m->second.rayPoint + m->second.dIst * m->second.unitRayDir;
1813 }
1814 }
1815
1816 const double q = Tools::volumeLengthQuality(&coords[0]);
1817 if (!std::isnormal(q)) {
1818 min_q = -2;
1819 break;
1820 }
1821 min_q = std::min(min_q, q);
1822 }
1823 return min_q / q0;
1824 };
1825
1826 Range checked_verts;
1827 auto add_trim_vert = [&](const EntityHandle v, const EntityHandle e) {
1828 if (!(trimNewVertices.find(v) != trimNewVertices.end())) {
1829 auto &r = edgesToTrim.at(e);
1830 VectorDouble3 ray_point = get_point_coords(v);
1831 VectorDouble3 new_pos = r.rayPoint + r.dIst * r.unitRayDir;
1832 VectorDouble3 unit_ray_dir = ray_point - new_pos;
1833 Range adj_tets;
1834 CHKERR moab.get_adjacencies(&v, 1, 3, false, adj_tets);
1835 adj_tets = intersect(adj_tets, cutNewVolumes);
1836 double q = get_quality_change(adj_tets, v, new_pos);
1838 VectorDouble3 ray_dir = new_pos - ray_point;
1839 double dist = norm_2(ray_dir);
1840 verticesOnTrimEdges[v].dIst = 1;
1841 verticesOnTrimEdges[v].lEngth = dist;
1842 verticesOnTrimEdges[v].unitRayDir = ray_dir;
1843 verticesOnTrimEdges[v].rayPoint = ray_point;
1844 trimNewVertices.insert(v);
1845 }
1846 checked_verts.insert(v);
1847 }
1848 };
1849
1850 auto add_no_move_trim = [&](const EntityHandle v, const EntityHandle e) {
1851 if (!(trimNewVertices.find(v) != trimNewVertices.end())) {
1852 auto &m = edgesToTrim.at(e);
1854 verticesOnTrimEdges[v].rayPoint = get_point_coords(v);
1855 verticesOnTrimEdges[v].lEngth = 0;
1856 verticesOnTrimEdges[v].dIst = 0;
1857 trimNewVertices.insert(v);
1858 checked_verts.insert(v);
1859 }
1860 };
1861
1862 // Iterate over all vertives close to surface front and check if those can
1863 // be moved
1864
1865 // filer nodes which distance is in given tolerance
1866 auto hi = verts_map.get<0>().upper_bound(tol_trim_close);
1867 verts_map.get<0>().erase(hi, verts_map.end());
1868
1869 auto remove_verts = [&](Range nodes) {
1870 for (auto n : nodes) {
1871 auto r = verts_map.get<1>().equal_range(n);
1872 verts_map.get<1>().erase(r.first, r.second);
1873 }
1874 };
1875
1876 auto trim_verts = [&](const Range verts, const bool move) {
1877 VertMapMultiIndex verts_map_tmp;
1878 for (auto p = corners.pair_begin(); p != corners.pair_end(); ++p) {
1879 auto lo = verts_map.get<1>().lower_bound(p->first);
1880 auto hi = verts_map.get<1>().upper_bound(p->second);
1881 verts_map_tmp.insert(lo, hi);
1882 }
1883 if (move) {
1884 for (auto &m : verts_map_tmp.get<0>())
1885 add_trim_vert(m.v, m.e);
1886 } else {
1887 for (auto &m : verts_map_tmp.get<0>())
1888 add_no_move_trim(m.v, m.e);
1889 }
1890 };
1891
1892 auto trim_edges = [&](const Range ents, const bool move) {
1893 VertMapMultiIndex verts_map_tmp;
1894 for (auto p = ents.pair_begin(); p != ents.pair_end(); ++p) {
1895 auto lo = verts_map.get<2>().lower_bound(p->first);
1896 auto hi = verts_map.get<2>().upper_bound(p->second);
1897 verts_map_tmp.insert(lo, hi);
1898 verts_map.get<2>().erase(lo, hi);
1899 }
1900 if (move) {
1901 for (auto &m : verts_map_tmp.get<0>())
1902 add_trim_vert(m.v, m.e);
1903 } else {
1904 for (auto &m : verts_map_tmp.get<0>())
1905 add_no_move_trim(m.v, m.e);
1906 }
1907 };
1908
1909 auto intersect_self = [&](Range &a, const Range b) { a = intersect(a, b); };
1910
1911 map<std::string, Range> range_maps;
1912
1913 // Edges on surface skin
1914 CHKERR skin.find_skin(0, cutNewSurfaces, false, range_maps["surface_skin"]);
1915 // Edges as trimmed
1916 intersect_self(range_maps["surface_skin"], trimEdges);
1917 // Edges are on fixed edges
1918 range_maps["fixed_edges_on_surface_skin"] =
1919 intersect(range_maps["surface_skin"], fix_edges);
1920
1921 // Edges adjacent to fixed edges
1922 CHKERR moab.get_adjacencies(fixed_edges_verts, 1, false,
1923 range_maps["fixed_edges_verts_edges"],
1924 moab::Interface::UNION);
1925 intersect_self(range_maps["fixed_edges_verts_edges"], trimEdges);
1926
1927 CHKERR moab.get_connectivity(
1928 range_maps["fixed_edges_verts_edges"],
1929 range_maps["fixed_edges_verts_edges_verts_on_the_skin"], false);
1930 intersect_self(range_maps["fixed_edges_verts_edges_verts_on_the_skin"],
1931 tets_skin_verts);
1932
1933 // Do not move nodes at the corners
1934 trim_verts(corners, false);
1935 remove_verts(corners);
1936
1937 // Trim edges on the body skin move
1938 trim_edges(range_maps["fixed_edges_on_surface_skin"], true);
1939 remove_verts(range_maps["fixed_edges_on_surface_skin_verts"]);
1940
1941 // Trim edges nodes on the skin but edge not on the skin do not
1942 trim_verts(range_maps["fixed_edges_verts_edges_verts_on_the_skin"], false);
1943 remove_verts(range_maps["fixed_edges_verts_edges_verts_on_the_skin"]);
1944
1945 // Trim edges on the cut skin move
1946 trim_edges(range_maps["surface_skin"], true);
1947 trim_verts(tets_skin_verts, false);
1948 remove_verts(tets_skin_verts);
1949
1950 for (auto &m : verts_map.get<0>())
1951 add_trim_vert(m.v, m.e);
1952
1953 for (auto m : verticesOnTrimEdges) {
1954 EntityHandle v = m.first;
1955 Range adj;
1956 CHKERR moab.get_adjacencies(&v, 1, 1, false, adj);
1957 for (auto e : adj) {
1958 edgesToTrim.erase(e);
1959 trimEdges.erase(e);
1960 }
1961 }
1962
1963 if (debug)
1964 if (!trimNewVertices.empty())
1965 CHKERR SaveData(moab)("trim_close_vertices.vtk", trimNewVertices);
1966
1967 if (debug)
1968 if (!trimEdges.empty())
1969 CHKERR SaveData(moab)("trim_edges.vtk", trimEdges);
1970
1972}
1973
1975 const bool debug) {
1976 CoreInterface &m_field = cOre;
1977 moab::Interface &moab = m_field.get_moab();
1978 MeshRefinement *refiner;
1979 auto ref_ents_ptr = m_field.get_ref_ents();
1981
1982 CHKERR m_field.getInterface(refiner);
1985
1986 trimNewVolumes.clear();
1987 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
1988 bit, BitRefLevel().set(), MBTET, trimNewVolumes);
1989
1990 for (map<EntityHandle, TreeData>::iterator mit = edgesToTrim.begin();
1991 mit != edgesToTrim.end(); mit++) {
1992 auto vit = ref_ents_ptr->get<Composite_ParentEnt_And_EntType_mi_tag>().find(
1993 boost::make_tuple(MBVERTEX, mit->first));
1994 if (vit ==
1995 ref_ents_ptr->get<Composite_ParentEnt_And_EntType_mi_tag>().end())
1996 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1997 "No vertex on trim edges, that make no sense");
1998
1999 const boost::shared_ptr<RefEntity> &ref_ent = *vit;
2000 if ((ref_ent->getBitRefLevel() & bit).any()) {
2001 EntityHandle vert = ref_ent->getEnt();
2002 trimNewVertices.insert(vert);
2003 verticesOnTrimEdges[vert] = mit->second;
2004 }
2005 }
2006
2007 // Get faces which are trimmed
2008 trimNewSurfaces.clear();
2009 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
2010 bit, bit, MBTRI, trimNewSurfaces);
2011
2012 Range trim_new_surfaces_nodes;
2013 CHKERR moab.get_connectivity(trimNewSurfaces, trim_new_surfaces_nodes, true);
2014 trim_new_surfaces_nodes = subtract(trim_new_surfaces_nodes, trimNewVertices);
2015 trim_new_surfaces_nodes = subtract(trim_new_surfaces_nodes, cutNewVertices);
2016 Range faces_not_on_surface;
2017 CHKERR moab.get_adjacencies(trim_new_surfaces_nodes, 2, false,
2018 faces_not_on_surface, moab::Interface::UNION);
2019 trimNewSurfaces = subtract(trimNewSurfaces, faces_not_on_surface);
2020
2021 // Get surfaces which are not trimmed and add them to surface
2022 Range all_surfaces_on_bit_level;
2023 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
2024 bit, BitRefLevel().set(), MBTRI, all_surfaces_on_bit_level);
2025 all_surfaces_on_bit_level =
2026 intersect(all_surfaces_on_bit_level, cutNewSurfaces);
2027
2028 trimNewSurfaces = unite(trimNewSurfaces, all_surfaces_on_bit_level);
2029
2031}
2032
2034 Range *corner_nodes,
2035 const double tol, Tag th,
2036 const bool debug) {
2037 CoreInterface &m_field = cOre;
2038 moab::Interface &moab = m_field.get_moab();
2039 Skinner skin(&moab);
2040
2041 auto get_adj = [&moab](const Range r, const int dim) {
2042 Range a;
2043 if (dim)
2044 CHKERR moab.get_adjacencies(r, dim, false, a, moab::Interface::UNION);
2045 else
2046 CHKERR moab.get_connectivity(r, a, true);
2047 return a;
2048 };
2049
2050 auto get_skin = [&skin](const auto r) {
2051 Range s;
2052 CHKERR skin.find_skin(0, r, false, s);
2053 return s;
2054 };
2055
2057
2058 auto trim_tets_skin = get_skin(trimNewVolumes);
2059 auto trim_tets_skin_edges = get_adj(trim_tets_skin, 1);
2060 auto trim_surface_edges = get_adj(trimNewSurfaces, 1);
2061
2062 auto contarain_edges =
2063 intersect(trim_surface_edges, get_adj(constrainSurface, 1));
2064 if (fixed_edges)
2065 contarain_edges.merge(
2066 intersect(fixed_edges->subset_by_type(MBEDGE), trim_surface_edges));
2067
2068 Range barrier_vertices(trimNewVertices);
2069 if (corner_nodes) {
2070 // Add nodes which re adjacent to corner nodes
2071 barrier_vertices.merge(get_adj(
2072 intersect(get_adj(intersect(*corner_nodes, barrier_vertices), 2),
2074 0));
2075 }
2076
2077 auto get_nodes_with_one_node_on_fixed_edge_other_not = [&]() {
2078 const auto fixed_edges_on_surface =
2079 intersect(*fixed_edges, trim_surface_edges);
2080 const auto skin_fixed_edges_on_surface = get_skin(fixed_edges_on_surface);
2081 const auto barrier_nodes = subtract(skin_fixed_edges_on_surface,
2082 get_adj(get_skin(trimNewSurfaces), 0));
2083 return barrier_nodes;
2084 };
2085
2086 if (fixed_edges)
2087 barrier_vertices.merge(get_nodes_with_one_node_on_fixed_edge_other_not());
2088
2089 auto add_close_surface_barrier = [&]() {
2091
2092 // Ten algorytm dedykuje destylarni Bowmore. Jeżeli coś pójdzie nie tak, coś
2093 // pęknie inaczej niż trzeba, to pełną odpowiedzialność ponosi Beam Suntory
2094 // UK Ltd.
2095
2096 CHKERR buildTree();
2097
2098 auto test_edges =
2099 subtract(trim_surface_edges, get_adj(constrainSurface, 1));
2100 if (fixed_edges)
2101 test_edges = subtract(test_edges, *fixed_edges);
2102 auto trim_surface_nodes = get_adj(test_edges, 0);
2103 trim_surface_nodes = subtract(trim_surface_nodes, barrier_vertices);
2104 if (fixed_edges)
2105 trim_surface_nodes =
2106 subtract(trim_surface_nodes, get_adj(*fixed_edges, 0));
2107
2108 Range trim_skin;
2109 trim_skin.merge(contarain_edges);
2110 if (fRont.empty())
2111 trim_skin.merge(get_skin(sUrface));
2112 else
2113 trim_skin.merge(fRont);
2114 if (fixed_edges)
2115 trim_skin.merge(intersect(*fixed_edges, trim_surface_edges));
2116
2117 if (debug && !trim_skin.empty())
2118 CHKERR SaveData(m_field.get_moab())("trim_skin.vtk", trim_skin);
2119
2120 VectorDouble coords(3 * trim_surface_nodes.size());
2121 if (th)
2122 CHKERR moab.tag_get_data(th, trim_surface_nodes, &*coords.begin());
2123 else
2124 CHKERR moab.get_coords(trim_surface_nodes, &*coords.begin());
2125
2126 if (!trim_skin.empty()) {
2127
2128 std::vector<double> min_distances(trim_surface_nodes.size(), -1);
2129 CHKERR cOre.getInterface<Tools>()->findMinDistanceFromTheEdges(
2130 &*coords.begin(), trim_surface_nodes.size(), trim_skin,
2131 &*min_distances.begin());
2132
2133 auto coords_ptr = &*coords.begin();
2134 auto min_dist = &*min_distances.begin();
2135
2136 std::vector<EntityHandle> add_nodes;
2137 add_nodes.reserve(trim_surface_nodes.size());
2138
2139 for (auto n : trim_surface_nodes) {
2140
2141 if ((*min_dist) > std::numeric_limits<double>::epsilon()) {
2142 VectorDouble3 point_out(3);
2143
2144 EntityHandle facets_out;
2145 CHKERR treeSurfPtr->closest_to_location(coords_ptr, rootSetSurf,
2146 &point_out[0], facets_out);
2147
2148 VectorDouble3 delta = point_out - getVectorAdaptor(coords_ptr, 3);
2149 VectorDouble3 normal(3);
2150 CHKERR m_field.getInterface<Tools>()->getTriNormal(facets_out,
2151 &normal[0]);
2152
2153 normal /= norm_2(normal);
2154 delta -= inner_prod(normal, delta) * normal;
2155
2156 double dist = norm_2(delta);
2157 if (dist < tol * (*min_dist))
2158 add_nodes.emplace_back(n);
2159 }
2160
2161 coords_ptr += 3;
2162 min_dist += 1;
2163 }
2164
2165 barrier_vertices.insert_list(add_nodes.begin(), add_nodes.end());
2166 }
2167
2169 };
2170
2171 auto remove_faces_on_skin = [&]() {
2173 Range skin_faces = intersect(trimNewSurfaces, trim_tets_skin);
2174 if (!skin_faces.empty()) {
2175 barrier_vertices.merge(get_adj(skin_faces, 0));
2176 for (auto f : skin_faces)
2177 trimNewSurfaces.erase(f);
2178 }
2180 };
2181
2182 auto get_trim_free_edges = [&]() {
2183 // get current surface skin
2184 Range trim_surf_skin;
2185 CHKERR skin.find_skin(0, trimNewSurfaces, false, trim_surf_skin);
2186 trim_surf_skin = subtract(trim_surf_skin, trim_tets_skin_edges);
2187
2188 Range trim_surf_skin_verts;
2189 CHKERR moab.get_connectivity(trim_surf_skin, trim_surf_skin_verts, true);
2190 for (auto e : barrier_vertices)
2191 trim_surf_skin_verts.erase(e);
2192
2193 Range free_edges;
2194 CHKERR moab.get_adjacencies(trim_surf_skin_verts, 1, false, free_edges,
2195 moab::Interface::UNION);
2196 free_edges =
2197 subtract(intersect(free_edges, trim_surf_skin), contarain_edges);
2198
2199 return free_edges;
2200 };
2201
2202 CHKERR remove_faces_on_skin();
2203 CHKERR add_close_surface_barrier();
2204
2205 if (debug && !barrier_vertices.empty())
2206 CHKERR SaveData(m_field.get_moab())("barrier_vertices.vtk",
2207 barrier_vertices);
2208
2209 int nn = 0;
2210 Range out_edges;
2211 while (!(out_edges = get_trim_free_edges()).empty()) {
2212
2213 if (debug && !trimNewSurfaces.empty())
2214 CHKERR SaveData(m_field.get_moab())(
2215 "trimNewSurfaces_" + boost::lexical_cast<std::string>(nn) + ".vtk",
2217
2218 if (debug && !out_edges.empty())
2219 CHKERR SaveData(m_field.get_moab())(
2220 "trimNewSurfacesOutsideVerts_" +
2221 boost::lexical_cast<std::string>(nn) + ".vtk",
2222 out_edges);
2223
2224 Range outside_faces;
2225 CHKERR moab.get_adjacencies(out_edges, 2, false, outside_faces,
2226 moab::Interface::UNION);
2227 trimNewSurfaces = subtract(trimNewSurfaces, outside_faces);
2228 ++nn;
2229 }
2230
2231 if (debug && !trimNewSurfaces.empty())
2232 CHKERR SaveData(m_field.get_moab())(
2233 "trimNewSurfaces_" + boost::lexical_cast<std::string>(nn) + ".vtk",
2235
2237}
2238
2241 Range &ents) {
2242 CoreInterface &m_field = cOre;
2243 moab::Interface &moab = m_field.get_moab();
2244 PrismInterface *interface;
2246 CHKERR m_field.getInterface(interface);
2247 // Remove tris on surface front
2248 {
2249 Range front_tris;
2250 EntityHandle meshset;
2251 CHKERR moab.create_meshset(MESHSET_SET, meshset);
2252 CHKERR moab.add_entities(meshset, ents);
2254 meshset, split_bit, true, front_tris);
2255 CHKERR moab.delete_entities(&meshset, 1);
2256 ents = subtract(ents, front_tris);
2257 }
2258 // Remove entities on skin
2259 Range tets;
2260 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
2261 split_bit, BitRefLevel().set(), MBTET, tets);
2262 // Remove entities on skin
2263 Skinner skin(&moab);
2264 Range tets_skin;
2265 CHKERR skin.find_skin(0, tets, false, tets_skin);
2266 ents = subtract(ents, tets_skin);
2267
2269}
2270
2272 const BitRefLevel bit,
2273 const Range &ents, Tag th) {
2274 CoreInterface &m_field = cOre;
2275 moab::Interface &moab = m_field.get_moab();
2276 PrismInterface *interface;
2278 CHKERR m_field.getInterface(interface);
2279 EntityHandle meshset_volume;
2280 CHKERR moab.create_meshset(MESHSET_SET, meshset_volume);
2281 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
2282 split_bit, BitRefLevel().set(), MBTET, meshset_volume);
2283 EntityHandle meshset_surface;
2284 CHKERR moab.create_meshset(MESHSET_SET, meshset_surface);
2285 CHKERR moab.add_entities(meshset_surface, ents);
2286 CHKERR interface->getSides(meshset_surface, split_bit, true);
2287 CHKERR interface->splitSides(meshset_volume, bit, meshset_surface, true,
2288 true);
2289 CHKERR moab.delete_entities(&meshset_volume, 1);
2290 CHKERR moab.delete_entities(&meshset_surface, 1);
2291 if (th) {
2292 Range prisms;
2293 ierr = m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
2294 bit, BitRefLevel().set(), MBPRISM, prisms);
2295 for (Range::iterator pit = prisms.begin(); pit != prisms.end(); pit++) {
2296 int num_nodes;
2297 const EntityHandle *conn;
2298 CHKERR moab.get_connectivity(*pit, conn, num_nodes, true);
2299 MatrixDouble data(3, 3);
2300 CHKERR moab.tag_get_data(th, conn, 3, &data(0, 0));
2301 CHKERR moab.tag_set_data(th, &conn[3], 3, &data(0, 0));
2302 }
2303 }
2305}
2306
2308 const int fraction_level, const Range &tets, const Range &surface,
2309 const Range &fixed_edges, const Range &corner_nodes, Range &edges_to_merge,
2310 Range &out_tets, Range &new_surf, Tag th, const bool update_meshsets,
2311 const BitRefLevel *bit_ptr, const bool debug) {
2312 CoreInterface &m_field = cOre;
2313 moab::Interface &moab = m_field.get_moab();
2315
2316 /**
2317 * \brief Merge nodes
2318 */
2319 struct MergeNodes {
2320 CoreInterface &mField;
2321 const bool onlyIfImproveQuality;
2322 Tag tH;
2323 bool updateMehsets;
2324
2325 MergeNodes(CoreInterface &m_field, const bool only_if_improve_quality,
2326 Tag th, bool update_mehsets)
2327 : mField(m_field), onlyIfImproveQuality(only_if_improve_quality),
2328 tH(th), updateMehsets(update_mehsets) {
2329 mField.getInterface(nodeMergerPtr);
2330 }
2331 NodeMergerInterface *nodeMergerPtr;
2332 MoFEMErrorCode mergeNodes(int line_search, EntityHandle father,
2333 EntityHandle mother, Range &proc_tets,
2334 bool add_child = true) {
2335 moab::Interface &moab = mField.get_moab();
2337 const EntityHandle conn[] = {father, mother};
2338 Range vert_tets;
2339 CHKERR moab.get_adjacencies(conn, 2, 3, false, vert_tets,
2340 moab::Interface::UNION);
2341 vert_tets = intersect(vert_tets, proc_tets);
2342 Range out_tets;
2343 CHKERR nodeMergerPtr->mergeNodes(father, mother, out_tets, &vert_tets,
2344 onlyIfImproveQuality, 0, line_search,
2345 tH);
2346
2347 if (add_child && nodeMergerPtr->getSuccessMerge()) {
2348
2349 Range::iterator lo, hi = proc_tets.begin();
2350 for (auto pt = vert_tets.pair_begin(); pt != vert_tets.pair_end();
2351 ++pt) {
2352 lo = proc_tets.lower_bound(hi, proc_tets.end(), pt->first);
2353 if (lo != proc_tets.end()) {
2354 hi = proc_tets.upper_bound(lo, proc_tets.end(), pt->second);
2355 proc_tets.erase(lo, hi);
2356 } else
2357 break;
2358 }
2359 proc_tets.merge(out_tets);
2360
2361 auto &parent_child_map = nodeMergerPtr->getParentChildMap();
2362
2363 struct ChangeChild {
2364 EntityHandle child;
2365 ChangeChild(const EntityHandle child) : child(child) {}
2366 void operator()(NodeMergerInterface::ParentChild &p) {
2367 p.cHild = child;
2368 }
2369 };
2370
2371 std::vector<decltype(parentsChildMap.get<0>().begin())> it_vec;
2372 it_vec.reserve(parentsChildMap.size());
2373
2374 for (auto &p : parent_child_map) {
2375
2376 it_vec.clear();
2377 for (auto it = parentsChildMap.get<0>().equal_range(p.pArent);
2378 it.first != it.second; ++it.first)
2379 it_vec.emplace_back(it.first);
2380
2381 for (auto it = parentsChildMap.get<1>().equal_range(p.pArent);
2382 it.first != it.second; ++it.first)
2383 it_vec.emplace_back(parentsChildMap.project<0>(it.first));
2384
2385 for (auto &it : it_vec)
2386 parentsChildMap.modify(it, ChangeChild(p.cHild));
2387 }
2388
2389 parentsChildMap.insert(parent_child_map.begin(),
2390 parent_child_map.end());
2391 }
2393 }
2394
2395 MoFEMErrorCode updateRangeByChilds(Range &new_surf, Range &edges_to_merge,
2396 Range &not_merged_edges,
2397 bool add_child) {
2398 moab::Interface &moab = mField.get_moab();
2400 if (add_child) {
2401
2402 std::vector<EntityHandle> parents_ents_vec(parentsChildMap.size());
2403 for (auto &it : parentsChildMap)
2404 parents_ents_vec.emplace_back(it.pArent);
2405 Range parent_ents;
2406 parent_ents.insert_list(parents_ents_vec.begin(),
2407 parents_ents_vec.end());
2408
2409 Range surf_parent_ents = intersect(new_surf, parent_ents);
2410 new_surf = subtract(new_surf, surf_parent_ents);
2411 Range child_surf_ents;
2412 CHKERR updateRangeByChilds(parentsChildMap, surf_parent_ents,
2413 child_surf_ents);
2414 new_surf.merge(child_surf_ents);
2415
2416 Range edges_to_merge_parent_ents =
2417 intersect(edges_to_merge, parent_ents);
2418 edges_to_merge = subtract(edges_to_merge, edges_to_merge_parent_ents);
2419 Range merged_child_edge_ents;
2420 CHKERR updateRangeByChilds(parentsChildMap, edges_to_merge_parent_ents,
2421 merged_child_edge_ents);
2422
2423 Range not_merged_edges_child_ents =
2424 intersect(not_merged_edges, parent_ents);
2425 not_merged_edges =
2426 subtract(not_merged_edges, not_merged_edges_child_ents);
2427 Range not_merged_child_edge_ents;
2428 CHKERR updateRangeByChilds(parentsChildMap, not_merged_edges_child_ents,
2429 not_merged_child_edge_ents);
2430
2431 merged_child_edge_ents =
2432 subtract(merged_child_edge_ents, not_merged_child_edge_ents);
2433 edges_to_merge.merge(merged_child_edge_ents);
2434 not_merged_edges.merge(not_merged_child_edge_ents);
2435
2436 if (updateMehsets) {
2438 (*mField.getInterface<MeshsetsManager>()), cubit_it)) {
2439 EntityHandle cubit_meshset = cubit_it->meshset;
2440 Range meshset_ents;
2441 CHKERR moab.get_entities_by_handle(cubit_meshset, meshset_ents,
2442 true);
2443 Range child_ents;
2444 CHKERR updateRangeByChilds(parentsChildMap, meshset_ents,
2445 child_ents);
2446 CHKERR moab.add_entities(cubit_meshset, child_ents);
2447 }
2448 }
2449 }
2450
2452 };
2453
2454 private:
2456 std::vector<EntityHandle> childsVec;
2457
2458 inline MoFEMErrorCode updateRangeByChilds(
2459 const NodeMergerInterface::ParentChildMap &parent_child_map,
2460 const Range &parents, Range &childs) {
2462 childsVec.clear();
2463 childsVec.reserve(parents.size());
2464 for (auto pit = parents.pair_begin(); pit != parents.pair_end(); pit++) {
2465 auto it = parent_child_map.lower_bound(pit->first);
2466 if (it != parent_child_map.end()) {
2467 for (auto hi_it = parent_child_map.upper_bound(pit->second);
2468 it != hi_it; ++it)
2469 childsVec.emplace_back(it->cHild);
2470 }
2471 }
2472 childs.insert_list(childsVec.begin(), childsVec.end());
2474 }
2475 };
2476
2477 /**
2478 * \brief Calculate edge length
2479 */
2480 struct LengthMap {
2481 Tag tH;
2482 CoreInterface &mField;
2483 moab::Interface &moab;
2484 const double maxLength;
2485 LengthMap(CoreInterface &m_field, Tag th, double max_length)
2486 : tH(th), mField(m_field), moab(m_field.get_moab()),
2487 maxLength(max_length) {
2488 gammaL = 1.;
2489 gammaQ = 3.;
2490 acceptedThrasholdMergeQuality = 0.5;
2491 }
2492
2493 double
2494 gammaL; ///< Controls importance of length when ranking edges for merge
2495 double
2496 gammaQ; ///< Controls importance of quality when ranking edges for merge
2497 double acceptedThrasholdMergeQuality; ///< Do not merge quality if quality
2498 ///< above accepted thrashold
2499
2500 MoFEMErrorCode operator()(const Range &tets, const Range &edges,
2501 LengthMapData_multi_index &length_map,
2502 double &ave) const {
2503 int num_nodes;
2504 const EntityHandle *conn;
2505 std::array<double, 12> coords;
2507 VectorAdaptor s0(3, ublas::shallow_array_adaptor<double>(3, &coords[0]));
2508 VectorAdaptor s1(3, ublas::shallow_array_adaptor<double>(3, &coords[3]));
2510
2511 struct NodeQuality {
2512 EntityHandle node;
2513 double quality;
2514 NodeQuality(const EntityHandle node) : node(node), quality(1) {}
2515 };
2516
2517 typedef multi_index_container<
2518 NodeQuality, indexed_by<
2519
2520 sequenced<>,
2521
2522 hashed_non_unique<tag<Ent_mi_tag>,
2523 member<NodeQuality, EntityHandle,
2524 &NodeQuality::node>>
2525
2526 >>
2527 NodeQuality_sequence;
2528
2529 NodeQuality_sequence node_quality_sequence;
2530
2531 Range edges_nodes;
2532 CHKERR moab.get_connectivity(tets, edges_nodes, false);
2533 Range edges_tets;
2534 CHKERR moab.get_adjacencies(edges, 3, false, edges_tets,
2535 moab::Interface::UNION);
2536 edges_tets = intersect(edges_tets, tets);
2537
2538 for (auto node : edges_nodes)
2539 node_quality_sequence.emplace_back(node);
2540
2541 for (auto tet : edges_tets) {
2542
2543 CHKERR moab.get_connectivity(tet, conn, num_nodes, true);
2544 if (tH)
2545 CHKERR moab.tag_get_data(tH, conn, num_nodes, coords.data());
2546 else
2547 CHKERR moab.get_coords(conn, num_nodes, coords.data());
2548
2549 const double q = Tools::volumeLengthQuality(coords.data());
2550
2551 for (auto n : {0, 1, 2, 3}) {
2552 auto it = node_quality_sequence.get<1>().find(conn[n]);
2553 if (it != node_quality_sequence.get<1>().end())
2554 const_cast<double &>(it->quality) = std::min(q, it->quality);
2555 }
2556 }
2557
2558 for (auto edge : edges) {
2559
2560 CHKERR moab.get_connectivity(edge, conn, num_nodes, true);
2561
2562 if (tH)
2563 CHKERR moab.tag_get_data(tH, conn, num_nodes, coords.data());
2564 else
2565 CHKERR moab.get_coords(conn, num_nodes, coords.data());
2566
2567 double q = 1;
2568 for (auto n : {0, 1}) {
2569 auto it = node_quality_sequence.get<1>().find(conn[n]);
2570 if (it != node_quality_sequence.get<1>().end())
2571 q = std::min(q, it->quality);
2572 }
2573
2574 if (q < acceptedThrasholdMergeQuality) {
2575 noalias(delta) = (s0 - s1) / maxLength;
2576 double dot = inner_prod(delta, delta);
2577 double val = pow(q, gammaQ) * pow(dot, 0.5 * gammaL);
2578 length_map.insert(LengthMapData(val, q, edge));
2579 }
2580 }
2581
2582 ave = 0;
2583 for (LengthMapData_multi_index::nth_index<0>::type::iterator mit =
2584 length_map.get<0>().begin();
2585 mit != length_map.get<0>().end(); mit++) {
2586 ave += mit->qUality;
2587 }
2588 ave /= length_map.size();
2590 }
2591 };
2592
2593 /**
2594 * \brief Topological relations
2595 */
2596 struct Toplogy {
2597
2598 CoreInterface &mField;
2599 Tag tH;
2600 Toplogy(CoreInterface &m_field, Tag th) : mField(m_field), tH(th) {}
2601
2602 enum TYPE {
2603 FREE = 0,
2604 SKIN = 1 << 0,
2605 SURFACE = 1 << 1,
2606 SURFACE_SKIN = 1 << 2,
2607 FRONT_ENDS = 1 << 3,
2608 FIX_EDGES = 1 << 4,
2609 FIX_CORNERS = 1 << 5
2610 };
2611
2612 typedef map<int, Range> SetsMap;
2613
2614 MoFEMErrorCode classifyVerts(const Range &surface, const Range &tets,
2615 const Range &fixed_edges,
2616 const Range &corner_nodes,
2617 const Range &constrain_surface,
2618 SetsMap &sets_map) const {
2619 moab::Interface &moab(mField.get_moab());
2620 Skinner skin(&moab);
2622
2623 sets_map[FIX_CORNERS].merge(corner_nodes);
2624 Range fixed_verts;
2625 CHKERR moab.get_connectivity(fixed_edges, fixed_verts, true);
2626 sets_map[FIX_EDGES].swap(fixed_verts);
2627
2628 Range tets_skin;
2629 CHKERR skin.find_skin(0, tets, false, tets_skin);
2630 Range tets_skin_edges;
2631 CHKERR moab.get_adjacencies(tets_skin, 1, false, tets_skin_edges,
2632 moab::Interface::UNION);
2633
2634 // constrain surface
2635 Range constrain_surface_verts;
2636 CHKERR moab.get_connectivity(constrain_surface, constrain_surface_verts,
2637 true);
2638 Range constrain_surface_edges;
2639 CHKERR moab.get_adjacencies(constrain_surface, 1, false,
2640 constrain_surface_edges,
2641 moab::Interface::UNION);
2642
2643 // surface skin
2644 Range surface_skin;
2645 CHKERR skin.find_skin(0, surface, false, surface_skin);
2646 Range front_in_the_body;
2647 front_in_the_body = subtract(surface_skin, tets_skin_edges);
2648 Range front_in_the_body_verts;
2649 CHKERR moab.get_connectivity(front_in_the_body, front_in_the_body_verts,
2650 true);
2651 Range front_ends;
2652 CHKERR skin.find_skin(0, front_in_the_body, false, front_ends);
2653 front_ends.merge(
2654 intersect(front_in_the_body_verts, constrain_surface_verts));
2655 sets_map[FRONT_ENDS].swap(front_ends);
2656
2657 Range surface_skin_verts;
2658 CHKERR moab.get_connectivity(surface_skin, surface_skin_verts, true);
2659 sets_map[SURFACE_SKIN].swap(surface_skin_verts);
2660
2661 // surface
2662 Range surface_verts;
2663 CHKERR moab.get_connectivity(surface, surface_verts, true);
2664 sets_map[SURFACE_SKIN].merge(
2665 intersect(constrain_surface_verts, surface_verts));
2666 sets_map[SURFACE].swap(surface_verts);
2667
2668 // skin
2669 Range tets_skin_verts;
2670 CHKERR moab.get_connectivity(tets_skin, tets_skin_verts, true);
2671 sets_map[SKIN].swap(tets_skin_verts);
2672 sets_map[SKIN].merge(constrain_surface_verts);
2673
2674 Range tets_verts;
2675 CHKERR moab.get_connectivity(tets, tets_verts, true);
2676 sets_map[FREE].swap(tets_verts);
2677
2679 }
2680
2681 MoFEMErrorCode getProcTets(const Range &tets, const Range &edges_to_merge,
2682 Range &proc_tets) const {
2683 moab::Interface &moab(mField.get_moab());
2685 Range edges_to_merge_verts;
2686 CHKERR moab.get_connectivity(edges_to_merge, edges_to_merge_verts, true);
2687 Range edges_to_merge_verts_tets;
2688 CHKERR moab.get_adjacencies(edges_to_merge_verts, 3, false,
2689 edges_to_merge_verts_tets,
2690 moab::Interface::UNION);
2691 edges_to_merge_verts_tets = intersect(edges_to_merge_verts_tets, tets);
2692 proc_tets.swap(edges_to_merge_verts_tets);
2694 }
2695
2696 MoFEMErrorCode removeBadEdges(const Range &surface, const Range &tets,
2697 const Range &fixed_edges,
2698 const Range &corner_nodes,
2699 const Range &constrain_surface,
2700 Range &edges_to_merge,
2701 Range &not_merged_edges) {
2702 moab::Interface &moab(mField.get_moab());
2704
2705 // find skin
2706 Skinner skin(&moab);
2707 Range tets_skin;
2708 CHKERR skin.find_skin(0, tets, false, tets_skin);
2709 Range surface_skin;
2710 CHKERR skin.find_skin(0, surface, false, surface_skin);
2711
2712 // constrain surface
2713 Range constrain_surface_verts;
2714 CHKERR moab.get_connectivity(constrain_surface, constrain_surface_verts,
2715 true);
2716 Range constrain_surface_edges;
2717 CHKERR moab.get_adjacencies(constrain_surface, 1, false,
2718 constrain_surface_edges,
2719 moab::Interface::UNION);
2720
2721 // end nodes
2722 Range tets_skin_edges;
2723 CHKERR moab.get_adjacencies(tets_skin, 1, false, tets_skin_edges,
2724 moab::Interface::UNION);
2725
2726 Range surface_front;
2727 surface_front = subtract(surface_skin, tets_skin_edges);
2728 Range surface_front_nodes;
2729 CHKERR moab.get_connectivity(surface_front, surface_front_nodes, true);
2730
2731 Range ends_nodes;
2732 CHKERR skin.find_skin(0, surface_front, false, ends_nodes);
2733 ends_nodes.merge(intersect(surface_front_nodes, constrain_surface_verts));
2734
2735 // remove bad merges
2736 surface_skin.merge(constrain_surface);
2737 tets_skin_edges.merge(constrain_surface_edges);
2738
2739 // get surface and body skin verts
2740 Range surface_edges;
2741 CHKERR moab.get_adjacencies(surface, 1, false, surface_edges,
2742 moab::Interface::UNION);
2743 // get nodes on the surface
2744 Range surface_edges_verts;
2745 CHKERR moab.get_connectivity(surface_edges, surface_edges_verts, true);
2746 // get vertices on the body skin
2747 Range tets_skin_edges_verts;
2748 CHKERR moab.get_connectivity(tets_skin_edges, tets_skin_edges_verts,
2749 true);
2750
2751 Range edges_to_remove;
2752
2753 // remove edges self connected to body skin
2754 {
2755 Range ents_nodes_and_edges;
2756 ents_nodes_and_edges.merge(tets_skin_edges_verts);
2757 ents_nodes_and_edges.merge(tets_skin_edges);
2758 CHKERR removeSelfConectingEdges(ents_nodes_and_edges, edges_to_remove,
2759 false);
2760 }
2761 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2762 not_merged_edges.merge(edges_to_remove);
2763
2764 // remove edges self connected to surface
2765 {
2766 Range ents_nodes_and_edges;
2767 ents_nodes_and_edges.merge(surface_edges_verts);
2768 ents_nodes_and_edges.merge(surface_edges);
2769 ents_nodes_and_edges.merge(tets_skin_edges_verts);
2770 ents_nodes_and_edges.merge(tets_skin_edges);
2771 CHKERR removeSelfConectingEdges(ents_nodes_and_edges, edges_to_remove,
2772 false);
2773 }
2774 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2775 not_merged_edges.merge(edges_to_remove);
2776
2777 // remove edges adjacent corner_nodes execpt those on fixed edges
2778 Range fixed_edges_nodes;
2779 CHKERR moab.get_connectivity(fixed_edges, fixed_edges_nodes, true);
2780 {
2781 Range ents_nodes_and_edges;
2782 ents_nodes_and_edges.merge(fixed_edges_nodes);
2783 ents_nodes_and_edges.merge(ends_nodes);
2784 ents_nodes_and_edges.merge(corner_nodes);
2785 ents_nodes_and_edges.merge(fixed_edges);
2786 CHKERR removeSelfConectingEdges(ents_nodes_and_edges, edges_to_remove,
2787 false);
2788 }
2789 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2790 not_merged_edges.merge(edges_to_remove);
2791
2792 // remove edges self connected to surface
2793 CHKERR removeSelfConectingEdges(surface_edges, edges_to_remove, false);
2794 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2795 not_merged_edges.merge(edges_to_remove);
2796
2797 // remove edges self contented on surface skin
2798 {
2799 Range ents_nodes_and_edges;
2800 ents_nodes_and_edges.merge(surface_skin);
2801 ents_nodes_and_edges.merge(fixed_edges_nodes);
2802 CHKERR removeSelfConectingEdges(ents_nodes_and_edges, edges_to_remove,
2803 false);
2804 }
2805 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2806 not_merged_edges.merge(edges_to_remove);
2807
2808 // remove edges connecting crack front and fixed edges, except those short
2809 {
2810 Range ents_nodes_and_edges;
2811 ents_nodes_and_edges.merge(surface_skin.subset_by_type(MBEDGE));
2812 ents_nodes_and_edges.merge(fixed_edges.subset_by_type(MBEDGE));
2813 CHKERR removeSelfConectingEdges(ents_nodes_and_edges, edges_to_remove,
2814 false);
2815 }
2816 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2817 not_merged_edges.merge(edges_to_remove);
2818
2819 // remove crack front nodes connected to the surface, except those short
2820 {
2821 Range ents_nodes_and_edges;
2822 ents_nodes_and_edges.merge(surface_front_nodes);
2823 ents_nodes_and_edges.merge(surface_front);
2824 ents_nodes_and_edges.merge(tets_skin_edges_verts);
2825 ents_nodes_and_edges.merge(tets_skin_edges);
2826 CHKERR removeSelfConectingEdges(ents_nodes_and_edges, edges_to_remove,
2827 false);
2828 }
2829 edges_to_merge = subtract(edges_to_merge, edges_to_remove);
2830 not_merged_edges.merge(edges_to_remove);
2831
2833 }
2834
2835 private:
2836 MoFEMErrorCode removeSelfConectingEdges(const Range &ents,
2837 Range &edges_to_remove,
2838 bool debug) const {
2839 moab::Interface &moab(mField.get_moab());
2841 // get nodes
2842 Range ents_nodes = ents.subset_by_type(MBVERTEX);
2843 if (ents_nodes.empty()) {
2844 CHKERR moab.get_connectivity(ents, ents_nodes, true);
2845 }
2846 // edges adj. to nodes
2847 Range ents_nodes_edges;
2848 CHKERR moab.get_adjacencies(ents_nodes, 1, false, ents_nodes_edges,
2849 moab::Interface::UNION);
2850 // nodes of adj. edges
2851 Range ents_nodes_edges_nodes;
2852 CHKERR moab.get_connectivity(ents_nodes_edges, ents_nodes_edges_nodes,
2853 true);
2854 // hanging nodes
2855 ents_nodes_edges_nodes = subtract(ents_nodes_edges_nodes, ents_nodes);
2856 Range ents_nodes_edges_nodes_edges;
2857 CHKERR moab.get_adjacencies(ents_nodes_edges_nodes, 1, false,
2858 ents_nodes_edges_nodes_edges,
2859 moab::Interface::UNION);
2860 // remove edges adj. to hanging edges
2861 ents_nodes_edges =
2862 subtract(ents_nodes_edges, ents_nodes_edges_nodes_edges);
2863 ents_nodes_edges =
2864 subtract(ents_nodes_edges, ents.subset_by_type(MBEDGE));
2865
2866 edges_to_remove.swap(ents_nodes_edges);
2867 if (debug) {
2868 CHKERR SaveData(moab)("edges_to_remove.vtk", edges_to_remove);
2869 }
2871 }
2872 };
2873
2874 Range not_merged_edges;
2875 CHKERR Toplogy(m_field, th)
2876 .removeBadEdges(surface, tets, fixed_edges, corner_nodes,
2877 constrainSurface, edges_to_merge, not_merged_edges);
2878 Toplogy::SetsMap sets_map;
2879 CHKERR Toplogy(m_field, th)
2880 .classifyVerts(surface, tets, fixed_edges, corner_nodes, constrainSurface,
2881 sets_map);
2882 if (debug) {
2883 for (Toplogy::SetsMap::reverse_iterator sit = sets_map.rbegin();
2884 sit != sets_map.rend(); sit++) {
2885 std::string name = "classification_verts_" +
2886 boost::lexical_cast<std::string>(sit->first) + ".vtk";
2887 if (!sit->second.empty())
2888 CHKERR SaveData(moab)(name, sit->second);
2889 }
2890 }
2891 Range proc_tets;
2892 CHKERR Toplogy(m_field, th).getProcTets(tets, edges_to_merge, proc_tets);
2893 out_tets = subtract(tets, proc_tets);
2894
2895 if (bit_ptr) {
2896 Range all_out_ents = out_tets;
2897 for (int dd = 2; dd >= 0; dd--) {
2898 CHKERR moab.get_adjacencies(out_tets, dd, false, all_out_ents,
2899 moab::Interface::UNION);
2900 }
2901 CHKERR m_field.getInterface<BitRefManager>()->addBitRefLevel(all_out_ents,
2902 *bit_ptr);
2903 }
2904
2905 int nb_nodes_merged = 0;
2906 LengthMapData_multi_index length_map;
2907 new_surf = surface;
2908
2909 auto save_merge_step = [&](const int pp, const Range collapsed_edges) {
2911 Range adj_faces;
2912 CHKERR moab.get_adjacencies(proc_tets, 2, false, adj_faces,
2913 moab::Interface::UNION);
2914 std::string name;
2915 name = "node_merge_step_" + boost::lexical_cast<std::string>(pp) + ".vtk";
2916 CHKERR SaveData(moab)(
2917 name, unite(intersect(new_surf, adj_faces), collapsed_edges));
2918 name =
2919 "edges_to_merge_step_" + boost::lexical_cast<std::string>(pp) + ".vtk";
2920 CHKERR SaveData(moab)(
2921 name, unite(intersect(new_surf, adj_faces), edges_to_merge));
2923 };
2924
2925 if (debug)
2926 CHKERR save_merge_step(0, Range());
2927
2928 double ave0 = 0, ave = 0, min = 0, min_p = 0, min_pp;
2929 for (int pp = 0; pp != nbMaxMergingCycles; ++pp) {
2930
2931 int nb_nodes_merged_p = nb_nodes_merged;
2932 length_map.clear();
2933 min_pp = min_p;
2934 min_p = min;
2935 CHKERR LengthMap(m_field, th, aveLength)(proc_tets, edges_to_merge,
2936 length_map, ave);
2937
2938 if(!length_map.empty())
2939 min = length_map.get<2>().begin()->qUality;
2940 if (pp == 0) {
2941 ave0 = ave;
2942 }
2943
2944 int nn = 0;
2945 Range collapsed_edges;
2946 MergeNodes merge_nodes(m_field, true, th, update_meshsets);
2947
2948 for (auto mit = length_map.get<0>().begin();
2949 mit != length_map.get<0>().end(); mit++, nn++) {
2950
2951 if (!mit->skip) {
2952
2953 auto get_father_and_mother =
2954 [&](EntityHandle &father, EntityHandle &mother, int &line_search) {
2956 int num_nodes;
2957 const EntityHandle *conn;
2958 CHKERR moab.get_connectivity(mit->eDge, conn, num_nodes, true);
2959 std::array<int, 2> conn_type = {0, 0};
2960 for (int nn = 0; nn != 2; nn++) {
2961 for (Toplogy::SetsMap::reverse_iterator sit = sets_map.rbegin();
2962 sit != sets_map.rend(); sit++) {
2963 if (sit->second.find(conn[nn]) != sit->second.end()) {
2964 conn_type[nn] |= sit->first;
2965 }
2966 }
2967 }
2968 int type_father, type_mother;
2969 if (conn_type[0] > conn_type[1]) {
2970 father = conn[0];
2971 mother = conn[1];
2972 type_father = conn_type[0];
2973 type_mother = conn_type[1];
2974 } else {
2975 father = conn[1];
2976 mother = conn[0];
2977 type_father = conn_type[1];
2978 type_mother = conn_type[0];
2979 }
2980 if (type_father == type_mother) {
2981 line_search = lineSearchSteps;
2982 }
2984 };
2985
2986 int line_search = 0;
2987 EntityHandle father, mother;
2988 CHKERR get_father_and_mother(father, mother, line_search);
2989 CHKERR merge_nodes.mergeNodes(line_search, father, mother, proc_tets);
2991 const EntityHandle father_and_mother[] = {father, mother};
2992 Range adj_tets;
2993 CHKERR moab.get_adjacencies(father_and_mother, 1, 3, false, adj_tets);
2994 Range adj_tets_nodes;
2995 CHKERR moab.get_connectivity(adj_tets, adj_tets_nodes, true);
2996 Range adj_edges;
2997 CHKERR moab.get_adjacencies(adj_tets_nodes, 1, false, adj_edges,
2998 moab::Interface::UNION);
2999 for (auto ait : adj_edges) {
3000 auto miit = length_map.get<1>().find(ait);
3001 if (miit != length_map.get<1>().end())
3002 (const_cast<LengthMapData &>(*miit)).skip = true;
3003 }
3004 nb_nodes_merged++;
3005 collapsed_edges.insert(mit->eDge);
3006 }
3007
3008 if (nn > static_cast<int>(length_map.size() / fraction_level))
3009 break;
3010 if (mit->qUality > ave)
3011 break;
3012 }
3013 }
3014
3015 CHKERR merge_nodes.updateRangeByChilds(new_surf, edges_to_merge,
3016 not_merged_edges, true);
3017
3018 MOFEM_LOG_C("WORLD", Sev::verbose,
3019 "(%d) Number of nodes merged %d ave q %3.4e min q %3.4e", pp,
3020 nb_nodes_merged, ave, min);
3021
3022 if (debug)
3023 CHKERR save_merge_step(pp + 1, collapsed_edges);
3024
3025 if (nb_nodes_merged == nb_nodes_merged_p)
3026 break;
3027 if (min > 1e-2 && min == min_pp)
3028 break;
3029 if (min > ave0)
3030 break;
3031
3032 Range adj_edges;
3033 CHKERR moab.get_adjacencies(proc_tets, 1, false, adj_edges,
3034 moab::Interface::UNION);
3035 edges_to_merge = intersect(edges_to_merge, adj_edges);
3036 CHKERR Toplogy(m_field, th)
3037 .removeBadEdges(new_surf, proc_tets, fixed_edges, corner_nodes,
3038 constrainSurface, edges_to_merge, not_merged_edges);
3039 }
3040
3041 auto reconstruct_refined_ents = [&]() {
3045 };
3046
3047 // Add function which reconstruct core multi-index. Node merging is messy
3048 // process and entity parent could be changed without notification to
3049 // multi-index. TODO Issue has to be tracked down better, however in principle
3050 // is better not to modify multi-index each time parent is changed, that makes
3051 // code slow. Is better to do it in the bulk as below.
3052 CHKERR reconstruct_refined_ents();
3053
3054 if (bit_ptr)
3055 CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevel(proc_tets,
3056 *bit_ptr);
3057
3058 out_tets.merge(proc_tets);
3059 Range adj_faces;
3060 CHKERR moab.get_adjacencies(out_tets, 2, false, adj_faces,
3061 moab::Interface::UNION);
3062 new_surf = intersect(new_surf, adj_faces);
3063
3065}
3066
3068 const int fraction_level, const BitRefLevel cut_bit,
3069 const BitRefLevel trim_bit, const BitRefLevel bit, const Range &surface,
3070 const Range &fixed_edges, const Range &corner_nodes, Tag th,
3071 const bool update_meshsets, const bool debug) {
3072 CoreInterface &m_field = cOre;
3074 Range tets_level;
3075 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
3076 trim_bit, BitRefLevel().set(), MBTET, tets_level);
3077
3078 Range edges_to_merge;
3079 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByRefLevel(
3080 trim_bit, cut_bit | trim_bit, edges_to_merge);
3081
3082 // get all entities not in database
3083 Range all_ents_not_in_database_before;
3084 CHKERR cOre.getInterface<BitRefManager>()->getAllEntitiesNotInDatabase(
3085 all_ents_not_in_database_before);
3086
3087 edges_to_merge = edges_to_merge.subset_by_type(MBEDGE);
3088 if (debug)
3089 CHKERR SaveData(m_field.get_moab())("edges_to_merge.vtk", edges_to_merge);
3090
3091 Range out_new_tets, out_new_surf;
3092 CHKERR mergeBadEdges(fraction_level, tets_level, surface, fixed_edges,
3093 corner_nodes, edges_to_merge, out_new_tets, out_new_surf,
3094 th, update_meshsets, &bit, debug);
3095
3096 // get all entities not in database after merge
3097 Range all_ents_not_in_database_after;
3098 CHKERR cOre.getInterface<BitRefManager>()->getAllEntitiesNotInDatabase(
3099 all_ents_not_in_database_after);
3100
3101 // delete hanging entities
3102 all_ents_not_in_database_after =
3103 subtract(all_ents_not_in_database_after, all_ents_not_in_database_before);
3104 Range meshsets;
3105 CHKERR m_field.get_moab().get_entities_by_type(0, MBENTITYSET, meshsets,
3106 false);
3107 for (auto m : meshsets)
3108 CHKERR m_field.get_moab().remove_entities(m,
3109 all_ents_not_in_database_after);
3110
3111 m_field.get_moab().delete_entities(all_ents_not_in_database_after);
3112
3113 mergedVolumes.swap(out_new_tets);
3114 mergedSurfaces.swap(out_new_surf);
3115
3117}
3118
3120 CoreInterface &m_field = cOre;
3121 moab::Interface &moab = m_field.get_moab();
3123 Range nodes;
3124 if (bit.none())
3125 CHKERR moab.get_entities_by_type(0, MBVERTEX, nodes);
3126 else
3127 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
3128 bit, BitRefLevel().set(), MBVERTEX, nodes);
3129 std::vector<double> coords(3 * nodes.size());
3130 CHKERR moab.get_coords(nodes, &coords[0]);
3131 CHKERR moab.tag_set_data(th, nodes, &coords[0]);
3133}
3134
3136 const BitRefLevel mask) {
3137 CoreInterface &m_field = cOre;
3138 moab::Interface &moab = m_field.get_moab();
3140 Range nodes;
3141 if (bit.none())
3142 CHKERR moab.get_entities_by_type(0, MBVERTEX, nodes);
3143 else
3144 CHKERR m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
3145 bit, mask, MBVERTEX, nodes);
3146 std::vector<double> coords(3 * nodes.size());
3147 CHKERR moab.tag_get_data(th, nodes, &coords[0]);
3148 CHKERR moab.set_coords(nodes, &coords[0]);
3150}
3151
3153 CoreInterface &m_field = cOre;
3154 moab::Interface &moab = m_field.get_moab();
3156 CHKERR SaveData(moab)(prefix + "out_vol.vtk", vOlume);
3157 CHKERR SaveData(moab)(prefix + "out_surface.vtk", sUrface);
3158 CHKERR SaveData(moab)(prefix + "out_cut_edges.vtk", cutEdges);
3159 CHKERR SaveData(moab)(prefix + "out_cut_new_volumes.vtk", cutNewVolumes);
3160 CHKERR SaveData(moab)(prefix + "out_cut_new_surfaces.vtk", cutNewSurfaces);
3161 CHKERR SaveData(moab)(prefix + "out_cut_zero_distance_ents.vtk",
3164}
3165
3167 moab::Interface &moab = cOre.getInterface<CoreInterface>()->get_moab();
3169 CHKERR SaveData(moab)("out_trim_surface.vtk", sUrface);
3170 CHKERR SaveData(moab)("out_trim_new_volumes.vtk", trimNewVolumes);
3171 CHKERR SaveData(moab)("out_trim_new_surfaces.vtk", trimNewSurfaces);
3172 CHKERR SaveData(moab)("out_trim_edges.vtk", trimEdges);
3174}
3175
3176} // namespace MoFEM
#define CutMeshFunctionBegin
static auto get_skin(MoFEM::Interface &m_field, Range body_ents)
#define MOFEM_LOG_C(channel, severity, format,...)
constexpr double a
@ QUIET
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define THROW_MESSAGE(msg)
Throw MoFEM exception.
static const bool debug
virtual const RefEntity_multiIndex * get_ref_ents() const =0
Get the ref ents object.
MoFEMErrorCode addBitRefLevel(const Range &ents, const BitRefLevel bit, int verb=QUIET) const
Add bit ref level to ref entity.
MoFEMErrorCode updateRangeByChildren(const Range &parent, Range &child, MoFEMTypes bh=MF_ZERO)
Update range by childrens.
#define MOFEM_LOG(channel, severity)
Log.
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
double tol
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorBoundedArray< double, 3 > VectorDouble3
Definition Types.hpp:92
VectorBoundedArray< double, 6 > VectorDouble6
Definition Types.hpp:95
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition Types.hpp:115
VectorBoundedArray< double, 12 > VectorDouble12
Definition Types.hpp:97
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
MatrixShallowArrayAdaptor< double > MatrixAdaptor
Matrix adaptor.
Definition Types.hpp:132
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition Templates.hpp:31
static const bool debug
MoFEMErrorCode reconstructMultiIndex(const MI &mi, MO &&mo=Modify_change_nothing())
Template used to reconstruct multi-index.
constexpr double t
plate stiffness
Definition plate.cpp:58
static constexpr double delta
FTensor::Index< 'm', 3 > m
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition Core.hpp:82
Interface to cut meshes.
MoFEMErrorCode findEdgesToTrim(Range *fixed_edges, Range *corner_nodes, Tag th=NULL, const double tol=1e-4, const bool debug=false)
Find edges to trimEdges.
MoFEMErrorCode copySurface(const Range surface, Tag th=NULL, double *shift=NULL, double *origin=NULL, double *transform=NULL, const std::string save_mesh="")
copy surface entities
MoFEMErrorCode setFront(const Range surface)
MoFEMErrorCode findLevelSetVolumes(Tag th, Range &vol_edges, const bool remove_adj_prims_edges, int verb=QUIET, const bool debug=false, const std::string edges_file_name=string())
Find level set volumes.
MoFEMErrorCode moveMidNodesOnTrimmedEdges(Tag th=NULL)
move trimmed edges mid nodes
MoFEMErrorCode trimOnly(const BitRefLevel trim_bit, Tag th, const double tol_trim_close, Range *fixed_edges=NULL, Range *corner_nodes=NULL, const bool update_meshsets=false, const bool debug=false)
Trim mesh only.
MoFEMErrorCode trimSurface(Range *fixed_edge, Range *corner_nodes, const double tol=1e-4, Tag th=NULL, const bool debug=false)
Trim surface from faces beyond front.
double aveLength
Average edge length.
MoFEMErrorCode trimEdgesInTheMiddle(const BitRefLevel bit, const bool debug=false)
trim edges
MoFEMErrorCode setCoords(Tag th, const BitRefLevel bit=BitRefLevel(), const BitRefLevel mask=BitRefLevel().set())
set coords from tag
const Range & projectZeroDistanceEnts() const
MoFEMErrorCode findEdgesToCut(Range vol, int verb=QUIET, const bool debug=false)
find edges to cut
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
MoFEMErrorCode cutEdgesInMiddle(const BitRefLevel bit, Range &cut_vols, Range &cut_surf, Range &cut_verts, const bool debug=false)
cut edges
MoFEMErrorCode setSurface(const Range surface)
set surface entities
map< EntityHandle, TreeData > edgesToTrim
map< EntityHandle, TreeData > verticesOnTrimEdges
MoFEMErrorCode buildTree()
build tree
double maxLength
Maximal edge length.
MoFEMErrorCode splitSides(const BitRefLevel split_bit, const BitRefLevel bit, const Range &ents, Tag th=NULL)
split sides
MoFEMErrorCode mergeBadEdges(const int fraction_level, const Range &tets, const Range &surface, const Range &fixed_edges, const Range &corner_nodes, Range &merged_nodes, Range &out_tets, Range &new_surf, Tag th, const bool update_meshsets=false, const BitRefLevel *bit_ptr=NULL, const bool debug=false)
Merge edges.
MoFEMErrorCode mergeVolumes(const Range volume)
merge volume entities
MoFEMErrorCode moveMidNodesOnCutEdges(Tag th=NULL)
projecting of mid edge nodes on new mesh on surface
map< EntityHandle, TreeData > verticesOnCutEdges
multi_index_container< LengthMapData, indexed_by< ordered_non_unique< member< LengthMapData, double, &LengthMapData::lEngth > >, hashed_unique< member< LengthMapData, EntityHandle, &LengthMapData::eDge > >, ordered_non_unique< member< LengthMapData, double, &LengthMapData::qUality > > > > LengthMapData_multi_index
MoFEMErrorCode saveCutEdges(const std::string prefix="")
MoFEMErrorCode snapSurfaceSkinToEdges(const Range fixed_edges, const double rel_tol, const double abs_tol, Tag th=nullptr, const bool debug=false)
MoFEMErrorCode removePathologicalFrontTris(const BitRefLevel split_bit, Range &ents)
Remove pathological elements on surface internal front.
MoFEMErrorCode setConstrainSurface(const Range surf)
Set the constrain surface object.
MoFEMErrorCode cutAndTrim(int &first_bit, Tag th, const double tol_geometry, const double tol_cut_close, const double tol_trim_close, Range *fixed_edges=NULL, Range *corner_nodes=NULL, const bool update_meshsets=false, const bool debug=false)
Cut and trim.
const Range & getNewTrimSurfaces() const
MoFEMErrorCode cutOnly(Range vol, const BitRefLevel cut_bit, Tag th, const double tol_geometry, const double tol_cut_close, Range *fixed_edges=NULL, Range *corner_nodes=NULL, const bool update_meshsets=false, const bool debug=false)
Cut mesh only.
MoFEMErrorCode setVolume(const Range volume)
set volume entities
MoFEMErrorCode createFrontLevelSets(Range vol, Tag th=nullptr, int verb=QUIET, const bool debug=false)
Calculate distance from mesh nodes to surface front.
MoFEMErrorCode setTagData(Tag th, const BitRefLevel bit=BitRefLevel())
set coords to tag
CutMeshInterface(const MoFEM::Core &core)
MoFEMErrorCode cutTrimAndMerge(int &first_bit, const int fraction_level, Tag th, const double tol_geometry, const double tol_cut_close, const double tol_trim_close, Range &fixed_edges, Range &corner_nodes, const bool update_meshsets=false, const bool debug=false)
Cut, trim and merge.
map< EntityHandle, TreeData > edgesToCut
MoFEMErrorCode makeFront(const bool debug=false)
Create front from the surface.
boost::shared_ptr< OrientedBoxTreeTool > treeSurfPtr
MoFEMErrorCode mergeSurface(const Range surface)
merge surface entities
MoFEMErrorCode refineMesh(const int init_bit_level, const int surf_levels, const int front_levels, Range *fixed_edges=nullptr, int verb=QUIET, const bool debug=false)
Refine and set level sets.
MoFEMErrorCode snapSurfaceToEdges(const Range surface_edges, const Range fixed_edges, const double rel_tol, const double abs_tol, Tag th=nullptr, const bool debug=false)
MoFEMErrorCode createSurfaceLevelSets(int verb=QUIET, const bool debug=false)
Calculate distance from mesh nodes to cut surface.
Mesh refinement interface.
MoFEMErrorCode refineTets(const EntityHandle meshset, const BitRefLevel &bit, int verb=QUIET, const bool debug=false)
refine TET in the meshset
MoFEMErrorCode addVerticesInTheMiddleOfEdges(const EntityHandle meshset, const BitRefLevel &bit, const bool recursive=false, int verb=QUIET, EntityHandle start_v=0)
make vertices in the middle of edges in meshset and add them to refinement levels defined by bit
Interface for managing meshsets containing materials and boundary conditions.
Merge node by collapsing edge between them.
multi_index_container< ParentChild, indexed_by< ordered_unique< member< ParentChild, EntityHandle, &ParentChild::pArent > >, ordered_non_unique< member< ParentChild, EntityHandle, &ParentChild::cHild > > > > ParentChildMap
bool getSuccessMerge()
Return true if successful merge.
MoFEMErrorCode mergeNodes(EntityHandle father, EntityHandle mother, Range &out_tets, Range *tets_ptr=NULL, const bool only_if_improve_quality=false, const double move=0, const int line_search=0, Tag th=NULL, const int verb=0)
merge nodes which sharing edge
ParentChildMap & getParentChildMap()
Get map of parent cand child.
Create interface from given surface and insert flat prisms in-between.
MoFEMErrorCode findFacesWithThreeNodesOnInternalSurfaceSkin(const EntityHandle sideset, const BitRefLevel mesh_bit_level, const bool recursive, Range &faces_with_three_nodes_on_front, int verb=QUIET)
Find triangles which have three nodes on internal surface skin.
MoFEMErrorCode getSides(const int msId, const CubitBCType cubit_bc_type, const BitRefLevel mesh_bit_level, const bool recursive, int verb=QUIET)
Store tetrahedra from each side of the interface separately in two child meshsets of the parent meshs...
MoFEMErrorCode splitSides(const EntityHandle meshset, const BitRefLevel &bit, const int msId, const CubitBCType cubit_bc_type, const bool add_interface_entities, const bool recursive=false, int verb=QUIET)
Split nodes and other entities of tetrahedra on both sides of the interface and insert flat prisms in...
Auxiliary tools.
Definition Tools.hpp:19
static SEGMENT_MIN_DISTANCE minDistanceFromSegments(const double *w_ptr, const double *v_ptr, const double *k_ptr, const double *l_ptr, double *const tvw_ptr=nullptr, double *const tlk_ptr=nullptr)
Find points on two segments in closest distance.
Definition Tools.cpp:472
static double volumeLengthQuality(const double *coords)
Calculate tetrahedron volume length quality.
Definition Tools.cpp:15
base class for all interface classes
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.