tetDecomposer.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "tetDecomposer.H"
30 #include "meshTools.H"
31 #include "polyMesh.H"
32 #include "polyTopoChange.H"
33 #include "mapPolyMesh.H"
34 #include "OFstream.H"
35 #include "edgeHashes.H"
36 #include "syncTools.H"
37 #include "triPointRef.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(tetDecomposer, 0);
44 }
45 
46 const Foam::Enum
47 <
49 >
51 ({
52  { decompositionType::FACE_CENTRE_TRIS, "faceCentre" },
53  { decompositionType::FACE_DIAG_TRIS, "faceDiagonal" },
54  { decompositionType::PYRAMID, "pyramid" },
55 });
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 void Foam::tetDecomposer::modifyFace
61 (
62  polyTopoChange& meshMod,
63  const face& f,
64  const label facei,
65  const label own,
66  const label nei,
67  const label patchi,
68  const label zoneI,
69  const bool zoneFlip
70 ) const
71 {
72  // First usage of face. Modify.
73  if (nei == -1 || own < nei)
74  {
75  meshMod.modifyFace
76  (
77  f, // modified face
78  facei, // label of face
79  own, // owner
80  nei, // neighbour
81  false, // face flip
82  patchi, // patch for face
83  zoneI, // zone for face
84  zoneFlip // face flip in zone
85  );
86  }
87  else
88  {
89  meshMod.modifyFace
90  (
91  f.reverseFace(), // modified face
92  facei, // label of face
93  nei, // owner
94  own, // neighbour
95  true, // face flip
96  patchi, // patch for face
97  zoneI, // zone for face
98  !zoneFlip // face flip in zone
99  );
100  }
101 }
102 
103 
104 void Foam::tetDecomposer::addFace
105 (
106  polyTopoChange& meshMod,
107  const face& f,
108  const label own,
109  const label nei,
110  const label masterPointID,
111  const label masterEdgeID,
112  const label masterFaceID,
113  const label patchi,
114  const label zoneI,
115  const bool zoneFlip
116 ) const
117 {
118  // Second or more usage of face. Add.
119  if (nei == -1 || own < nei)
120  {
121  meshMod.addFace
122  (
123  f, // modified face
124  own, // owner
125  nei, // neighbour
126  masterPointID, // master point
127  masterEdgeID, // master edge
128  masterFaceID, // master face
129  false, // face flip
130  patchi, // patch for face
131  zoneI, // zone for face
132  zoneFlip // face flip in zone
133  );
134  }
135  else
136  {
137  meshMod.addFace
138  (
139  f.reverseFace(), // modified face
140  nei, // owner
141  own, // neighbour
142  masterPointID, // master point
143  masterEdgeID, // master edge
144  masterFaceID, // master face
145  true, // face flip
146  patchi, // patch for face
147  zoneI, // zone for face
148  !zoneFlip // face flip in zone
149  );
150  }
151 }
152 
153 
154 // Work out triangle index given the starting vertex in the face
155 Foam::label Foam::tetDecomposer::triIndex(const label facei, const label fp)
156 const
157 {
158  const face& f = mesh_.faces()[facei];
159  const label fp0 = max(0, mesh_.tetBasePtIs()[facei]);
160 
161  // Work out triangle index on this face
162  label thisTrii;
163  if (fp == fp0)
164  {
165  thisTrii = 0;
166  }
167  else if (fp == f.rcIndex(fp0))
168  {
169  thisTrii = f.size()-3;
170  }
171  else
172  {
173  thisTrii = (fp-fp0-1) % (f.size()-2);
174  }
175  return thisTrii;
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
180 
181 Foam::tetDecomposer::tetDecomposer(const polyMesh& mesh)
182 :
183  mesh_(mesh)
184 {}
185 
186 
187 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
188 
190 (
191  const decompositionType decomposeType,
192  const bitSet& decomposeCell,
193  polyTopoChange& meshMod
194 )
195 {
196  cellToPoint_.setSize(mesh_.nCells(), -1);
197  forAll(mesh_.cellCentres(), celli)
198  {
199  if (decomposeCell[celli])
200  {
201  // Any point on the cell
202  label masterPointi = mesh_.faces()[mesh_.cells()[celli][0]][0];
203 
204  cellToPoint_[celli] = meshMod.addPoint
205  (
206  mesh_.cellCentres()[celli],
207  masterPointi,
208  -1,
209  true
210  );
211  }
212  }
213 
214 
215  // Determine for every face whether it borders a cell that is decomposed
216  bitSet decomposeFace(mesh_.nFaces());
217  {
218  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
219  {
220  label own = mesh_.faceOwner()[facei];
221  label nei = mesh_.faceNeighbour()[facei];
222  if (decomposeCell[own] || decomposeCell[nei])
223  {
224  decomposeFace[facei] = true;
225  }
226  }
227 
228  boolList neiDecomposeCell(mesh_.nBoundaryFaces());
229  forAll(neiDecomposeCell, bFacei)
230  {
231  label facei = mesh_.nInternalFaces()+bFacei;
232  label own = mesh_.faceOwner()[facei];
233  neiDecomposeCell[bFacei] = decomposeCell[own];
234  }
235  syncTools::swapBoundaryFaceList(mesh_, neiDecomposeCell);
236 
237  for
238  (
239  label facei = mesh_.nInternalFaces();
240  facei < mesh_.nFaces();
241  facei++
242  )
243  {
244  label own = mesh_.faceOwner()[facei];
245  label bFacei = facei-mesh_.nInternalFaces();
246  if (decomposeCell[own] || neiDecomposeCell[bFacei])
247  {
248  decomposeFace[facei] = true;
249  }
250  }
251  }
252 
253 
254  // Add face centre points
255  if (decomposeType == FACE_CENTRE_TRIS)
256  {
257  faceToPoint_.setSize(mesh_.nFaces(), -1);
258  forAll(mesh_.faceCentres(), facei)
259  {
260  if (decomposeFace[facei])
261  {
262  // Any point on the face
263  const label masterPointi = mesh_.faces()[facei][0];
264 
265  faceToPoint_[facei] = meshMod.addPoint
266  (
267  mesh_.faceCentres()[facei],
268  masterPointi,
269  -1,
270  true
271  );
272  }
273  }
274  }
275 
276 
277  // Per face, per point (faceCentre) or triangle (faceDiag) the (existing
278  // or added) cell on either side
279  faceOwnerCells_.setSize(mesh_.nFaces());
280  faceNeighbourCells_.setSize(mesh_.nFaces());
281 
282  if (decomposeType == FACE_CENTRE_TRIS)
283  {
284  forAll(faceOwnerCells_, facei)
285  {
286  const face& f = mesh_.faces()[facei];
287  faceOwnerCells_[facei].setSize(f.size(), -1);
288  faceNeighbourCells_[facei].setSize(f.size(), -1);
289  }
290  }
291  else if (decomposeType == FACE_DIAG_TRIS)
292  {
293  // Force construction of diagonal decomposition
294  (void)mesh_.tetBasePtIs();
295 
296  forAll(faceOwnerCells_, facei)
297  {
298  const face& f = mesh_.faces()[facei];
299  faceOwnerCells_[facei].setSize(f.size()-2, -1);
300  faceNeighbourCells_[facei].setSize(f.size()-2, -1);
301  }
302  }
303  else
304  {
305  forAll(faceOwnerCells_, facei)
306  {
307  faceOwnerCells_[facei].setSize(1, -1);
308  faceNeighbourCells_[facei].setSize(1, -1);
309  }
310  }
311 
312 
313  // Add internal cells. Note: done in same order as pyramid triangle
314  // creation later to maintain same ordering.
315  forAll(mesh_.cells(), celli)
316  {
317  const cell& cFaces = mesh_.cells()[celli];
318 
319  // Whether cell has already been modified (all other cells get added)
320  bool modifiedCell = false;
321 
322  forAll(cFaces, cFacei)
323  {
324  label facei = cFaces[cFacei];
325  const face& f = mesh_.faces()[facei];
326 
327  // Get reference to either owner or neighbour
328  labelList& added =
329  (
330  (mesh_.faceOwner()[facei] == celli)
331  ? faceOwnerCells_[facei]
332  : faceNeighbourCells_[facei]
333  );
334 
335  if (decomposeCell[celli])
336  {
337  if (decomposeType == FACE_CENTRE_TRIS)
338  {
339  forAll(f, fp)
340  {
341  if (!modifiedCell)
342  {
343  // Reuse cell itself
344  added[fp] = celli;
345  modifiedCell = true;
346  }
347  else
348  {
349  added[fp] = meshMod.addCell
350  (
351  -1, // masterPoint
352  -1, // masterEdge
353  -1, // masterFace
354  celli, // masterCell
355  mesh_.cellZones().whichZone(celli)
356  );
357  }
358  }
359  }
360  else if (decomposeType == FACE_DIAG_TRIS)
361  {
362  for (label triI = 0; triI < f.size()-2; triI++)
363  {
364  if (!modifiedCell)
365  {
366  // Reuse cell itself
367  added[triI] = celli;
368  modifiedCell = true;
369  }
370  else
371  {
372  added[triI] = meshMod.addCell
373  (
374  -1, // masterPoint
375  -1, // masterEdge
376  -1, // masterFace
377  celli, // masterCell
378  mesh_.cellZones().whichZone(celli)
379  );
380  //Pout<< "For cell:" << celli
381  // << " at:" << mesh_.cellCentres()[celli]
382  // << " face:" << facei
383  // << " at:" << mesh_.faceCentres()[facei]
384  // << " tri:" << triI
385  // << " added cell:" << added[triI] << endl;
386  }
387  }
388  }
389  else // if (decomposeType == PYRAMID)
390  {
391  // Pyramidal decomposition.
392  // Assign same cell to all slots
393  if (!modifiedCell)
394  {
395  // Reuse cell itself
396  added = celli;
397  modifiedCell = true;
398  }
399  else
400  {
401  added = meshMod.addCell
402  (
403  -1, // masterPoint
404  -1, // masterEdge
405  -1, // masterFace
406  celli, // masterCell
407  mesh_.cellZones().whichZone(celli)
408  );
409  }
410  }
411  }
412  else
413  {
414  // All vertices/triangles address to original cell
415  added = celli;
416  }
417  }
418  }
419 
420 
421 
422  // Add triangle faces
423  face triangle(3);
424 
425  forAll(mesh_.faces(), facei)
426  {
427  label own = mesh_.faceOwner()[facei];
428  const labelList& addedOwn = faceOwnerCells_[facei];
429  const labelList& addedNei = faceNeighbourCells_[facei];
430  const face& f = mesh_.faces()[facei];
431 
432  label patchi = -1;
433  if (facei >= mesh_.nInternalFaces())
434  {
435  patchi = mesh_.boundaryMesh().whichPatch(facei);
436  }
437 
438  label zoneI = mesh_.faceZones().whichZone(facei);
439  bool zoneFlip = false;
440  if (zoneI != -1)
441  {
442  const faceZone& fz = mesh_.faceZones()[zoneI];
443  zoneFlip = fz.flipMap()[fz.whichFace(facei)];
444  }
445 
446  //Pout<< "Face:" << facei << " at:" << mesh_.faceCentres()[facei]
447  // << endl;
448 
449  if (decomposeType == FACE_CENTRE_TRIS && decomposeFace[facei])
450  {
451  forAll(f, fp)
452  {
453  // 1. Front triangle (decomposition of face itself)
454  // (between owner and neighbour cell)
455  {
456  triangle[0] = f[fp];
457  triangle[1] = f[f.fcIndex(fp)];
458  triangle[2] = faceToPoint_[facei];
459 
460  //Pout<< " triangle:" << triangle
461  // << " points:"
462  // << UIndirectList<point>(meshMod.points(), triangle)
463  // << " between:" << addedOwn[fp]
464  // << " and:" << addedNei[fp] << endl;
465 
466 
467  if (fp == 0)
468  {
469  modifyFace
470  (
471  meshMod,
472  triangle,
473  facei,
474  addedOwn[fp],
475  addedNei[fp],
476  patchi,
477  zoneI,
478  zoneFlip
479  );
480  }
481  else
482  {
483  addFace
484  (
485  meshMod,
486  triangle,
487  addedOwn[fp],
488  addedNei[fp],
489  -1, //point
490  -1, //edge
491  facei, //face
492  patchi,
493  zoneI,
494  zoneFlip
495  );
496  }
497  }
498 
499 
500  // 2. Within owner cell - to cell centre
501  if (decomposeCell[own])
502  {
503  label newOwn = addedOwn[f.rcIndex(fp)];
504  label newNei = addedOwn[fp];
505 
506  triangle[0] = f[fp];
507  triangle[1] = cellToPoint_[own];
508  triangle[2] = faceToPoint_[facei];
509 
510  addFace
511  (
512  meshMod,
513  triangle,
514  newOwn,
515  newNei,
516  f[fp], //point
517  -1, //edge
518  -1, //face
519  -1, //patchi
520  -1, //zone
521  false
522  );
523  }
524  // 2b. Within neighbour cell - to cell centre
525  if
526  (
527  facei < mesh_.nInternalFaces()
528  && decomposeCell[mesh_.faceNeighbour()[facei]]
529  )
530  {
531  label newOwn = addedNei[f.rcIndex(fp)];
532  label newNei = addedNei[fp];
533 
534  triangle[0] = f[fp];
535  triangle[1] = faceToPoint_[facei];
536  triangle[2] = cellToPoint_[mesh_.faceNeighbour()[facei]];
537 
538  addFace
539  (
540  meshMod,
541  triangle,
542  newOwn,
543  newNei,
544  f[fp], //point
545  -1, //edge
546  -1, //face
547  -1, //patchi
548  -1, //zone
549  false
550  );
551  }
552  }
553  }
554  else if (decomposeType == FACE_DIAG_TRIS && decomposeFace[facei])
555  {
556  label fp0 = max(mesh_.tetBasePtIs()[facei], 0);
557  label fp = f.fcIndex(fp0);
558 
559  for (label triI = 0; triI < f.size()-2; triI++)
560  {
561  label nextTri = triI+1;
562  if (nextTri >= f.size()-2)
563  {
564  nextTri -= f.size()-2;
565  }
566  label nextFp = f.fcIndex(fp);
567 
568 
569  // Triangle triI consisiting of f[fp0], f[fp], f[nextFp]
570 
571 
572  // 1. Front triangle (decomposition of face itself)
573  // (between owner and neighbour cell)
574  {
575  triangle[0] = f[fp0];
576  triangle[1] = f[fp];
577  triangle[2] = f[nextFp];
578 
579  if (triI == 0)
580  {
581  modifyFace
582  (
583  meshMod,
584  triangle,
585  facei,
586  addedOwn[triI],
587  addedNei[triI],
588  patchi,
589  zoneI,
590  zoneFlip
591  );
592  }
593  else
594  {
595  addFace
596  (
597  meshMod,
598  triangle,
599  addedOwn[triI],
600  addedNei[triI],
601  -1, //point
602  -1, //edge
603  facei, //face
604  patchi,
605  zoneI,
606  zoneFlip
607  );
608  }
609  }
610 
611 
612  // 2. Within owner cell - diagonal to cell centre
613  if (triI < f.size()-3)
614  {
615  if (decomposeCell[own])
616  {
617  label newOwn = addedOwn[triI];
618  label newNei = addedOwn[nextTri];
619 
620  triangle[0] = f[fp0];
621  triangle[1] = f[nextFp];
622  triangle[2] = cellToPoint_[own];
623 
624  addFace
625  (
626  meshMod,
627  triangle,
628  newOwn,
629  newNei,
630  f[fp], //point
631  -1, //edge
632  -1, //face
633  -1, //patchi
634  -1, //zone
635  false
636  );
637  }
638 
639  // 2b. Within neighbour cell - to cell centre
640  if
641  (
642  facei < mesh_.nInternalFaces()
643  && decomposeCell[mesh_.faceNeighbour()[facei]]
644  )
645  {
646  label newOwn = addedNei[triI];
647  label newNei = addedNei[nextTri];
648 
649  triangle[0] = f[nextFp];
650  triangle[1] = f[fp0];
651  triangle[2] =
652  cellToPoint_[mesh_.faceNeighbour()[facei]];
653 
654  addFace
655  (
656  meshMod,
657  triangle,
658  newOwn,
659  newNei,
660  f[fp], //point
661  -1, //edge
662  -1, //face
663  -1, //patchi
664  -1, //zone
665  false
666  );
667  }
668  }
669 
670 
671  fp = nextFp;
672  }
673  }
674  else
675  {
676  // No decomposition. Use zero'th slot.
677  modifyFace
678  (
679  meshMod,
680  f,
681  facei,
682  addedOwn[0],
683  addedNei[0],
684  patchi,
685  zoneI,
686  zoneFlip
687  );
688  }
689  }
690 
691 
692 
693  // Add triangles for all edges.
694  EdgeMap<label> edgeToFace;
695 
696  forAll(mesh_.cells(), celli)
697  {
698  const cell& cFaces = mesh_.cells()[celli];
699 
700  edgeToFace.clear();
701 
702  forAll(cFaces, cFacei)
703  {
704  label facei = cFaces[cFacei];
705 
706  if (decomposeCell[celli])
707  {
708  const face& f = mesh_.faces()[facei];
709  //const labelList& fEdges = mesh_.faceEdges()[facei];
710  forAll(f, fp)
711  {
712  label p0 = f[fp];
713  label p1 = f[f.fcIndex(fp)];
714  const edge e(p0, p1);
715 
716  EdgeMap<label>::const_iterator edgeFnd = edgeToFace.find(e);
717  if (edgeFnd == edgeToFace.end())
718  {
719  edgeToFace.insert(e, facei);
720  }
721  else
722  {
723  // Found the other face on the edge.
724  label otherFacei = edgeFnd();
725  const face& otherF = mesh_.faces()[otherFacei];
726 
727  // Found the other face on the edge. Note that since
728  // we are looping in the same order the tets added for
729  // otherFacei will be before those of facei
730 
731  label otherFp = otherF.find(p0);
732  if (otherF.nextLabel(otherFp) == p1)
733  {
734  // ok. otherFp is first vertex of edge.
735  }
736  else if (otherF.prevLabel(otherFp) == p1)
737  {
738  otherFp = otherF.rcIndex(otherFp);
739  }
740  else
741  {
743  << "problem." << abort(FatalError);
744  }
745 
746 
747  // Triangle from edge to cell centre
748  if (mesh_.faceOwner()[facei] == celli)
749  {
750  triangle[0] = p0;
751  triangle[1] = p1;
752  triangle[2] = cellToPoint_[celli];
753  }
754  else
755  {
756  triangle[0] = p1;
757  triangle[1] = p0;
758  triangle[2] = cellToPoint_[celli];
759  }
760 
761  // Determine tets on either side
762  label thisTet, otherTet;
763 
764  if (decomposeType == FACE_CENTRE_TRIS)
765  {
766  if (mesh_.faceOwner()[facei] == celli)
767  {
768  thisTet = faceOwnerCells_[facei][fp];
769  }
770  else
771  {
772  thisTet = faceNeighbourCells_[facei][fp];
773  }
774 
775  if (mesh_.faceOwner()[otherFacei] == celli)
776  {
777  otherTet =
778  faceOwnerCells_[otherFacei][otherFp];
779  }
780  else
781  {
782  otherTet =
783  faceNeighbourCells_[otherFacei][otherFp];
784  }
785  }
786  else if (decomposeType == FACE_DIAG_TRIS)
787  {
788  label thisTriI = triIndex(facei, fp);
789  if (mesh_.faceOwner()[facei] == celli)
790  {
791  thisTet = faceOwnerCells_[facei][thisTriI];
792  }
793  else
794  {
795  thisTet = faceNeighbourCells_[facei][thisTriI];
796  }
797 
798  label otherTriI = triIndex(otherFacei, otherFp);
799  if (mesh_.faceOwner()[otherFacei] == celli)
800  {
801  otherTet =
802  faceOwnerCells_[otherFacei][otherTriI];
803  }
804  else
805  {
806  otherTet =
807  faceNeighbourCells_[otherFacei][otherTriI];
808  }
809  }
810  else
811  {
812  if (mesh_.faceOwner()[facei] == celli)
813  {
814  thisTet = faceOwnerCells_[facei][0];
815  }
816  else
817  {
818  thisTet = faceNeighbourCells_[facei][0];
819  }
820  if (mesh_.faceOwner()[otherFacei] == celli)
821  {
822  otherTet = faceOwnerCells_[otherFacei][0];
823  }
824  else
825  {
826  otherTet =
827  faceNeighbourCells_[otherFacei][0];
828  }
829  }
830 
831  addFace
832  (
833  meshMod,
834  triangle,
835  otherTet,
836  thisTet,
837  -1, //masterPoint
838  -1, //fEdges[fp], //masterEdge
839  facei, //masterFace
840  -1, //patchi
841  -1, //zone
842  false
843  );
844  }
845  }
846  }
847  }
848  }
849 }
850 
851 
853 {
854  inplaceRenumber(map.reversePointMap(), cellToPoint_);
855  inplaceRenumber(map.reversePointMap(), faceToPoint_);
856 
857  forAll(faceOwnerCells_, facei)
858  {
859  inplaceRenumber(map.reverseCellMap(), faceOwnerCells_[facei]);
860  }
861  forAll(faceNeighbourCells_, facei)
862  {
863  inplaceRenumber(map.reverseCellMap(), faceNeighbourCells_[facei]);
864  }
865 }
866 
867 
868 // ************************************************************************* //
Foam::polyTopoChange::addCell
label addCell(const label masterPointID, const label masterEdgeID, const label masterFaceID, const label masterCellID, const label zoneID)
Add cell. Return new cell label.
Definition: polyTopoChange.C:2895
meshTools.H
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::faceZone::flipMap
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:272
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::tetDecomposer::updateMesh
void updateMesh(const mapPolyMesh &)
Force recalculation of locally stored data on topological change.
Definition: tetDecomposer.C:852
Foam::HashTable::end
iterator end() noexcept
iterator to signal the end (for any HashTable)
Definition: HashTableIterI.H:256
Foam::tetDecomposer::decompositionTypeNames
static const Enum< decompositionType > decompositionTypeNames
Definition: tetDecomposer.H:76
mapPolyMesh.H
polyTopoChange.H
Foam::polyTopoChange::addPoint
label addPoint(const point &pt, const label masterPointID, const label zoneID, const bool inCell)
Add point. Return new point label.
Definition: polyTopoChange.C:2610
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:99
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::tetDecomposer::setRefinement
void setRefinement(const decompositionType decomposeType, const bitSet &decomposeCell, polyTopoChange &meshMod)
Insert all changes into meshMod to convert the polyMesh into.
Definition: tetDecomposer.C:190
Foam::HashTable::insert
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
triPointRef.H
polyMesh.H
syncTools.H
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:61
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::triangle
A triangle primitive used to calculate face normals and swept volumes.
Definition: triangle.H:59
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::tetDecomposer::decompositionType
decompositionType
Definition: tetDecomposer.H:69
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:372
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::HashTable::find
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:114
Foam::EdgeMap< label >
Foam::face::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:181
Foam::mapPolyMesh::reverseCellMap
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:532
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:175
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
edgeHashes.H
Foam::HashTable::clear
void clear()
Clear all entries from table.
Definition: HashTable.C:630
f
labelList f(nPoints)
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
Foam::List< bool >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
tetDecomposer.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
p0
const volScalarField & p0
Definition: EEqn.H:36
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54