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