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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "meshTools.H"
27 #include "polyMesh.H"
28 #include "hexMatcher.H"
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
33 (
34  const vector& n,
35  const vectorField& faceNormals,
36  const labelList& faceLabels
37 )
38 {
39  forAll(faceLabels, i)
40  {
41  if ((faceNormals[faceLabels[i]] & n) < SMALL)
42  {
43  // Found normal in different direction from n.
44  return false;
45  }
46  }
47 
48  return true;
49 }
50 
51 
53 {
54  vectorField octantNormal(8);
55  octantNormal[mXmYmZ] = vector(-1, -1, -1);
56  octantNormal[pXmYmZ] = vector(1, -1, -1);
57  octantNormal[mXpYmZ] = vector(-1, 1, -1);
58  octantNormal[pXpYmZ] = vector(1, 1, -1);
59 
60  octantNormal[mXmYpZ] = vector(-1, -1, 1);
61  octantNormal[pXmYpZ] = vector(1, -1, 1);
62  octantNormal[mXpYpZ] = vector(-1, 1, 1);
63  octantNormal[pXpYpZ] = vector(1, 1, 1);
64 
65  octantNormal /= mag(octantNormal);
66 
67 
68  vectorField pn(pp.nPoints());
69 
70  const vectorField& faceNormals = pp.faceNormals();
71  const vectorField& pointNormals = pp.pointNormals();
72  const labelListList& pointFaces = pp.pointFaces();
73 
74  forAll(pointFaces, pointI)
75  {
76  const labelList& pFaces = pointFaces[pointI];
77 
78  if (visNormal(pointNormals[pointI], faceNormals, pFaces))
79  {
80  pn[pointI] = pointNormals[pointI];
81  }
82  else
83  {
85  << "Average point normal not visible for point:"
86  << pp.meshPoints()[pointI] << endl;
87 
88  label visOctant =
90  | pXmYmZMask
91  | mXpYmZMask
92  | pXpYmZMask
93  | mXmYpZMask
94  | pXmYpZMask
95  | mXpYpZMask
96  | pXpYpZMask;
97 
98  forAll(pFaces, i)
99  {
100  const vector& n = faceNormals[pFaces[i]];
101 
102  if (n.x() > SMALL)
103  {
104  // All -x octants become invisible
105  visOctant &= ~mXmYmZMask;
106  visOctant &= ~mXmYpZMask;
107  visOctant &= ~mXpYmZMask;
108  visOctant &= ~mXpYpZMask;
109  }
110  else if (n.x() < -SMALL)
111  {
112  // All +x octants become invisible
113  visOctant &= ~pXmYmZMask;
114  visOctant &= ~pXmYpZMask;
115  visOctant &= ~pXpYmZMask;
116  visOctant &= ~pXpYpZMask;
117  }
118 
119  if (n.y() > SMALL)
120  {
121  visOctant &= ~mXmYmZMask;
122  visOctant &= ~mXmYpZMask;
123  visOctant &= ~pXmYmZMask;
124  visOctant &= ~pXmYpZMask;
125  }
126  else if (n.y() < -SMALL)
127  {
128  visOctant &= ~mXpYmZMask;
129  visOctant &= ~mXpYpZMask;
130  visOctant &= ~pXpYmZMask;
131  visOctant &= ~pXpYpZMask;
132  }
133 
134  if (n.z() > SMALL)
135  {
136  visOctant &= ~mXmYmZMask;
137  visOctant &= ~mXpYmZMask;
138  visOctant &= ~pXmYmZMask;
139  visOctant &= ~pXpYmZMask;
140  }
141  else if (n.z() < -SMALL)
142  {
143  visOctant &= ~mXmYpZMask;
144  visOctant &= ~mXpYpZMask;
145  visOctant &= ~pXmYpZMask;
146  visOctant &= ~pXpYpZMask;
147  }
148  }
149 
150  label visI = -1;
151 
152  label mask = 1;
153 
154  for (label octant = 0; octant < 8; octant++)
155  {
156  if (visOctant & mask)
157  {
158  // Take first visible octant found. Note: could use
159  // combination of visible octants here.
160  visI = octant;
161 
162  break;
163  }
164  mask <<= 1;
165  }
166 
167  if (visI != -1)
168  {
169  // Take a vector in this octant.
170  pn[pointI] = octantNormal[visI];
171  }
172  else
173  {
174  pn[pointI] = vector::zero;
175 
177  << "No visible octant for point:" << pp.meshPoints()[pointI]
178  << " cooord:" << pp.points()[pp.meshPoints()[pointI]] << nl
179  << "Normal set to " << pn[pointI] << endl;
180  }
181  }
182  }
183 
184  return pn;
185 }
186 
187 
189 (
190  const primitiveMesh& mesh,
191  const label edgeI
192 )
193 {
194  vector eVec = mesh.edges()[edgeI].vec(mesh.points());
195 
196  eVec /= mag(eVec);
197 
198  return eVec;
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 triad& t,
216  const point& pt
217 )
218 {
219  forAll(t, dirI)
220  {
221  writeOBJ(os, pt, pt + t[dirI]);
222  }
223 }
224 
225 
227 (
228  Ostream& os,
229  const point& p1,
230  const point& p2,
231  label& count
232 )
233 {
234  os << "v" << ' ' << p1.x() << ' ' << p1.y() << ' ' << p1.z() << endl;
235  os << "v" << ' ' << p2.x() << ' ' << p2.y() << ' ' << p2.z() << endl;
236 
237  os << "l" << " " << (count + 1) << " " << (count + 2) << endl;
238 
239  count += 2;
240 }
241 
242 
244 (
245  Ostream& os,
246  const point& p1,
247  const point& p2
248 )
249 {
250  os << "v" << ' ' << p1.x() << ' ' << p1.y() << ' ' << p1.z() << endl;
251 
252  os << "vn"
253  << ' ' << p2.x() - p1.x()
254  << ' ' << p2.y() - p1.y()
255  << ' ' << p2.z() - p1.z() << endl;
256 }
257 
258 
260 (
261  Ostream& os,
262  const cellList& cells,
263  const faceList& faces,
264  const pointField& points,
265  const labelList& cellLabels
266 )
267 {
268  labelHashSet usedFaces(4*cellLabels.size());
269 
270  forAll(cellLabels, i)
271  {
272  const cell& cFaces = cells[cellLabels[i]];
273 
274  forAll(cFaces, j)
275  {
276  usedFaces.insert(cFaces[j]);
277  }
278  }
279 
280  writeOBJ(os, faces, points, usedFaces.toc());
281 }
282 
283 
285 (
286  const primitiveMesh& mesh,
287  const label cellI,
288  const label edgeI
289 )
290 {
291  return findIndex(mesh.edgeCells(edgeI), cellI) != -1;
292 }
293 
294 
296 (
297  const primitiveMesh& mesh,
298  const label faceI,
299  const label edgeI
300 )
301 {
302  return findIndex(mesh.faceEdges(faceI), edgeI) != -1;
303 }
304 
305 
307 (
308  const primitiveMesh& mesh,
309  const label cellI,
310  const label faceI
311 )
312 {
313  if (mesh.isInternalFace(faceI))
314  {
315  if
316  (
317  (mesh.faceOwner()[faceI] == cellI)
318  || (mesh.faceNeighbour()[faceI] == cellI)
319  )
320  {
321  return true;
322  }
323  }
324  else
325  {
326  if (mesh.faceOwner()[faceI] == cellI)
327  {
328  return true;
329  }
330  }
331  return false;
332 }
333 
334 
336 (
337  const edgeList& edges,
338  const labelList& candidates,
339  const label v0,
340  const label v1
341 )
342 {
343  forAll(candidates, i)
344  {
345  label edgeI = candidates[i];
346 
347  const edge& e = edges[edgeI];
348 
349  if ((e[0] == v0 && e[1] == v1) || (e[0] == v1 && e[1] == v0))
350  {
351  return edgeI;
352  }
353  }
354  return -1;
355 }
356 
357 
359 (
360  const primitiveMesh& mesh,
361  const label v0,
362  const label v1
363 )
364 {
365  const edgeList& edges = mesh.edges();
366 
367  const labelList& v0Edges = mesh.pointEdges()[v0];
368 
369  forAll(v0Edges, i)
370  {
371  label edgeI = v0Edges[i];
372 
373  const edge& e = edges[edgeI];
374 
375  if ((e.start() == v1) || (e.end() == v1))
376  {
377  return edgeI;
378  }
379  }
380  return -1;
381 }
382 
383 
385 (
386  const primitiveMesh& mesh,
387  const label f0,
388  const label f1
389 )
390 {
391  const labelList& f0Edges = mesh.faceEdges()[f0];
392  const labelList& f1Edges = mesh.faceEdges()[f1];
393 
394  forAll(f0Edges, f0EdgeI)
395  {
396  label edge0 = f0Edges[f0EdgeI];
397 
398  forAll(f1Edges, f1EdgeI)
399  {
400  label edge1 = f1Edges[f1EdgeI];
401 
402  if (edge0 == edge1)
403  {
404  return edge0;
405  }
406  }
407  }
409  << "Faces " << f0 << " and " << f1 << " do not share an edge"
410  << abort(FatalError);
411 
412  return -1;
413 
414 }
415 
416 
418 (
419  const primitiveMesh& mesh,
420  const label cell0I,
421  const label cell1I
422 )
423 {
424  const cell& cFaces = mesh.cells()[cell0I];
425 
426  forAll(cFaces, cFaceI)
427  {
428  label faceI = cFaces[cFaceI];
429 
430  if
431  (
432  mesh.isInternalFace(faceI)
433  && (
434  mesh.faceOwner()[faceI] == cell1I
435  || mesh.faceNeighbour()[faceI] == cell1I
436  )
437  )
438  {
439  return faceI;
440  }
441  }
442 
443 
445  << "No common face for"
446  << " cell0I:" << cell0I << " faces:" << cFaces
447  << " cell1I:" << cell1I << " faces:"
448  << mesh.cells()[cell1I]
449  << abort(FatalError);
450 
451  return -1;
452 }
453 
454 
456 (
457  const primitiveMesh& mesh,
458  const label cellI,
459  const label edgeI,
460  label& face0,
461  label& face1
462 )
463 {
464  const labelList& eFaces = mesh.edgeFaces(edgeI);
465 
466  face0 = -1;
467  face1 = -1;
468 
469  forAll(eFaces, eFaceI)
470  {
471  label faceI = eFaces[eFaceI];
472 
473  if (faceOnCell(mesh, cellI, faceI))
474  {
475  if (face0 == -1)
476  {
477  face0 = faceI;
478  }
479  else
480  {
481  face1 = faceI;
482 
483  return;
484  }
485  }
486  }
487 
488  if ((face0 == -1) || (face1 == -1))
489  {
491  << "Can not find faces using edge " << mesh.edges()[edgeI]
492  << " on cell " << cellI << abort(FatalError);
493  }
494 }
495 
496 
498 (
499  const primitiveMesh& mesh,
500  const labelList& edgeLabels,
501  const label thisEdgeI,
502  const label thisVertI
503 )
504 {
505  forAll(edgeLabels, edgeLabelI)
506  {
507  label edgeI = edgeLabels[edgeLabelI];
508 
509  if (edgeI != thisEdgeI)
510  {
511  const edge& e = mesh.edges()[edgeI];
512 
513  if ((e.start() == thisVertI) || (e.end() == thisVertI))
514  {
515  return edgeI;
516  }
517  }
518  }
519 
521  << "Can not find edge in "
522  << UIndirectList<edge>(mesh.edges(), edgeLabels)()
523  << " connected to edge "
524  << thisEdgeI << " with vertices " << mesh.edges()[thisEdgeI]
525  << " on side " << thisVertI << abort(FatalError);
526 
527  return -1;
528 }
529 
530 
532 (
533  const primitiveMesh& mesh,
534  const label cellI,
535  const label faceI,
536  const label edgeI
537 )
538 {
539  label face0;
540  label face1;
541 
542  getEdgeFaces(mesh, cellI, edgeI, face0, face1);
543 
544  if (face0 == faceI)
545  {
546  return face1;
547  }
548  else
549  {
550  return face0;
551  }
552 }
553 
554 
556 (
557  const primitiveMesh& mesh,
558  const label otherCellI,
559  const label faceI
560 )
561 {
562  if (!mesh.isInternalFace(faceI))
563  {
565  << "Face " << faceI << " is not internal"
566  << abort(FatalError);
567  }
568 
569  label newCellI = mesh.faceOwner()[faceI];
570 
571  if (newCellI == otherCellI)
572  {
573  newCellI = mesh.faceNeighbour()[faceI];
574  }
575  return newCellI;
576 }
577 
578 
580 (
581  const primitiveMesh& mesh,
582  const label faceI,
583  const label startEdgeI,
584  const label startVertI,
585  const label nEdges
586 )
587 {
588  const labelList& fEdges = mesh.faceEdges(faceI);
589 
590  label edgeI = startEdgeI;
591 
592  label vertI = startVertI;
593 
594  for (label iter = 0; iter < nEdges; iter++)
595  {
596  edgeI = otherEdge(mesh, fEdges, edgeI, vertI);
597 
598  vertI = mesh.edges()[edgeI].otherVertex(vertI);
599  }
600 
601  return edgeI;
602 }
603 
604 
606 (
607  const polyMesh& mesh,
608  point& pt
609 )
610 {
611  const Vector<label>& dirs = mesh.geometricD();
612 
613  const point& min = mesh.bounds().min();
614  const point& max = mesh.bounds().max();
615 
616  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
617  {
618  if (dirs[cmpt] == -1)
619  {
620  pt[cmpt] = 0.5*(min[cmpt] + max[cmpt]);
621  }
622  }
623 }
624 
625 
627 (
628  const polyMesh& mesh,
629  pointField& pts
630 )
631 {
632  const Vector<label>& dirs = mesh.geometricD();
633 
634  const point& min = mesh.bounds().min();
635  const point& max = mesh.bounds().max();
636 
637  bool isConstrained = false;
638  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
639  {
640  if (dirs[cmpt] == -1)
641  {
642  isConstrained = true;
643  break;
644  }
645  }
646 
647  if (isConstrained)
648  {
649  forAll(pts, i)
650  {
651  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
652  {
653  if (dirs[cmpt] == -1)
654  {
655  pts[i][cmpt] = 0.5*(min[cmpt] + max[cmpt]);
656  }
657  }
658  }
659  }
660 }
661 
662 
664 (
665  const polyMesh& mesh,
666  const Vector<label>& dirs,
667  vector& d
668 )
669 {
670  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
671  {
672  if (dirs[cmpt] == -1)
673  {
674  d[cmpt] = 0.0;
675  }
676  }
677 }
678 
679 
681 (
682  const polyMesh& mesh,
683  const Vector<label>& dirs,
684  vectorField& d
685 )
686 {
687  bool isConstrained = false;
688  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
689  {
690  if (dirs[cmpt] == -1)
691  {
692  isConstrained = true;
693  break;
694  }
695  }
696 
697  if (isConstrained)
698  {
699  forAll(d, i)
700  {
701  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
702  {
703  if (dirs[cmpt] == -1)
704  {
705  d[i][cmpt] = 0.0;
706  }
707  }
708  }
709  }
710 }
711 
712 
714 (
715  const primitiveMesh& mesh,
716  const label cellI,
717  const label e0,
718  label& e1,
719  label& e2,
720  label& e3
721 )
722 {
723  // Go to any face using e0
724  label faceI = meshTools::otherFace(mesh, cellI, -1, e0);
725 
726  // Opposite edge on face
727  e1 = meshTools::walkFace(mesh, faceI, e0, mesh.edges()[e0].end(), 2);
728 
729  faceI = meshTools::otherFace(mesh, cellI, faceI, e1);
730 
731  e2 = meshTools::walkFace(mesh, faceI, e1, mesh.edges()[e1].end(), 2);
732 
733  faceI = meshTools::otherFace(mesh, cellI, faceI, e2);
734 
735  e3 = meshTools::walkFace(mesh, faceI, e2, mesh.edges()[e2].end(), 2);
736 }
737 
738 
740 (
741  const primitiveMesh& mesh,
742  const label cellI,
743  const label startEdgeI
744 )
745 {
746  if (!hexMatcher().isA(mesh, cellI))
747  {
749  << "Not a hex : cell:" << cellI << abort(FatalError);
750  }
751 
752 
753  vector avgVec(normEdgeVec(mesh, startEdgeI));
754 
755  label edgeI = startEdgeI;
756 
757  label faceI = -1;
758 
759  for (label i = 0; i < 3; i++)
760  {
761  // Step to next face, next edge
762  faceI = meshTools::otherFace(mesh, cellI, faceI, edgeI);
763 
764  vector eVec(normEdgeVec(mesh, edgeI));
765 
766  if ((eVec & avgVec) > 0)
767  {
768  avgVec += eVec;
769  }
770  else
771  {
772  avgVec -= eVec;
773  }
774 
775  label vertI = mesh.edges()[edgeI].end();
776 
777  edgeI = meshTools::walkFace(mesh, faceI, edgeI, vertI, 2);
778  }
779 
780  avgVec /= mag(avgVec) + VSMALL;
781 
782  return avgVec;
783 }
784 
785 
787 (
788  const primitiveMesh& mesh,
789  const label cellI,
790  const vector& cutDir
791 )
792 {
793  if (!hexMatcher().isA(mesh, cellI))
794  {
796  << "Not a hex : cell:" << cellI << abort(FatalError);
797  }
798 
799  const labelList& cEdges = mesh.cellEdges()[cellI];
800 
801  labelHashSet doneEdges(2*cEdges.size());
802 
803  scalar maxCos = -GREAT;
804  label maxEdgeI = -1;
805 
806  for (label i = 0; i < 4; i++)
807  {
808  forAll(cEdges, cEdgeI)
809  {
810  label e0 = cEdges[cEdgeI];
811 
812  if (!doneEdges.found(e0))
813  {
814  vector avgDir(edgeToCutDir(mesh, cellI, e0));
815 
816  scalar cosAngle = mag(avgDir & cutDir);
817 
818  if (cosAngle > maxCos)
819  {
820  maxCos = cosAngle;
821  maxEdgeI = e0;
822  }
823 
824  // Mark off edges in cEdges.
825  label e1, e2, e3;
826  getParallelEdges(mesh, cellI, e0, e1, e2, e3);
827 
828  doneEdges.insert(e0);
829  doneEdges.insert(e1);
830  doneEdges.insert(e2);
831  doneEdges.insert(e3);
832  }
833  }
834  }
835 
836  forAll(cEdges, cEdgeI)
837  {
838  if (!doneEdges.found(cEdges[cEdgeI]))
839  {
841  << "Cell:" << cellI << " edges:" << cEdges << endl
842  << "Edge:" << cEdges[cEdgeI] << " not yet handled"
843  << abort(FatalError);
844  }
845  }
846 
847  if (maxEdgeI == -1)
848  {
850  << "Problem : did not find edge aligned with " << cutDir
851  << " on cell " << cellI << abort(FatalError);
852  }
853 
854  return maxEdgeI;
855 }
856 
857 
858 // ************************************************************************* //
Foam::meshTools::mXpYmZMask
static const label mXpYmZMask
Definition: meshTools.H:75
Foam::meshTools::constrainDirection
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:664
meshTools.H
Foam::PrimitivePatch::pointFaces
const labelListList & pointFaces() const
Return point-face addressing.
Definition: PrimitivePatchTemplate.C:352
Foam::meshTools::mXmYpZ
static const label mXmYpZ
Definition: meshTools.H:68
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::meshTools::mXmYpZMask
static const label mXmYpZMask
Definition: meshTools.H:78
Foam::HashTable::toc
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
writeOBJ
void writeOBJ(Ostream &os, label &vertI, const tetPoints &tet)
Definition: Test-tetTetOverlap.C:38
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Foam::meshTools::edgeOnCell
bool edgeOnCell(const primitiveMesh &, const label cellI, const label edgeI)
Is edge used by cell.
Definition: meshTools.C:285
Foam::meshTools::pXpYmZ
static const label pXpYmZ
Definition: meshTools.H:66
hexMatcher.H
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::meshTools::pXpYpZMask
static const label pXpYpZMask
Definition: meshTools.H:81
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::meshTools::mXpYpZ
static const label mXpYpZ
Definition: meshTools.H:70
Foam::meshTools::pXmYmZ
static const label pXmYmZ
Definition: meshTools.H:64
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::meshTools::pXmYpZMask
static const label pXmYpZMask
Definition: meshTools.H:79
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::meshTools::faceOnCell
bool faceOnCell(const primitiveMesh &, const label cellI, const label faceI)
Is face used by cell.
Definition: meshTools.C:307
Foam::meshTools::getParallelEdges
void getParallelEdges(const primitiveMesh &, const label cellI, const label e0, label &, label &, label &)
Given edge on hex find other 'parallel', non-connected edges.
Definition: meshTools.C:714
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("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]=cellShape(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
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::meshTools::pXmYmZMask
static const label pXmYmZMask
Definition: meshTools.H:74
Foam::hexMatcher
A cellMatcher for hex cells.
Definition: hexMatcher.H:51
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::meshTools::findEdge
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:336
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::meshTools::pXmYpZ
static const label pXmYpZ
Definition: meshTools.H:69
f1
scalar f1
Definition: createFields.H:28
Foam::meshTools::edgeToCutDir
vector edgeToCutDir(const primitiveMesh &, const label cellI, const label edgeI)
Given edge on hex find all 'parallel' (i.e. non-connected)
Definition: meshTools.C:740
Foam::meshTools::calcBoxPointNormals
vectorField calcBoxPointNormals(const primitivePatch &pp)
Calculate point normals on a 'box' mesh (all edges aligned with.
Definition: meshTools.C:52
Foam::meshTools::walkFace
label walkFace(const primitiveMesh &, 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:580
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatchTemplate.H:293
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
Foam::meshTools::getSharedEdge
label getSharedEdge(const primitiveMesh &, const label f0, const label f1)
Return edge shared by two faces. Throws error if no edge found.
Definition: meshTools.C:385
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::meshTools::otherEdge
label otherEdge(const primitiveMesh &, const labelList &edgeLabels, const label edgeI, const label vertI)
Return label of other edge (out of candidates edgeLabels)
Definition: meshTools.C:498
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::meshTools::mXmYmZMask
static const label mXmYmZMask
Definition: meshTools.H:73
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Foam::meshTools::edgeOnFace
bool edgeOnFace(const primitiveMesh &, const label faceI, const label edgeI)
Is edge used by face.
Definition: meshTools.C:296
Foam::meshTools::mXpYmZ
static const label mXpYmZ
Definition: meshTools.H:65
Foam::meshTools::getEdgeFaces
void getEdgeFaces(const primitiveMesh &, 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:456
Foam::triad
Representation of a 3D Cartesian coordinate system as a Vector of vectors.
Definition: triad.H:65
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::PrimitivePatch::pointNormals
const Field< PointType > & pointNormals() const
Return point normals for patch.
Definition: PrimitivePatchTemplate.C:540
Foam::meshTools::cutDirToEdge
label cutDirToEdge(const primitiveMesh &, const label cellI, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:787
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::meshTools::normEdgeVec
vector normEdgeVec(const primitiveMesh &, const label edgeI)
Normalized edge vector.
Definition: meshTools.C:189
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face normals for patch.
Definition: PrimitivePatchTemplate.C:520
Foam::meshTools::mXmYmZ
static const label mXmYmZ
Definition: meshTools.H:63
Foam::meshTools::otherFace
label otherFace(const primitiveMesh &, const label cellI, const label faceI, const label edgeI)
Return face on cell using edgeI but not faceI. Throws error.
Definition: meshTools.C:532
Foam::Vector< scalar >
Foam::meshTools::getSharedFace
label getSharedFace(const primitiveMesh &, const label cell0, const label cell1)
Return face shared by two cells. Throws error if none found.
Definition: meshTools.C:418
Foam::isA
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::meshTools::constrainToMeshCentre
void constrainToMeshCentre(const polyMesh &mesh, point &pt)
Set the constrained components of position to mesh centre.
Definition: meshTools.C:606
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::meshTools::pXpYmZMask
static const label pXpYmZMask
Definition: meshTools.H:76
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::meshTools::visNormal
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:33
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::meshTools::mXpYpZMask
static const label mXpYpZMask
Definition: meshTools.H:80
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::meshTools::otherCell
label otherCell(const primitiveMesh &, const label cellI, const label faceI)
Return cell on other side of face. Throws error.
Definition: meshTools.C:556
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::meshTools::pXpYpZ
static const label pXpYpZ
Definition: meshTools.H:71
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79