Loading...
Searching...
No Matches
removePoints.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "BiIndirectList.H"
30#include "removePoints.H"
31#include "PstreamReduceOps.H"
32#include "polyMesh.H"
33#include "polyTopoChange.H"
34#include "polyRemovePoint.H"
35#include "polyAddPoint.H"
36#include "polyModifyFace.H"
37#include "syncTools.H"
38#include "faceSet.H"
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
43namespace Foam
45
47
48// Combine-reduce operator to combine data on faces. Takes care
49// of reverse orientation on coupled face.
50template<class T, template<class> class CombineOp>
51class faceEqOp
52{
53public:
54
55 void operator()(List<T>& x, const UList<T>& y) const
56 {
57 if (y.size() > 0)
58 {
59 if (x.size() == 0)
60 {
61 x = y;
62 }
63 else
64 {
65 label j = 0;
66 forAll(x, i)
67 {
68 CombineOp<T>()(x[i], y[j]);
69 j = y.rcIndex(j);
70 }
71 }
72 }
73 }
74};
75
76}
77
78
79// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
80
81void Foam::removePoints::modifyFace
82(
83 const label facei,
84 const face& newFace,
85 polyTopoChange& meshMod
86) const
87{
88 // Get other face data.
89 label patchi = -1;
90 label owner = mesh_.faceOwner()[facei];
91 label neighbour = -1;
92
93 if (mesh_.isInternalFace(facei))
94 {
95 neighbour = mesh_.faceNeighbour()[facei];
96 }
97 else
98 {
99 patchi = mesh_.boundaryMesh().whichPatch(facei);
100 }
101
102 label zoneID = mesh_.faceZones().whichZone(facei);
103
104 bool zoneFlip = false;
105
106 if (zoneID >= 0)
107 {
108 const faceZone& fZone = mesh_.faceZones()[zoneID];
109
110 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
111 }
112
113 meshMod.setAction
114 (
115 polyModifyFace
116 (
117 newFace, // modified face
118 facei, // label of face being modified
119 owner, // owner
120 neighbour, // neighbour
121 false, // face flip
122 patchi, // patch for face
123 false, // remove from zone
124 zoneID, // zone for face
125 zoneFlip // face flip in zone
127 );
128}
129
130
131// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
132
133Foam::removePoints::removePoints
134(
135 const polyMesh& mesh,
136 const bool undoable
137)
138:
139 mesh_(mesh),
140 undoable_(undoable)
141{}
142
143
144// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145
147(
148 const scalar minCos,
149 boolList& pointCanBeDeleted
150) const
151{
152 // Containers to store two edges per point:
153 // -1 : not filled
154 // >= 0 : edge label
155 // -2 : more than two edges for point
156 labelList edge0(mesh_.nPoints(), -1);
157 labelList edge1(mesh_.nPoints(), -1);
158
159 const edgeList& edges = mesh_.edges();
160
161 forAll(edges, edgeI)
162 {
163 const edge& e = edges[edgeI];
164
165 forAll(e, eI)
166 {
167 label pointi = e[eI];
168
169 if (edge0[pointi] == -2)
170 {
171 // Already too many edges
172 }
173 else if (edge0[pointi] == -1)
174 {
175 // Store first edge using point
176 edge0[pointi] = edgeI;
177 }
178 else
179 {
180 // Already one edge using point. Check second container.
181
182 if (edge1[pointi] == -1)
183 {
184 // Store second edge using point
185 edge1[pointi] = edgeI;
186 }
187 else
188 {
189 // Third edge using point. Mark.
190 edge0[pointi] = -2;
191 edge1[pointi] = -2;
192 }
193 }
194 }
195 }
196
197
198 // Check the ones used by only 2 edges that these are sufficiently in line.
199 const pointField& points = mesh_.points();
200
201 pointCanBeDeleted.setSize(mesh_.nPoints());
202 pointCanBeDeleted = false;
203 //label nDeleted = 0;
204
205 forAll(edge0, pointi)
206 {
207 if (edge0[pointi] >= 0 && edge1[pointi] >= 0)
208 {
209 // Point used by two edges exactly
210
211 const edge& e0 = edges[edge0[pointi]];
212 const edge& e1 = edges[edge1[pointi]];
213
214 label common = e0.commonVertex(e1);
215 label vLeft = e0.otherVertex(common);
216 label vRight = e1.otherVertex(common);
217
218 const vector e0Vec = normalised(points[common] - points[vLeft]);
219 const vector e1Vec = normalised(points[vRight] - points[common]);
220
221 if ((e0Vec & e1Vec) > minCos)
222 {
223 pointCanBeDeleted[pointi] = true;
224 //nDeleted++;
225 }
226 }
227 else if (edge0[pointi] == -1)
228 {
229 // point not used at all
230 pointCanBeDeleted[pointi] = true;
231 //nDeleted++;
232 }
233 }
234 edge0.clear();
235 edge1.clear();
236
237
238 // Protect any points on faces that would collapse down to nothing
239 // No particular intelligence so might protect too many points
240 forAll(mesh_.faces(), facei)
241 {
242 const face& f = mesh_.faces()[facei];
243
244 label nCollapse = 0;
245 forAll(f, fp)
246 {
247 if (pointCanBeDeleted[f[fp]])
248 {
249 nCollapse++;
250 }
251 }
252
253 if ((f.size() - nCollapse) < 3)
254 {
255 // Just unmark enough points
256 forAll(f, fp)
257 {
258 if (pointCanBeDeleted[f[fp]])
259 {
260 pointCanBeDeleted[f[fp]] = false;
261 --nCollapse;
262 if (nCollapse == 0)
263 {
264 break;
265 }
266 }
267 }
268 }
269 }
270
271
272 // Point can be deleted only if all processors want to delete it
274 (
275 mesh_,
276 pointCanBeDeleted,
278 true // null value
279 );
280
281 label nDeleted = 0;
282 forAll(pointCanBeDeleted, pointi)
283 {
284 if (pointCanBeDeleted[pointi])
285 {
286 nDeleted++;
288 }
289
290 return returnReduce(nDeleted, sumOp<label>());
291}
292
293
295(
296 const boolList& pointCanBeDeleted,
297 polyTopoChange& meshMod
298)
299{
300 // Count deleted points
301 label nDeleted = 0;
302 forAll(pointCanBeDeleted, pointi)
303 {
304 if (pointCanBeDeleted[pointi])
305 {
306 nDeleted++;
307 }
308 }
309
310 // Faces (in mesh face labels) affected by points removed. Will hopefully
311 // be only a few.
312 labelHashSet facesAffected(4*nDeleted);
313
314
315 // Undo: from global mesh point to index in savedPoint_
316 Map<label> pointToSaved;
317
318 // Size undo storage
319 if (undoable_)
320 {
321 savedPoints_.resize(nDeleted);
322 pointToSaved.reserve(nDeleted);
323 }
324
325
326 // Remove points
327 // ~~~~~~~~~~~~~
328
329 nDeleted = 0;
330
331 forAll(pointCanBeDeleted, pointi)
332 {
333 if (pointCanBeDeleted[pointi])
334 {
335 if (undoable_)
336 {
337 pointToSaved.insert(pointi, nDeleted);
338 savedPoints_[nDeleted++] = mesh_.points()[pointi];
339 }
340 meshMod.setAction(polyRemovePoint(pointi));
341
342 // Store faces affected
343 const labelList& pFaces = mesh_.pointFaces()[pointi];
344
345 facesAffected.insert(pFaces);
346 }
347 }
348
349
350
351 // Update faces
352 // ~~~~~~~~~~~~
353
354
355 if (undoable_)
356 {
357 savedFaceLabels_.setSize(facesAffected.size());
358 savedFaces_.setSize(facesAffected.size());
359 }
360 label nSaved = 0;
361
362 for (const label facei : facesAffected)
363 {
364 const face& f = mesh_.faces()[facei];
365
366 face newFace(f.size());
367
368 label newI = 0;
369
370 forAll(f, fp)
371 {
372 label pointi = f[fp];
373
374 if (!pointCanBeDeleted[pointi])
375 {
376 newFace[newI++] = pointi;
377 }
378 }
379 newFace.setSize(newI);
380
381 // Actually change the face to the new vertices
382 modifyFace(facei, newFace, meshMod);
383
384 // Save the face. Negative indices are into savedPoints_
385 if (undoable_)
386 {
387 savedFaceLabels_[nSaved] = facei;
388
389 face& savedFace = savedFaces_[nSaved++];
390 savedFace.setSize(f.size());
391
392 forAll(f, fp)
393 {
394 label pointi = f[fp];
395
396 if (pointCanBeDeleted[pointi])
397 {
398 savedFace[fp] = -pointToSaved[pointi]-1;
399 }
400 else
401 {
402 savedFace[fp] = pointi;
403 }
404 }
405 }
406 }
407
408 if (undoable_)
409 {
410 // DEBUG: Compare the stored faces with the current ones.
411 if (debug)
412 {
413 forAll(savedFaceLabels_, saveI)
414 {
415 // Points from the mesh
416 List<point> meshPoints
417 (
419 (
420 mesh_.points(),
421 mesh_.faces()[savedFaceLabels_[saveI]] // mesh face
422 )
423 );
424
425 // Points from the saved face
426 List<point> keptPoints
427 (
429 (
430 mesh_.points(),
431 savedPoints_,
432 savedFaces_[saveI] // saved face
433 )()
434 );
435
436 if (meshPoints != keptPoints)
437 {
439 << "facei:" << savedFaceLabels_[saveI] << nl
440 << "meshPoints:" << meshPoints << nl
441 << "keptPoints:" << keptPoints << nl
442 << abort(FatalError);
444 }
445 }
446 }
447}
448
449
451{
452 if (undoable_)
453 {
454 forAll(savedFaceLabels_, localI)
455 {
456 if (savedFaceLabels_[localI] >= 0)
457 {
458 label newFacei = map.reverseFaceMap()[savedFaceLabels_[localI]];
459
460 if (newFacei == -1)
461 {
463 << "Old face " << savedFaceLabels_[localI]
464 << " seems to have disappeared."
465 << abort(FatalError);
466 }
467 savedFaceLabels_[localI] = newFacei;
468 }
469 }
470
471
472 // Renumber mesh vertices (indices >=0). Leave saved vertices
473 // (<0) intact.
474 forAll(savedFaces_, i)
475 {
476 face& f = savedFaces_[i];
477
478 forAll(f, fp)
479 {
480 label pointi = f[fp];
481
482 if (pointi >= 0)
483 {
484 f[fp] = map.reversePointMap()[pointi];
485
486 if (f[fp] == -1)
487 {
489 << "Old point " << pointi
490 << " seems to have disappeared."
491 << abort(FatalError);
492 }
493 }
494 }
495 }
496
497
498 // DEBUG: Compare the stored faces with the current ones.
499 if (debug)
500 {
501 forAll(savedFaceLabels_, saveI)
502 {
503 if (savedFaceLabels_[saveI] >= 0)
504 {
505 const face& f = mesh_.faces()[savedFaceLabels_[saveI]];
506
507 // Get kept points of saved faces.
508 const face& savedFace = savedFaces_[saveI];
509
510 face keptFace(savedFace.size());
511 label keptFp = 0;
512
513 forAll(savedFace, fp)
514 {
515 label pointi = savedFace[fp];
516
517 if (pointi >= 0)
518 {
519 keptFace[keptFp++] = pointi;
520 }
521 }
522 keptFace.setSize(keptFp);
523
524 // Compare as faces (since f might have rotated and
525 // face::operator== takes care of this)
526 if (keptFace != f)
527 {
529 << "facei:" << savedFaceLabels_[saveI] << nl
530 << "face:" << f << nl
531 << "keptFace:" << keptFace << nl
532 << "saved points:"
534 (
535 mesh_.points(),
536 savedPoints_,
537 savedFace
538 )() << nl
539 << abort(FatalError);
540 }
541 }
542 }
543 }
545}
546
547
548// Given list of faces to undo picks up the local indices of the faces
549// to restore. Additionally it also picks up all the faces that use
550// any of the deleted points.
552(
553 const labelList& undoFaces,
554 labelList& localFaces,
555 labelList& localPoints
556) const
557{
558 if (!undoable_)
559 {
561 << "removePoints not constructed with"
562 << " unrefinement capability."
563 << abort(FatalError);
564 }
565
566 if (debug)
567 {
568 // Check if synced.
569 faceSet undoFacesSet(mesh_, "undoFacesSet", undoFaces);
570 label sz = undoFacesSet.size();
571
572 undoFacesSet.sync(mesh_);
573 if (sz != undoFacesSet.size())
574 {
576 << "undoFaces not synchronised across coupled faces." << endl
577 << "Size before sync:" << sz
578 << " after sync:" << undoFacesSet.size()
579 << abort(FatalError);
580 }
581 }
582
583
584 // Problem: input undoFaces are synced. But e.g.
585 // two faces, A (uncoupled) and B(coupled) share a deleted point. A gets
586 // passed in to be restored. Now picking up the deleted point and the
587 // original faces using it picks up face B. But not its coupled neighbour!
588 // Problem is that we cannot easily synchronise the deleted points
589 // (e.g. syncPointList) since they're not in the mesh anymore - only the
590 // faces are. So we have to collect the points-to-restore as indices
591 // in the faces (which is information we can synchronise)
592
593
594
595 // Mark points-to-restore
596 labelHashSet localPointsSet(undoFaces.size());
597
598 {
599 // Create set of faces to be restored
600 labelHashSet undoFacesSet(undoFaces);
601
602 forAll(savedFaceLabels_, saveI)
603 {
604 if (savedFaceLabels_[saveI] < 0)
605 {
607 << "Illegal face label " << savedFaceLabels_[saveI]
608 << " at index " << saveI
609 << abort(FatalError);
610 }
611
612 if (undoFacesSet.found(savedFaceLabels_[saveI]))
613 {
614 const face& savedFace = savedFaces_[saveI];
615
616 forAll(savedFace, fp)
617 {
618 if (savedFace[fp] < 0)
619 {
620 label savedPointi = -savedFace[fp]-1;
621
622 if (savedPoints_[savedPointi] == vector::max)
623 {
625 << "Trying to restore point " << savedPointi
626 << " from mesh face " << savedFaceLabels_[saveI]
627 << " saved face:" << savedFace
628 << " which has already been undone."
629 << abort(FatalError);
630 }
631
632 localPointsSet.insert(savedPointi);
633 }
634 }
635 }
636 }
637
638
639 // Per boundary face, per index in face whether the point needs
640 // restoring. Note that this is over all saved faces, not just over
641 // the ones in undoFaces.
642
643 boolListList faceVertexRestore(mesh_.nBoundaryFaces());
644
645 // Populate with my local points-to-restore.
646 forAll(savedFaces_, saveI)
647 {
648 label bFacei = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
649
650 if (bFacei >= 0)
651 {
652 const face& savedFace = savedFaces_[saveI];
653
654 boolList& fRestore = faceVertexRestore[bFacei];
655
656 fRestore.setSize(savedFace.size());
657 fRestore = false;
658
659 forAll(savedFace, fp)
660 {
661 if (savedFace[fp] < 0)
662 {
663 label savedPointi = -savedFace[fp]-1;
664
665 if (localPointsSet.found(savedPointi))
666 {
667 fRestore[fp] = true;
668 }
669 }
670 }
671 }
672 }
673
675 (
676 mesh_,
677 faceVertexRestore,
678 faceEqOp<bool, orEqOp>(), // special operator to handle faces
679 Foam::dummyTransform() // no transformation
680 );
681
682 // So now if any of the points-to-restore is used by any coupled face
683 // anywhere the corresponding index in faceVertexRestore will be set.
684
685 // Now combine the localPointSet and the (sychronised)
686 // boundary-points-to-restore.
687
688 forAll(savedFaces_, saveI)
689 {
690 label bFacei = savedFaceLabels_[saveI] - mesh_.nInternalFaces();
691
692 if (bFacei >= 0)
693 {
694 const boolList& fRestore = faceVertexRestore[bFacei];
695
696 const face& savedFace = savedFaces_[saveI];
697
698 forAll(fRestore, fp)
699 {
700 // Does neighbour require point restored?
701 if (fRestore[fp])
702 {
703 if (savedFace[fp] >= 0)
704 {
706 << "Problem: on coupled face:"
707 << savedFaceLabels_[saveI]
708 << " fc:"
709 << mesh_.faceCentres()[savedFaceLabels_[saveI]]
710 << endl
711 << " my neighbour tries to restore the vertex"
712 << " at index " << fp
713 << " whereas my saved face:" << savedFace
714 << " does not indicate a deleted vertex"
715 << " at that position."
716 << abort(FatalError);
717 }
718
719 label savedPointi = -savedFace[fp]-1;
720
721 localPointsSet.insert(savedPointi);
722 }
723 }
724 }
725 }
726 }
727
728 localPoints = localPointsSet.toc();
729
730
731 // Collect all saved faces using any localPointsSet
732 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
733
734 labelHashSet localFacesSet(2*undoFaces.size());
735
736 forAll(savedFaces_, saveI)
737 {
738 const face& savedFace = savedFaces_[saveI];
739
740 forAll(savedFace, fp)
741 {
742 if (savedFace[fp] < 0)
743 {
744 label savedPointi = -savedFace[fp]-1;
745
746 if (localPointsSet.found(savedPointi))
747 {
748 localFacesSet.insert(saveI);
749 }
750 }
751 }
752 }
753 localFaces = localFacesSet.toc();
754
755
756 // Note that at this point the localFaces to restore can contain points
757 // that are not going to be restored! The localFaces though will
758 // be guaranteed to be all the ones affected by the restoration of the
759 // localPoints.
760}
761
762
764(
765 const labelList& localFaces,
766 const labelList& localPoints,
767 polyTopoChange& meshMod
768)
769{
770 if (!undoable_)
771 {
773 << "removePoints not constructed with"
774 << " unrefinement capability."
775 << abort(FatalError);
776 }
777
778
779 // Per savedPoint -1 or the restored point label
780 labelList addedPoints(savedPoints_.size(), -1);
781
782 forAll(localPoints, i)
783 {
784 label localI = localPoints[i];
785
786 if (savedPoints_[localI] == vector::max)
787 {
789 << "Saved point " << localI << " already restored!"
790 << abort(FatalError);
791 }
792
793 addedPoints[localI] = meshMod.setAction
794 (
796 (
797 savedPoints_[localI], // point
798 -1, // master point
799 -1, // zone for point
800 true // supports a cell
801 )
802 );
803
804 // Mark the restored points so they are not restored again.
805 savedPoints_[localI] = vector::max;
806 }
807
808 forAll(localFaces, i)
809 {
810 label saveI = localFaces[i];
811
812 // Modify indices into saved points (so < 0) to point to the
813 // added points.
814 face& savedFace = savedFaces_[saveI];
815
816 face newFace(savedFace.size());
817 label newFp = 0;
818
819 bool hasSavedPoints = false;
820
821 forAll(savedFace, fp)
822 {
823 if (savedFace[fp] < 0)
824 {
825 label addedPointi = addedPoints[-savedFace[fp]-1];
826
827 if (addedPointi != -1)
828 {
829 savedFace[fp] = addedPointi;
830 newFace[newFp++] = addedPointi;
831 }
832 else
833 {
834 hasSavedPoints = true;
835 }
836 }
837 else
838 {
839 newFace[newFp++] = savedFace[fp];
840 }
841 }
842 newFace.setSize(newFp);
843
844 modifyFace(savedFaceLabels_[saveI], newFace, meshMod);
845
846 if (!hasSavedPoints)
847 {
848 // Face fully restored. Mark for compaction later on
849 savedFaceLabels_[saveI] = -1;
850 savedFaces_[saveI].clear();
851 }
852 }
853
854
855 // Compact face labels
856 label newSaveI = 0;
857
858 forAll(savedFaceLabels_, saveI)
859 {
860 if (savedFaceLabels_[saveI] != -1)
861 {
862 if (newSaveI != saveI)
863 {
864 savedFaceLabels_[newSaveI] = savedFaceLabels_[saveI];
865 savedFaces_[newSaveI].transfer(savedFaces_[saveI]);
866 }
867 newSaveI++;
868 }
869 }
870
871 savedFaceLabels_.setSize(newSaveI);
872 savedFaces_.setSize(newSaveI);
873
874
875 // Check that all faces have been restored that use any restored points
876 if (debug)
877 {
878 forAll(savedFaceLabels_, saveI)
879 {
880 const face& savedFace = savedFaces_[saveI];
881
882 forAll(savedFace, fp)
883 {
884 if (savedFace[fp] < 0)
885 {
886 label addedPointi = addedPoints[-savedFace[fp]-1];
887
888 if (addedPointi != -1)
889 {
891 << "Face:" << savedFaceLabels_[saveI]
892 << " savedVerts:" << savedFace
893 << " uses restored point:" << -savedFace[fp]-1
894 << " with new pointlabel:" << addedPointi
895 << abort(FatalError);
896 }
897 }
898 }
899 }
900 }
901}
902
903
904// ************************************************************************* //
scalar y
Inter-processor communication reduction functions.
Indexes into negList (negative index) or posList (zero or positive index).
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
A List with indirect addressing. Like IndirectList but does not store addressing.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form max
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
label commonVertex(const edge &other) const
Return vertex common with other edge or -1 on failure.
Definition edgeI.H:158
label otherVertex(const label vertex) const
Given one vertex label, return the other one.
Definition edgeI.H:174
void operator()(List< T > &x, const UList< T > &y) const
A list of face labels.
Definition faceSet.H:50
virtual void sync(const polyMesh &mesh)
Sync faceSet across coupled patches.
Definition faceSet.C:152
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reverseFaceMap() const noexcept
Reverse face map.
const labelList & reversePointMap() const noexcept
Reverse point map.
Class containing data for point addition.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition polyMesh.C:1107
Class containing data for point removal.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Removes selected points from mesh and updates faces using these points.
label countPointUsage(const scalar minCos, boolList &pointCanBeDeleted) const
Mark in pointCanBeDeleted the points that can be deleted.
void setUnrefinement(const labelList &localFaces, const labelList &localPoints, polyTopoChange &)
Restore selected faces and vertices.
void setRefinement(const boolList &, polyTopoChange &)
Play commands into polyTopoChange to remove points. Gets.
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
void getUnrefimentSet(const labelList &undoFaces, labelList &localFaces, labelList &localPoints) const
Given set of faces to restore calculates a consistent set of.
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=UPstream::parRun())
Synchronize values on boundary faces only.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Dummy transform to be used with syncTools.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< List< bool > > boolListList
List of boolList.
Definition boolList.H:36
errorManip< error > abort(error &err)
Definition errorManip.H:139
List< bool > boolList
A List of bools.
Definition List.H:60
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299