meshTools.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) 2015-2020 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 "meshTools.H"
30#include "polyMesh.H"
31#include "hexMatcher.H"
32#include "treeBoundBox.H"
33
34// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35
37(
38 const vector& n,
40 const labelList& faceLabels
41)
42{
43 forAll(faceLabels, i)
44 {
45 if ((faceNormals[faceLabels[i]] & n) < SMALL)
46 {
47 // Found normal in different direction from n.
48 return false;
49 }
50 }
51
52 return true;
53}
54
55
57{
58 vectorField octantNormal(8);
59 octantNormal[mXmYmZ] = vector(-1, -1, -1);
60 octantNormal[pXmYmZ] = vector(1, -1, -1);
61 octantNormal[mXpYmZ] = vector(-1, 1, -1);
62 octantNormal[pXpYmZ] = vector(1, 1, -1);
63
64 octantNormal[mXmYpZ] = vector(-1, -1, 1);
65 octantNormal[pXmYpZ] = vector(1, -1, 1);
66 octantNormal[mXpYpZ] = vector(-1, 1, 1);
67 octantNormal[pXpYpZ] = vector(1, 1, 1);
68
69 octantNormal /= mag(octantNormal);
70
71
72 vectorField pn(pp.nPoints());
73
75 const vectorField& pointNormals = pp.pointNormals();
76 const labelListList& pointFaces = pp.pointFaces();
77
78 forAll(pointFaces, pointi)
79 {
80 const labelList& pFaces = pointFaces[pointi];
81
82 if (visNormal(pointNormals[pointi], faceNormals, pFaces))
83 {
84 pn[pointi] = pointNormals[pointi];
85 }
86 else
87 {
89 << "Average point normal not visible for point:"
90 << pp.meshPoints()[pointi] << endl;
91
92 label visOctant =
93 mXmYmZMask
94 | pXmYmZMask
95 | mXpYmZMask
96 | pXpYmZMask
97 | mXmYpZMask
98 | pXmYpZMask
99 | mXpYpZMask
100 | pXpYpZMask;
101
102 forAll(pFaces, i)
103 {
104 const vector& n = faceNormals[pFaces[i]];
105
106 if (n.x() > SMALL)
107 {
108 // All -x octants become invisible
109 visOctant &= ~mXmYmZMask;
110 visOctant &= ~mXmYpZMask;
111 visOctant &= ~mXpYmZMask;
112 visOctant &= ~mXpYpZMask;
113 }
114 else if (n.x() < -SMALL)
115 {
116 // All +x octants become invisible
117 visOctant &= ~pXmYmZMask;
118 visOctant &= ~pXmYpZMask;
119 visOctant &= ~pXpYmZMask;
120 visOctant &= ~pXpYpZMask;
121 }
122
123 if (n.y() > SMALL)
124 {
125 visOctant &= ~mXmYmZMask;
126 visOctant &= ~mXmYpZMask;
127 visOctant &= ~pXmYmZMask;
128 visOctant &= ~pXmYpZMask;
129 }
130 else if (n.y() < -SMALL)
131 {
132 visOctant &= ~mXpYmZMask;
133 visOctant &= ~mXpYpZMask;
134 visOctant &= ~pXpYmZMask;
135 visOctant &= ~pXpYpZMask;
136 }
137
138 if (n.z() > SMALL)
139 {
140 visOctant &= ~mXmYmZMask;
141 visOctant &= ~mXpYmZMask;
142 visOctant &= ~pXmYmZMask;
143 visOctant &= ~pXpYmZMask;
144 }
145 else if (n.z() < -SMALL)
146 {
147 visOctant &= ~mXmYpZMask;
148 visOctant &= ~mXpYpZMask;
149 visOctant &= ~pXmYpZMask;
150 visOctant &= ~pXpYpZMask;
151 }
152 }
153
154 label visI = -1;
155
156 label mask = 1;
157
158 for (label octant = 0; octant < 8; octant++)
159 {
160 if (visOctant & mask)
161 {
162 // Take first visible octant found. Note: could use
163 // combination of visible octants here.
164 visI = octant;
165
166 break;
167 }
168 mask <<= 1;
169 }
170
171 if (visI != -1)
172 {
173 // Take a vector in this octant.
174 pn[pointi] = octantNormal[visI];
175 }
176 else
177 {
178 pn[pointi] = Zero;
179
181 << "No visible octant for point:" << pp.meshPoints()[pointi]
182 << " cooord:" << pp.points()[pp.meshPoints()[pointi]] << nl
183 << "Normal set to " << pn[pointi] << endl;
184 }
185 }
186 }
187
188 return pn;
189}
190
191
193(
194 const primitiveMesh& mesh,
195 const label edgeI
196)
197{
198 return mesh.edges()[edgeI].unitVec(mesh.points());
199}
200
201
203(
204 Ostream& os,
205 const point& pt
206)
207{
208 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
209}
210
211
213(
214 Ostream& os,
215 const UList<point>& pts
216)
217{
218 forAll(pts, i)
219 {
220 const point& pt = pts[i];
221 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
222 }
223
224}
225
226
228(
229 Ostream& os,
230 const triad& t,
231 const point& origin
232)
233{
234 forAll(t, dirI)
235 {
236 writeOBJ(os, origin, origin + t[dirI]);
237 }
238}
239
240
242(
243 Ostream& os,
244 const point& p1,
245 const point& p2,
246 label& count
247)
248{
249 os << "v " << p1.x() << ' ' << p1.y() << ' ' << p1.z() << nl;
250 os << "v " << p2.x() << ' ' << p2.y() << ' ' << p2.z() << nl;
251 os << "l " << (count + 1) << " " << (count + 2) << endl;
252
253 count += 2;
254}
255
256
258(
259 Ostream& os,
260 const point& p1,
261 const point& p2
262)
263{
264 os << "v " << p1.x() << ' ' << p1.y() << ' ' << p1.z() << nl;
265 os << "vn "
266 << (p2.x() - p1.x()) << ' '
267 << (p2.y() - p1.y()) << ' '
268 << (p2.z() - p1.z()) << endl;
269}
270
271
273(
274 Ostream& os,
275 const treeBoundBox& bb
276)
277{
278 writeOBJ(os, bb.points());
279
280 for (const edge& e : treeBoundBox::edges)
281 {
282 os << "l " << (e[0] + 1) << ' ' << (e[1] + 1) << nl;
283 }
284}
285
286
288(
289 Ostream& os,
290 const cellList& cells,
291 const faceList& faces,
292 const UList<point>& points,
293 const labelList& cellLabels
294)
295{
296 labelHashSet usedFaces(4*cellLabels.size());
297
298 for (const label celli : cellLabels)
299 {
300 usedFaces.insert(cells[celli]);
301 }
302
303 writeOBJ(os, faces, points, usedFaces.toc());
304}
305
306
308(
309 const primitiveMesh& mesh,
310 const label celli,
311 const label edgeI
312)
313{
314 return mesh.edgeCells(edgeI).found(celli);
315}
316
317
319(
320 const primitiveMesh& mesh,
321 const label facei,
322 const label edgeI
323)
324{
325 return mesh.faceEdges(facei).found(edgeI);
326}
327
328
330(
331 const primitiveMesh& mesh,
332 const label celli,
333 const label facei
334)
335{
336 if (mesh.isInternalFace(facei))
337 {
338 if
339 (
340 (mesh.faceOwner()[facei] == celli)
341 || (mesh.faceNeighbour()[facei] == celli)
342 )
343 {
344 return true;
345 }
346 }
347 else
348 {
349 if (mesh.faceOwner()[facei] == celli)
350 {
351 return true;
352 }
353 }
354 return false;
355}
356
357
359(
360 const edgeList& edges,
361 const labelList& candidates,
362 const label v0,
363 const label v1
364)
365{
366 forAll(candidates, i)
367 {
368 const label edgeI = candidates[i];
369
370 const edge& e = edges[edgeI];
371
372 if ((e[0] == v0 && e[1] == v1) || (e[0] == v1 && e[1] == v0))
373 {
374 return edgeI;
375 }
376 }
377 return -1;
378}
379
380
382(
383 const primitiveMesh& mesh,
384 const label v0,
385 const label v1
386)
387{
388 const edgeList& edges = mesh.edges();
389
390 const labelList& v0Edges = mesh.pointEdges()[v0];
391
392 forAll(v0Edges, i)
393 {
394 label edgeI = v0Edges[i];
395
396 const edge& e = edges[edgeI];
397
398 if ((e.start() == v1) || (e.end() == v1))
399 {
400 return edgeI;
401 }
402 }
403 return -1;
404}
405
406
408(
409 const primitiveMesh& mesh,
410 const label f0,
411 const label f1
412)
413{
414 const labelList& f0Edges = mesh.faceEdges()[f0];
415 const labelList& f1Edges = mesh.faceEdges()[f1];
416
417 forAll(f0Edges, f0EdgeI)
418 {
419 label edge0 = f0Edges[f0EdgeI];
420
421 forAll(f1Edges, f1EdgeI)
422 {
423 label edge1 = f1Edges[f1EdgeI];
424
425 if (edge0 == edge1)
426 {
427 return edge0;
428 }
429 }
430 }
432 << "Faces " << f0 << " and " << f1 << " do not share an edge"
433 << abort(FatalError);
434
435 return -1;
436
437}
438
439
441(
442 const primitiveMesh& mesh,
443 const label cell0I,
444 const label cell1I
445)
446{
447 const cell& cFaces = mesh.cells()[cell0I];
448
449 forAll(cFaces, cFacei)
450 {
451 label facei = cFaces[cFacei];
452
453 if
454 (
455 mesh.isInternalFace(facei)
456 && (
457 mesh.faceOwner()[facei] == cell1I
458 || mesh.faceNeighbour()[facei] == cell1I
459 )
460 )
461 {
462 return facei;
463 }
464 }
465
466
468 << "No common face for"
469 << " cell0I:" << cell0I << " faces:" << cFaces
470 << " cell1I:" << cell1I << " faces:"
471 << mesh.cells()[cell1I]
472 << abort(FatalError);
473
474 return -1;
475}
476
477
479(
480 const primitiveMesh& mesh,
481 const label celli,
482 const label edgeI,
483 label& face0,
484 label& face1
485)
486{
487 const labelList& eFaces = mesh.edgeFaces(edgeI);
488
489 face0 = -1;
490 face1 = -1;
491
492 forAll(eFaces, eFacei)
493 {
494 label facei = eFaces[eFacei];
495
496 if (faceOnCell(mesh, celli, facei))
497 {
498 if (face0 == -1)
499 {
500 face0 = facei;
501 }
502 else
503 {
504 face1 = facei;
505
506 return;
507 }
508 }
509 }
510
511 if ((face0 == -1) || (face1 == -1))
512 {
514 << "Can not find faces using edge " << mesh.edges()[edgeI]
515 << " on cell " << celli << abort(FatalError);
516 }
517}
518
519
521(
522 const primitiveMesh& mesh,
523 const labelList& edgeLabels,
524 const label thisEdgeI,
525 const label thisVertI
526)
527{
528 forAll(edgeLabels, edgeLabelI)
529 {
530 label edgeI = edgeLabels[edgeLabelI];
531
532 if (edgeI != thisEdgeI)
533 {
534 const edge& e = mesh.edges()[edgeI];
535
536 if ((e.start() == thisVertI) || (e.end() == thisVertI))
537 {
538 return edgeI;
539 }
540 }
541 }
542
544 << "Can not find edge in "
545 << UIndirectList<edge>(mesh.edges(), edgeLabels)
546 << " connected to edge "
547 << thisEdgeI << " with vertices " << mesh.edges()[thisEdgeI]
548 << " on side " << thisVertI << abort(FatalError);
549
550 return -1;
551}
552
553
555(
556 const primitiveMesh& mesh,
557 const label celli,
558 const label facei,
559 const label edgeI
560)
561{
562 label face0;
563 label face1;
564
565 getEdgeFaces(mesh, celli, edgeI, face0, face1);
566
567 if (face0 == facei)
568 {
569 return face1;
570 }
571 else
572 {
573 return face0;
574 }
575}
576
577
579(
580 const primitiveMesh& mesh,
581 const label otherCelli,
582 const label facei
583)
584{
585 if (!mesh.isInternalFace(facei))
586 {
588 << "Face " << facei << " is not internal"
589 << abort(FatalError);
590 }
591
592 label newCelli = mesh.faceOwner()[facei];
593
594 if (newCelli == otherCelli)
595 {
596 newCelli = mesh.faceNeighbour()[facei];
597 }
598 return newCelli;
599}
600
601
603(
604 const primitiveMesh& mesh,
605 const label facei,
606 const label startEdgeI,
607 const label startVertI,
608 const label nEdges
609)
610{
611 const labelList& fEdges = mesh.faceEdges(facei);
612
613 label edgeI = startEdgeI;
614
615 label vertI = startVertI;
616
617 for (label iter = 0; iter < nEdges; iter++)
618 {
619 edgeI = otherEdge(mesh, fEdges, edgeI, vertI);
620
621 vertI = mesh.edges()[edgeI].otherVertex(vertI);
622 }
623
624 return edgeI;
625}
626
627
629(
630 const polyMesh& mesh,
631 point& pt
632)
633{
634 const Vector<label>& dirs = mesh.geometricD();
635
636 const point& min = mesh.bounds().min();
637 const point& max = mesh.bounds().max();
638
639 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
640 {
641 if (dirs[cmpt] == -1)
642 {
643 pt[cmpt] = 0.5*(min[cmpt] + max[cmpt]);
644 }
645 }
646}
647
648
650(
651 const polyMesh& mesh,
652 pointField& pts
653)
654{
655 const Vector<label>& dirs = mesh.geometricD();
656
657 const point& min = mesh.bounds().min();
658 const point& max = mesh.bounds().max();
659
660 bool isConstrained = false;
661 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
662 {
663 if (dirs[cmpt] == -1)
664 {
665 isConstrained = true;
666 break;
667 }
668 }
669
670 if (isConstrained)
671 {
672 forAll(pts, i)
673 {
674 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
675 {
676 if (dirs[cmpt] == -1)
677 {
678 pts[i][cmpt] = 0.5*(min[cmpt] + max[cmpt]);
679 }
680 }
681 }
682 }
683}
684
685
687(
688 const polyMesh& mesh,
689 const Vector<label>& dirs,
690 vector& d
691)
692{
693 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
694 {
695 if (dirs[cmpt] == -1)
696 {
697 d[cmpt] = 0.0;
698 }
699 }
700}
701
702
704(
705 const polyMesh& mesh,
706 const Vector<label>& dirs,
707 vectorField& d
708)
709{
710 bool isConstrained = false;
711 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
712 {
713 if (dirs[cmpt] == -1)
714 {
715 isConstrained = true;
716 break;
717 }
718 }
719
720 if (isConstrained)
721 {
722 forAll(d, i)
723 {
724 for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
725 {
726 if (dirs[cmpt] == -1)
727 {
728 d[i][cmpt] = 0.0;
729 }
730 }
731 }
732 }
733}
734
735
737(
738 const primitiveMesh& mesh,
739 const label celli,
740 const label e0,
741 label& e1,
742 label& e2,
743 label& e3
744)
745{
746 // Go to any face using e0
747 label facei = meshTools::otherFace(mesh, celli, -1, e0);
748
749 // Opposite edge on face
750 e1 = meshTools::walkFace(mesh, facei, e0, mesh.edges()[e0].end(), 2);
751
752 facei = meshTools::otherFace(mesh, celli, facei, e1);
753
754 e2 = meshTools::walkFace(mesh, facei, e1, mesh.edges()[e1].end(), 2);
755
756 facei = meshTools::otherFace(mesh, celli, facei, e2);
757
758 e3 = meshTools::walkFace(mesh, facei, e2, mesh.edges()[e2].end(), 2);
759}
760
761
763(
764 const primitiveMesh& mesh,
765 const label celli,
766 const label startEdgeI
767)
768{
769 if (!hexMatcher::test(mesh, celli))
770 {
772 << "Not a hex : cell:" << celli << abort(FatalError);
773 }
774
775
776 vector avgVec(normEdgeVec(mesh, startEdgeI));
777
778 label edgeI = startEdgeI;
779
780 label facei = -1;
781
782 for (label i = 0; i < 3; i++)
783 {
784 // Step to next face, next edge
785 facei = meshTools::otherFace(mesh, celli, facei, edgeI);
786
787 vector eVec(normEdgeVec(mesh, edgeI));
788
789 if ((eVec & avgVec) > 0)
790 {
791 avgVec += eVec;
792 }
793 else
794 {
795 avgVec -= eVec;
796 }
797
798 label vertI = mesh.edges()[edgeI].end();
799
800 edgeI = meshTools::walkFace(mesh, facei, edgeI, vertI, 2);
801 }
802
803 avgVec.normalise();
804
805 return avgVec;
806}
807
808
810(
811 const primitiveMesh& mesh,
812 const label celli,
813 const vector& cutDir
814)
815{
816 if (!hexMatcher::test(mesh, celli))
817 {
819 << "Not a hex : cell:" << celli << abort(FatalError);
820 }
821
822 const labelList& cEdges = mesh.cellEdges()[celli];
823
824 labelHashSet doneEdges(2*cEdges.size());
825
826 scalar maxCos = -GREAT;
827 label maxEdgeI = -1;
828
829 for (label i = 0; i < 4; i++)
830 {
831 forAll(cEdges, cEdgeI)
832 {
833 label e0 = cEdges[cEdgeI];
834
835 if (!doneEdges.found(e0))
836 {
837 vector avgDir(edgeToCutDir(mesh, celli, e0));
838
839 scalar cosAngle = mag(avgDir & cutDir);
840
841 if (cosAngle > maxCos)
842 {
843 maxCos = cosAngle;
844 maxEdgeI = e0;
845 }
846
847 // Mark off edges in cEdges.
848 label e1, e2, e3;
849 getParallelEdges(mesh, celli, e0, e1, e2, e3);
850
851 doneEdges.insert(e0);
852 doneEdges.insert(e1);
853 doneEdges.insert(e2);
854 doneEdges.insert(e3);
855 }
856 }
857 }
858
859 forAll(cEdges, cEdgeI)
860 {
861 if (!doneEdges.found(cEdges[cEdgeI]))
862 {
864 << "Cell:" << celli << " edges:" << cEdges << endl
865 << "Edge:" << cEdges[cEdgeI] << " not yet handled"
866 << abort(FatalError);
867 }
868 }
869
870 if (maxEdgeI == -1)
871 {
873 << "Problem : did not find edge aligned with " << cutDir
874 << " on cell " << celli << abort(FatalError);
875 }
876
877 return maxEdgeI;
878}
879
880
881// ************************************************************************* //
label n
Y[inertIndex] max(0.0)
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition: HashTable.C:122
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of faces which address into the list of points.
label nPoints() const
Number of points supporting patch faces.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const Field< point_type > & pointNormals() const
Return point normals for patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
const labelListList & pointFaces() const
Return point-face addressing.
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
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:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition: VectorI.H:123
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
tmp< pointField > points() const
Vertex coordinates. In octant coding.
Definition: treeBoundBox.C:92
Representation of a 3D Cartesian coordinate system as a Vector of row vectors.
Definition: triad.H:67
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
const pointField & points
const cellShapeList & cells
surfaceVectorField faceNormals(mesh.Sf()/mesh.magSf())
#define WarningInFunction
Report a warning using Foam::Warning.
bool faceOnCell(const primitiveMesh &mesh, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:330
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:687
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:629
bool edgeOnCell(const primitiveMesh &mesh, const label celli, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:308
label otherCell(const primitiveMesh &mesh, const label celli, const label facei)
Return cell on other side of face. Throws error.
Definition: meshTools.C:579
vectorField calcBoxPointNormals(const primitivePatch &pp)
Calculate point normals on a 'box' mesh (all edges aligned with.
Definition: meshTools.C:56
label walkFace(const primitiveMesh &mesh, const label facei, const label startEdgeI, const label startVertI, const label nEdges)
Returns label of edge nEdges away from startEdge (in the direction.
Definition: meshTools.C:603
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:359
label otherEdge(const primitiveMesh &mesh, const labelList &edgeLabels, const label thisEdgeI, const label thisVertI)
Return label of other edge (out of candidates edgeLabels)
Definition: meshTools.C:521
label getSharedEdge(const primitiveMesh &mesh, const label f0, const label f1)
Return edge shared by two faces. Throws error if no edge found.
Definition: meshTools.C:408
label getSharedFace(const primitiveMesh &mesh, const label cell0, const label cell1)
Return face shared by two cells. Throws error if none found.
Definition: meshTools.C:441
vector normEdgeVec(const primitiveMesh &, const label edgeI)
Normalized edge vector.
Definition: meshTools.C:193
label otherFace(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:555
bool visNormal(const vector &n, const vectorField &faceNormals, const labelList &faceLabels)
Check if n is in same direction as normals of all faceLabels.
Definition: meshTools.C:37
vector edgeToCutDir(const primitiveMesh &mesh, const label celli, const label edgeI)
Given edge on hex find all 'parallel' (i.e. non-connected)
Definition: meshTools.C:763
bool edgeOnFace(const primitiveMesh &mesh, const label facei, const label edgeI)
Is edge used by face.
Definition: meshTools.C:319
void getParallelEdges(const primitiveMesh &mesh, const label celli, const label e0, label &e1, label &e2, label &e3)
Given edge on hex find other 'parallel', non-connected edges.
Definition: meshTools.C:737
void getEdgeFaces(const primitiveMesh &mesh, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:479
label cutDirToEdge(const primitiveMesh &mesh, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:810
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
error FatalError
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
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];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } 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};const 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< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333