addPatchCellLayer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
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 "addPatchCellLayer.H"
27 #include "polyMesh.H"
28 #include "polyTopoChange.H"
29 #include "meshTools.H"
30 #include "polyTopoChangeMap.H"
31 #include "syncTools.H"
32 #include "globalIndex.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::label Foam::addPatchCellLayer::nbrFace
45 (
46  const labelListList& edgeFaces,
47  const label edgeI,
48  const label facei
49 )
50 {
51  const labelList& eFaces = edgeFaces[edgeI];
52 
53  if (eFaces.size() == 2)
54  {
55  return (eFaces[0] != facei ? eFaces[0] : eFaces[1]);
56  }
57  else
58  {
59  return -1;
60  }
61 }
62 
63 
64 void Foam::addPatchCellLayer::addVertex
65 (
66  const label pointi,
67  face& f,
68  label& fp
69 )
70 {
71  if (fp == 0)
72  {
73  f[fp++] = pointi;
74  }
75  else
76  {
77  if (f[fp-1] != pointi && f[0] != pointi)
78  {
79  f[fp++] = pointi;
80  }
81  }
82 }
83 
84 
85 bool Foam::addPatchCellLayer::sameEdgeNeighbour
86 (
87  const indirectPrimitivePatch& pp,
88  const labelListList& globalEdgeFaces,
89  const boolList& doneEdge,
90  const label thisGlobalFacei,
91  const label nbrGlobalFacei,
92  const label edgeI
93 ) const
94 {
95  const edge& e = pp.edges()[edgeI];
96 
97  return
98  !doneEdge[edgeI] // not yet handled
99  && (
100  addedPoints_[e[0]].size() // is extruded
101  || addedPoints_[e[1]].size()
102  )
103  && (
104  nbrFace(globalEdgeFaces, edgeI, thisGlobalFacei)
105  == nbrGlobalFacei // is to same neighbour
106  );
107 }
108 
109 
110 Foam::labelPair Foam::addPatchCellLayer::getEdgeString
111 (
112  const indirectPrimitivePatch& pp,
113  const labelListList& globalEdgeFaces,
114  const boolList& doneEdge,
115  const label patchFacei,
116  const label globalFacei
117 ) const
118 {
119  const labelList& fEdges = pp.faceEdges()[patchFacei];
120 
121  label startFe = -1;
122  label endFe = -1;
123 
124  // Get edge that hasn't been done yet but needs extrusion
125  forAll(fEdges, fp)
126  {
127  label edgeI = fEdges[fp];
128  const edge& e = pp.edges()[edgeI];
129 
130  if
131  (
132  !doneEdge[edgeI]
133  && ( addedPoints_[e[0]].size() || addedPoints_[e[1]].size() )
134  )
135  {
136  startFe = fp;
137  break;
138  }
139  }
140 
141  if (startFe != -1)
142  {
143  // We found an edge that needs extruding but hasn't been done yet.
144  // Now find the face on the other side
145  label nbrGlobalFacei = nbrFace
146  (
147  globalEdgeFaces,
148  fEdges[startFe],
149  globalFacei
150  );
151 
152  if (nbrGlobalFacei == -1)
153  {
154  // Proper boundary edge. Only extrude single edge.
155  endFe = startFe;
156  }
157  else
158  {
159  // Search back for edge
160  // - which hasn't been handled yet
161  // - with same neighbour
162  // - that needs extrusion
163  while (true)
164  {
165  label prevFp = fEdges.rcIndex(startFe);
166 
167  if
168  (
169  !sameEdgeNeighbour
170  (
171  pp,
172  globalEdgeFaces,
173  doneEdge,
174  globalFacei,
175  nbrGlobalFacei,
176  fEdges[prevFp]
177  )
178  )
179  {
180  break;
181  }
182  startFe = prevFp;
183  }
184 
185  // Search forward for end of string
186  endFe = startFe;
187  while (true)
188  {
189  label nextFp = fEdges.fcIndex(endFe);
190 
191  if
192  (
193  !sameEdgeNeighbour
194  (
195  pp,
196  globalEdgeFaces,
197  doneEdge,
198  globalFacei,
199  nbrGlobalFacei,
200  fEdges[nextFp]
201  )
202  )
203  {
204  break;
205  }
206  endFe = nextFp;
207  }
208  }
209  }
210 
211  return labelPair(startFe, endFe);
212 }
213 
214 
215 Foam::label Foam::addPatchCellLayer::addSideFace
216 (
217  const polyMesh& mesh,
218  const indirectPrimitivePatch& pp,
219  const labelListList& addedCells, // per pp face the new extruded cell
220  const face& newFace,
221  const label newPatchID,
222 
223  const label ownFacei, // pp face that provides owner
224  const label nbrFacei,
225  const label layerI, // layer
226  const label numEdgeFaces, // number of layers for edge
227  const labelList& meshFaces, // precalculated edgeFaces
228  polyTopoChange& meshMod
229 ) const
230 {
231  label masterFacei = -1;
232  label addedFacei = -1;
233 
234  // Is patch edge external edge of indirectPrimitivePatch?
235  if (nbrFacei == -1)
236  {
237  // External edge so external face.
238 
239  const polyBoundaryMesh& patches = mesh.boundaryMesh();
240 
241  // Loop over all faces connected to edge and see if we can find a face
242  // that is otherPatchID
243 
244  // Get my mesh face
245  label meshFacei = pp.addressing()[ownFacei];
246 
247  forAll(meshFaces, k)
248  {
249  label facei = meshFaces[k];
250 
251  if
252  (
253  (facei != meshFacei)
254  && (patches.whichPatch(facei) == newPatchID)
255  )
256  {
257  // Found the patch face. Use it to map from
258  masterFacei = facei;
259  break;
260  }
261  }
262 
263  // Determine if different number of layer on owner and neighbour side
264  // (relevant only for coupled faces). See section for internal edge
265  // below.
266 
267  label layerOwn;
268 
269  if (addedCells[ownFacei].size() < numEdgeFaces)
270  {
271  label offset = numEdgeFaces - addedCells[ownFacei].size();
272  if (layerI <= offset)
273  {
274  layerOwn = 0;
275  }
276  else
277  {
278  layerOwn = layerI - offset;
279  }
280  }
281  else
282  {
283  layerOwn = layerI;
284  }
285 
286  addedFacei = meshMod.addFace
287  (
288  newFace, // face
289  addedCells[ownFacei][layerOwn], // owner
290  -1, // neighbour
291  masterFacei, // master face
292  false, // flux flip
293  newPatchID // patch for face
294  );
295  }
296  else
297  {
298  // When adding side faces we need to modify neighbour and owners
299  // in region where layer mesh is stopped. Determine which side
300  // has max number of faces and make sure layers match closest to
301  // original pp if there are different number of layers.
302 
303  label layerNbr;
304  label layerOwn;
305 
306  if (addedCells[ownFacei].size() > addedCells[nbrFacei].size())
307  {
308  label offset =
309  addedCells[ownFacei].size() - addedCells[nbrFacei].size();
310 
311  layerOwn = layerI;
312 
313  if (layerI <= offset)
314  {
315  layerNbr = 0;
316  }
317  else
318  {
319  layerNbr = layerI - offset;
320  }
321  }
322  else if (addedCells[nbrFacei].size() > addedCells[ownFacei].size())
323  {
324  label offset =
325  addedCells[nbrFacei].size() - addedCells[ownFacei].size();
326 
327  layerNbr = layerI;
328 
329  if (layerI <= offset)
330  {
331  layerOwn = 0;
332  }
333  else
334  {
335  layerOwn = layerI - offset;
336  }
337  }
338  else
339  {
340  // Same number of layers on both sides.
341  layerNbr = layerI;
342  layerOwn = layerI;
343  }
344 
345  addedFacei = meshMod.addFace
346  (
347  newFace, // face
348  addedCells[ownFacei][layerOwn], // owner
349  addedCells[nbrFacei][layerNbr], // neighbour
350  -1, // master face
351  false, // flux flip
352  -1 // patch for face
353  );
354  }
355 
356  return addedFacei;
357 }
358 
359 
360 Foam::label Foam::addPatchCellLayer::findProcPatch
361 (
362  const polyMesh& mesh,
363  const label nbrProcID
364 )
365 {
366  const polyBoundaryMesh& patches = mesh.boundaryMesh();
367 
369  {
371 
372  if
373  (
374  refCast<const processorPolyPatch>(patches[patchi]).neighbProcNo()
375  == nbrProcID
376  )
377  {
378  return patchi;
379  }
380  }
381  return -1;
382 }
383 
384 
385 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
386 
388 (
389  const polyMesh& mesh,
390  const bool addToMesh
391 )
392 :
393  addToMesh_(addToMesh),
394  addedPoints_(0),
395  pointZonesAddedPoints_(mesh.pointZones().size()),
396  faceZonesAddedFaces_(mesh.faceZones().size()),
397  cellZonesAddedCells_(mesh.cellZones().size()),
398  layerFaces_(0)
399 {}
400 
401 
402 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
403 
405 (
406  const polyMesh& mesh,
407  const labelListList& layerFaces
408 )
409 {
410  labelListList layerCells(layerFaces.size());
411 
412  forAll(layerFaces, patchFacei)
413  {
414  const labelList& faceLabels = layerFaces[patchFacei];
415 
416  if (faceLabels.size())
417  {
418  labelList& added = layerCells[patchFacei];
419  added.setSize(faceLabels.size()-1);
420 
421  for (label i = 0; i < faceLabels.size()-1; i++)
422  {
423  added[i] = mesh.faceNeighbour()[faceLabels[i]];
424  }
425  }
426  }
427  return layerCells;
428 }
429 
430 
432 (
433  const polyMesh& mesh
434 ) const
435 {
436  return addedCells(mesh, layerFaces_);
437 }
438 
439 
441 (
442  const polyMesh& mesh,
443  const globalIndex& globalFaces,
444  const indirectPrimitivePatch& pp
445 )
446 {
447  // Precalculate mesh edges for pp.edges.
448  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
449 
450  // From mesh edge to global face labels. Non-empty sublists only for
451  // pp edges.
452  labelListList globalEdgeFaces(mesh.nEdges());
453 
454  const labelListList& edgeFaces = pp.edgeFaces();
455 
456  forAll(edgeFaces, edgeI)
457  {
458  label meshEdgeI = meshEdges[edgeI];
459 
460  const labelList& eFaces = edgeFaces[edgeI];
461 
462  // Store face and processor as unique tag.
463  labelList& globalEFaces = globalEdgeFaces[meshEdgeI];
464  globalEFaces.setSize(eFaces.size());
465  forAll(eFaces, i)
466  {
467  globalEFaces[i] = globalFaces.toGlobal(pp.addressing()[eFaces[i]]);
468  }
469  }
470 
471  // Synchronise across coupled edges.
473  (
474  mesh,
475  globalEdgeFaces,
476  uniqueEqOp(),
477  labelList() // null value
478  );
479 
480  // Extract pp part
481  return labelListList(UIndirectList<labelList>(globalEdgeFaces, meshEdges));
482 }
483 
484 
486 (
487  const polyMesh& mesh,
488  const globalIndex& globalFaces,
489  const labelListList& globalEdgeFaces,
490  const indirectPrimitivePatch& pp,
491 
492  labelList& sidePatchID,
493  label& nPatches,
494  Map<label>& nbrProcToPatch,
495  Map<label>& patchToNbrProc
496 )
497 {
499 
500  // Precalculate mesh edges for pp.edges.
501  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
502 
503  sidePatchID.setSize(pp.nEdges());
504  sidePatchID = -1;
505 
506  // These also get determined but not (yet) exported:
507  // - whether face is created from other face or edge
508  // - what zone&orientation face should have
509 
510  labelList masterFacei(pp.nEdges(), -1);
511 
512  nPatches = patches.size();
513 
514  forAll(globalEdgeFaces, edgeI)
515  {
516  const labelList& eGlobalFaces = globalEdgeFaces[edgeI];
517  if
518  (
519  eGlobalFaces.size() == 2
520  && pp.edgeFaces()[edgeI].size() == 1
521  )
522  {
523  // Locally but not globally a boundary edge. Hence a coupled
524  // edge. Find the patch to use if on different
525  // processors.
526 
527  label f0 = eGlobalFaces[0];
528  label f1 = eGlobalFaces[1];
529 
530  label otherProci = -1;
531  if (globalFaces.isLocal(f0) && !globalFaces.isLocal(f1))
532  {
533  otherProci = globalFaces.whichProcID(f1);
534  }
535  else if (!globalFaces.isLocal(f0) && globalFaces.isLocal(f1))
536  {
537  otherProci = globalFaces.whichProcID(f0);
538  }
539 
540 
541  if (otherProci != -1)
542  {
543  sidePatchID[edgeI] = findProcPatch(mesh, otherProci);
544  if (sidePatchID[edgeI] == -1)
545  {
546  // Cannot find a patch to processor. See if already
547  // marked for addition
548  if (nbrProcToPatch.found(otherProci))
549  {
550  sidePatchID[edgeI] = nbrProcToPatch[otherProci];
551  }
552  else
553  {
554  sidePatchID[edgeI] = nPatches;
555  nbrProcToPatch.insert(otherProci, nPatches);
556  patchToNbrProc.insert(nPatches, otherProci);
557  nPatches++;
558  }
559  }
560  }
561  }
562  }
563 
564 
565 
566  // Determine face properties for all other boundary edges
567  // ------------------------------------------------------
568 
569  const labelListList& edgeFaces = pp.edgeFaces();
570 
571  DynamicList<label> dynMeshEdgeFaces;
572 
573  forAll(edgeFaces, edgeI)
574  {
575  if (edgeFaces[edgeI].size() == 1 && sidePatchID[edgeI] == -1)
576  {
577  // Proper, uncoupled patch edge.
578 
579  label myFacei = pp.addressing()[edgeFaces[edgeI][0]];
580 
581  // Pick up any boundary face on this edge and use its properties
582  label meshEdgeI = meshEdges[edgeI];
583  const labelList& meshFaces = mesh.edgeFaces
584  (
585  meshEdgeI,
586  dynMeshEdgeFaces
587  );
588 
589  forAll(meshFaces, k)
590  {
591  label facei = meshFaces[k];
592 
593  if (facei != myFacei && !mesh.isInternalFace(facei))
594  {
595  sidePatchID[edgeI] = mesh.boundaryMesh().whichPatch(facei);
596  masterFacei[edgeI] = facei;
597 
598  break;
599  }
600  }
601  }
602  }
603 
604 
605  // Now hopefully every boundary edge has a side patch. Check
606  if (debug)
607  {
608  forAll(edgeFaces, edgeI)
609  {
610  if (edgeFaces[edgeI].size() == 1 && sidePatchID[edgeI] == -1)
611  {
612  const edge& e = pp.edges()[edgeI];
614  << "Have no sidePatchID for edge " << edgeI << " points "
615  << pp.points()[pp.meshPoints()[e[0]]]
616  << pp.points()[pp.meshPoints()[e[1]]]
617  << endl;
618  }
619  }
620  }
621 
622 
623  // Now we have sidepatch see if we have patchface or edge to map from
624  forAll(edgeFaces, edgeI)
625  {
626  if
627  (
628  edgeFaces[edgeI].size() == 1
629  && sidePatchID[edgeI] != -1
630  && masterFacei[edgeI] == -1
631  )
632  {
633  // 1. Do we have a boundary face to map from
634 
635  label myFacei = pp.addressing()[edgeFaces[edgeI][0]];
636 
637  // Pick up any boundary face on this edge and use its properties
638  label meshEdgeI = meshEdges[edgeI];
639  const labelList& meshFaces = mesh.edgeFaces
640  (
641  meshEdgeI,
642  dynMeshEdgeFaces
643  );
644 
645  forAll(meshFaces, k)
646  {
647  const label facei = meshFaces[k];
648 
649  if
650  (
651  facei != myFacei
652  && !mesh.isInternalFace(facei)
653  && patches.whichPatch(facei) == sidePatchID[edgeI]
654  )
655  {
656  sidePatchID[edgeI] = mesh.boundaryMesh().whichPatch(facei);
657  masterFacei[edgeI] = facei;
658  break;
659  }
660  }
661  }
662  }
663 }
664 
665 
667 (
668  const polyMesh& mesh,
669  const globalIndex& globalFaces,
670  const labelListList& globalEdgeFaces,
671  const scalarField& expansionRatio,
672  const indirectPrimitivePatch& pp,
673  const labelList& sidePatchID,
674  const labelList& exposedPatchID,
675  const labelList& nFaceLayers,
676  const labelList& nPointLayers,
677  const vectorField& firstLayerDisp,
678  const labelList& faceCellZones,
679  polyTopoChange& meshMod
680 )
681 {
682  if (debug)
683  {
684  Pout<< "addPatchCellLayer::setRefinement : Adding up to "
685  << gMax(nPointLayers)
686  << " layers of cells to indirectPrimitivePatch with "
687  << pp.nPoints() << " points" << endl;
688  }
689 
690  if
691  (
692  pp.nPoints() != firstLayerDisp.size()
693  || pp.nPoints() != nPointLayers.size()
694  || pp.size() != nFaceLayers.size()
695  )
696  {
698  << "Size of new points is not same as number of points used by"
699  << " the face subset" << endl
700  << " patch.nPoints:" << pp.nPoints()
701  << " displacement:" << firstLayerDisp.size()
702  << " nPointLayers:" << nPointLayers.size() << nl
703  << " patch.nFaces:" << pp.size()
704  << " nFaceLayers:" << nFaceLayers.size()
705  << abort(FatalError);
706  }
707 
708  forAll(nPointLayers, i)
709  {
710  if (nPointLayers[i] < 0)
711  {
713  << "Illegal number of layers " << nPointLayers[i]
714  << " at patch point " << i << abort(FatalError);
715  }
716  }
717  forAll(nFaceLayers, i)
718  {
719  if (nFaceLayers[i] < 0)
720  {
722  << "Illegal number of layers " << nFaceLayers[i]
723  << " at patch face " << i << abort(FatalError);
724  }
725  }
726 
727  forAll(globalEdgeFaces, edgeI)
728  {
729  if (globalEdgeFaces[edgeI].size() > 2)
730  {
731  const edge& e = pp.edges()[edgeI];
732 
733  if (nPointLayers[e[0]] > 0 || nPointLayers[e[1]] > 0)
734  {
736  << "Trying to extrude edge "
737  << e.line(pp.localPoints())
738  << " which is non-manifold (has "
739  << globalEdgeFaces[edgeI].size()
740  << " faces using it)"
741  << abort(FatalError);
742  }
743  }
744  }
745 
746 
747  const labelList& meshPoints = pp.meshPoints();
748 
749  // Some storage for edge-face-addressing.
751 
752  // Precalculate mesh edges for pp.edges.
753  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
754 
755  if (debug)
756  {
757  // Check synchronisation
758  // ~~~~~~~~~~~~~~~~~~~~~
759 
760  {
761  labelList n(mesh.nPoints(), 0);
762  UIndirectList<label>(n, meshPoints) = nPointLayers;
764 
765  // Non-synced
766  forAll(meshPoints, i)
767  {
768  label meshPointi = meshPoints[i];
769 
770  if (n[meshPointi] != nPointLayers[i])
771  {
773  << "At mesh point:" << meshPointi
774  << " coordinate:" << mesh.points()[meshPointi]
775  << " specified nLayers:" << nPointLayers[i] << endl
776  << "On coupled point a different nLayers:"
777  << n[meshPointi] << " was specified."
778  << abort(FatalError);
779  }
780  }
781 
782 
783  // Check that nPointLayers equals the max layers of connected faces
784  // (or 0). Anything else makes no sense.
785  labelList nFromFace(mesh.nPoints(), 0);
786  forAll(nFaceLayers, i)
787  {
788  const face& f = pp[i];
789 
790  forAll(f, fp)
791  {
792  label pointi = f[fp];
793 
794  nFromFace[pointi] = max(nFromFace[pointi], nFaceLayers[i]);
795  }
796  }
798  (
799  mesh,
800  nFromFace,
801  maxEqOp<label>(),
802  label(0)
803  );
804 
805  forAll(nPointLayers, i)
806  {
807  label meshPointi = meshPoints[i];
808 
809  if
810  (
811  nPointLayers[i] > 0
812  && nPointLayers[i] != nFromFace[meshPointi]
813  )
814  {
816  << "At mesh point:" << meshPointi
817  << " coordinate:" << mesh.points()[meshPointi]
818  << " specified nLayers:" << nPointLayers[i] << endl
819  << "but the max nLayers of surrounding faces is:"
820  << nFromFace[meshPointi]
821  << abort(FatalError);
822  }
823  }
824  }
825 
826  {
828  UIndirectList<point>(d, meshPoints) = firstLayerDisp;
830  (
831  mesh,
832  d,
833  minEqOp<vector>(),
835  );
836 
837  forAll(meshPoints, i)
838  {
839  label meshPointi = meshPoints[i];
840 
841  if (mag(d[meshPointi] - firstLayerDisp[i]) > small)
842  {
844  << "At mesh point:" << meshPointi
845  << " coordinate:" << mesh.points()[meshPointi]
846  << " specified displacement:" << firstLayerDisp[i]
847  << endl
848  << "On coupled point a different displacement:"
849  << d[meshPointi] << " was specified."
850  << abort(FatalError);
851  }
852  }
853  }
854 
855  // Check that edges of pp (so ones that become boundary faces)
856  // connect to only one boundary face. Guarantees uniqueness of
857  // patch that they go into so if this is a coupled patch both
858  // sides decide the same.
859  // ~~~~~~~~~~~~~~~~~~~~~~
860 
861  for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
862  {
863  const edge& e = pp.edges()[edgeI];
864 
865  if (nPointLayers[e[0]] > 0 || nPointLayers[e[1]] > 0)
866  {
867  // Edge is to become a face
868 
869  const labelList& eFaces = pp.edgeFaces()[edgeI];
870 
871  // First check: pp should be single connected.
872  if (eFaces.size() != 1)
873  {
875  << "boundary-edge-to-be-extruded:"
876  << pp.points()[meshPoints[e[0]]]
877  << pp.points()[meshPoints[e[1]]]
878  << " has more than two faces using it:" << eFaces
879  << abort(FatalError);
880  }
881 
882  label myFacei = pp.addressing()[eFaces[0]];
883 
884  label meshEdgeI = meshEdges[edgeI];
885 
886  // Mesh faces using edge
887  const labelList& meshFaces = mesh.edgeFaces(meshEdgeI, ef);
888 
889  // Check that there is only one patchface using edge.
891 
892  label bFacei = -1;
893 
894  forAll(meshFaces, i)
895  {
896  label facei = meshFaces[i];
897 
898  if (facei != myFacei)
899  {
900  if (!mesh.isInternalFace(facei))
901  {
902  if (bFacei == -1)
903  {
904  bFacei = facei;
905  }
906  else
907  {
909  << "boundary-edge-to-be-extruded:"
910  << pp.points()[meshPoints[e[0]]]
911  << pp.points()[meshPoints[e[1]]]
912  << " has more than two boundary faces"
913  << " using it:"
914  << bFacei << " fc:"
915  << mesh.faceCentres()[bFacei]
916  << " patch:" << patches.whichPatch(bFacei)
917  << " and " << facei << " fc:"
918  << mesh.faceCentres()[facei]
919  << " patch:" << patches.whichPatch(facei)
920  << abort(FatalError);
921  }
922  }
923  }
924  }
925  }
926  }
927  }
928 
929 
931 
932  // Precalculated patchID for each patch face
933  labelList patchID(pp.size());
934 
935  forAll(pp, patchFacei)
936  {
937  label meshFacei = pp.addressing()[patchFacei];
938 
939  patchID[patchFacei] = patches.whichPatch(meshFacei);
940  }
941 
942 
943  // From master point (in patch point label) to added points (in mesh point
944  // label)
945  addedPoints_.setSize(pp.nPoints());
946 
947  // Mark points that do not get extruded by setting size of addedPoints_ to 0
948  label nTruncated = 0;
949 
950  forAll(nPointLayers, patchPointi)
951  {
952  if (nPointLayers[patchPointi] > 0)
953  {
954  addedPoints_[patchPointi].setSize(nPointLayers[patchPointi]);
955  }
956  else
957  {
958  nTruncated++;
959  }
960  }
961 
962  if (debug)
963  {
964  Pout<< "Not adding points at " << nTruncated << " out of "
965  << pp.nPoints() << " points" << endl;
966  }
967 
968 
969  //
970  // Create new points
971  //
972 
973  // If creating new mesh: copy existing patch points
974  labelList copiedPatchPoints;
975  if (!addToMesh_)
976  {
977  copiedPatchPoints.setSize(firstLayerDisp.size());
978  forAll(firstLayerDisp, patchPointi)
979  {
980  if (addedPoints_[patchPointi].size())
981  {
982  const label meshPointi = meshPoints[patchPointi];
983 
984  copiedPatchPoints[patchPointi] = meshMod.addPoint
985  (
986  mesh.points()[meshPointi], // point
987  -1, // master point
988  true // supports a cell
989  );
990 
991  forAll(mesh.pointZones(), zonei)
992  {
993  if (mesh.pointZones()[zonei].localIndex(meshPointi) != -1)
994  {
995  pointZonesAddedPoints_[zonei].insert
996  (
997  copiedPatchPoints[patchPointi]
998  );
999  }
1000  }
1001  }
1002  }
1003  }
1004 
1005 
1006  // Create points for additional layers
1007  forAll(firstLayerDisp, patchPointi)
1008  {
1009  if (addedPoints_[patchPointi].size())
1010  {
1011  const label meshPointi = meshPoints[patchPointi];
1012 
1013  point pt = mesh.points()[meshPointi];
1014  vector disp = firstLayerDisp[patchPointi];
1015 
1016  const labelList zones(mesh.pointZones().whichZones(meshPointi));
1017 
1018  forAll(addedPoints_[patchPointi], i)
1019  {
1020  pt += disp;
1021 
1022  const label addedVertI = meshMod.addPoint
1023  (
1024  pt, // point
1025  (addToMesh_ ? meshPointi : -1), // master point
1026  true // supports a cell
1027  );
1028 
1029  addedPoints_[patchPointi][i] = addedVertI;
1030 
1031  forAll(zones, zonei)
1032  {
1033  pointZonesAddedPoints_[zonei].insert(addedVertI);
1034  }
1035 
1036  disp *= expansionRatio[patchPointi];
1037  }
1038  }
1039  }
1040 
1041 
1042  //
1043  // Add cells to all boundaryFaces
1044  //
1045 
1046  labelListList addedCells(pp.size());
1047 
1048  forAll(pp, patchFacei)
1049  {
1050  if (nFaceLayers[patchFacei] > 0)
1051  {
1052  addedCells[patchFacei].setSize(nFaceLayers[patchFacei]);
1053 
1054  const label meshFacei = pp.addressing()[patchFacei];
1055 
1056  const labelList ownZones
1057  (
1058  mesh.cellZones().whichZones(mesh.faceOwner()[meshFacei])
1059  );
1060 
1061  for (label i = 0; i < nFaceLayers[patchFacei]; i++)
1062  {
1063  // Note: add from cell (owner of patch face) or from face?
1064  // for now add from cell so we can map easily.
1065  addedCells[patchFacei][i] = meshMod.addCell
1066  (
1067  (addToMesh_ ? mesh.faceOwner()[meshFacei] : -1)
1068  );
1069 
1070  forAll(ownZones, zonei)
1071  {
1072  cellZonesAddedCells_[ownZones[zonei]].insert
1073  (
1074  addedCells[patchFacei][i]
1075  );
1076  }
1077 
1078  if (faceCellZones.size() && faceCellZones[patchFacei] != -1)
1079  {
1080  cellZonesAddedCells_[faceCellZones[patchFacei]].insert
1081  (
1082  addedCells[patchFacei][i]
1083  );
1084  }
1085  }
1086  }
1087  }
1088 
1089 
1090  // Create faces on top of the original patch faces.
1091  // These faces are created from original patch faces outwards so in order
1092  // of increasing cell number. So orientation should be same as original
1093  // patch face for them to have owner<neighbour.
1094 
1095  layerFaces_.setSize(pp.size());
1096 
1097  forAll(pp.localFaces(), patchFacei)
1098  {
1099  label meshFacei = pp.addressing()[patchFacei];
1100 
1101  if (addedCells[patchFacei].size())
1102  {
1103  layerFaces_[patchFacei].setSize(addedCells[patchFacei].size() + 1);
1104 
1105  // Get duplicated vertices on the patch face.
1106  const face& f = pp.localFaces()[patchFacei];
1107 
1108  face newFace(f.size());
1109 
1110  forAll(addedCells[patchFacei], i)
1111  {
1112  forAll(f, fp)
1113  {
1114  if (addedPoints_[f[fp]].empty())
1115  {
1116  // Keep original point
1117  newFace[fp] =
1118  (
1119  addToMesh_
1120  ? meshPoints[f[fp]]
1121  : copiedPatchPoints[f[fp]]
1122  );
1123  }
1124  else
1125  {
1126  // Get new outside point
1127  label offset =
1128  addedPoints_[f[fp]].size()
1129  - addedCells[patchFacei].size();
1130  newFace[fp] = addedPoints_[f[fp]][i+offset];
1131  }
1132  }
1133 
1134 
1135  // Get new neighbour
1136  label nei;
1137  label patchi;
1138 
1139  if (i == addedCells[patchFacei].size()-1)
1140  {
1141  // Top layer so is patch face.
1142  nei = -1;
1143  patchi = patchID[patchFacei];
1144  }
1145  else
1146  {
1147  // Internal face between layer i and i+1
1148  nei = addedCells[patchFacei][i+1];
1149  patchi = -1;
1150  }
1151 
1152  layerFaces_[patchFacei][i+1] = meshMod.addFace
1153  (
1154  newFace, // face
1155  addedCells[patchFacei][i], // owner
1156  nei, // neighbour
1157  (addToMesh_ ? meshFacei : -1), // master face
1158  false, // flux flip
1159  patchi // patch for face
1160  );
1161  }
1162  }
1163  }
1164 
1165 
1166  //
1167  // Modify old patch faces to be on the inside
1168  //
1169 
1170  if (addToMesh_)
1171  {
1172  forAll(pp, patchFacei)
1173  {
1174  if (addedCells[patchFacei].size())
1175  {
1176  label meshFacei = pp.addressing()[patchFacei];
1177 
1178  layerFaces_[patchFacei][0] = meshFacei;
1179 
1180  meshMod.modifyFace
1181  (
1182  pp[patchFacei], // modified face
1183  meshFacei, // label of face
1184  mesh.faceOwner()[meshFacei], // owner
1185  addedCells[patchFacei][0], // neighbour
1186  false, // face flip
1187  -1 // patch for face
1188  );
1189  }
1190  }
1191  }
1192  else
1193  {
1194  // If creating new mesh: reverse original faces and put them
1195  // in the exposed patch ID.
1196  forAll(pp, patchFacei)
1197  {
1198  if (nFaceLayers[patchFacei] > 0)
1199  {
1200  // Reverse and renumber old patch face.
1201  face f(pp.localFaces()[patchFacei].reverseFace());
1202  forAll(f, fp)
1203  {
1204  f[fp] = copiedPatchPoints[f[fp]];
1205  }
1206 
1207  layerFaces_[patchFacei][0] = meshMod.addFace
1208  (
1209  f, // modified face
1210  addedCells[patchFacei][0], // owner
1211  -1, // neighbour
1212  -1, // masterFace
1213  true, // face flip
1214  exposedPatchID[patchFacei] // patch for face
1215  );
1216  }
1217  }
1218  }
1219 
1220 
1221  //
1222  // Create 'side' faces, one per edge that is being extended.
1223  //
1224 
1225  // Get number of layers per edge. This is 0 if edge is not extruded;
1226  // max of connected faces otherwise.
1227  labelList edgeLayers(pp.nEdges());
1228 
1229  forAll(edgeLayers, ppEdgei)
1230  {
1231  edgeLayers[ppEdgei] = 0;
1232 
1233  const edge& e = pp.edges()[ppEdgei];
1234 
1235  if (nPointLayers[e[0]] != 0 || nPointLayers[e[1]] != 0)
1236  {
1237  const labelList& eFaces = pp.edgeFaces()[ppEdgei];
1238 
1239  forAll(eFaces, eFacei)
1240  {
1241  edgeLayers[ppEdgei] = max
1242  (
1243  nFaceLayers[eFaces[eFacei]],
1244  edgeLayers[ppEdgei]
1245  );
1246  }
1247  }
1248  }
1249 
1251  (
1252  mesh,
1253  meshEdges,
1254  edgeLayers,
1255  maxEqOp<label>(),
1256  label(0)
1257  );
1258 
1259 
1260  // Determine the face zones that logically extrude from each patch edge
1261  List<labelList> ppEdgeFaceZones(pp.nEdges());
1262  {
1263  // Create a map from mesh face to patch face
1264  labelList meshFacePpFace(mesh.nFaces(), -1);
1265  UIndirectList<label>(meshFacePpFace, pp.addressing()) =
1266  identityMap(pp.size());
1267 
1268  // Do each patch edge in turn
1269  forAll(pp.edges(), ppEdgei)
1270  {
1271  const label meshEdgei = meshEdges[ppEdgei];
1272  const edge& meshEdge = mesh.edges()[meshEdgei];
1273 
1274  const labelList& edgeMeshFaces = mesh.edgeFaces(meshEdgei, ef);
1275 
1276  // Determine the face zones that connect to this edge
1277  DynamicList<label> edgeFaceZones;
1278  forAll(edgeMeshFaces, edgeFacei)
1279  {
1280  const label meshFacei = edgeMeshFaces[edgeFacei];
1281  const face& meshFace = mesh.faces()[meshFacei];
1282 
1283  // Ignore faces that are actually part of the patch
1284  if (meshFacePpFace[meshFacei] != -1) continue;
1285 
1286  const labelList faceFaceZones =
1287  mesh.faceZones().whichZones(meshFacei);
1288 
1289  forAll(faceFaceZones, faceFaceZonei)
1290  {
1291  const label faceZonei = faceFaceZones[faceFaceZonei];
1292  const faceZone& faceZone = mesh.faceZones()[faceZonei];
1293 
1294  // Ignore anything that isn't a surface
1295  if (!faceZone.oriented()) continue;
1296 
1297  const label faceZoneFacei = faceZone.localIndex(meshFacei);
1298 
1299  const label direction =
1300  (faceZone.flipMap()[faceZoneFacei] ? -1 : +1)
1301  *meshFace.edgeDirection(meshEdge);
1302 
1303  edgeFaceZones.append(direction*(1 + faceZonei));
1304  }
1305  }
1306 
1307  // Copy into the main list
1308  ppEdgeFaceZones[ppEdgei].transfer(edgeFaceZones);
1309  }
1310 
1311  // Synchronise
1312  if (Pstream::parRun())
1313  {
1315  (
1316  mesh,
1317  meshEdges,
1318  ppEdgeFaceZones,
1320  labelList()
1321  );
1322  }
1323 
1324  // Have differently oriented faces "cancel out"
1325  forAll(ppEdgeFaceZones, ppEdgei)
1326  {
1327  if (ppEdgeFaceZones[ppEdgei].empty()) continue;
1328 
1329  labelList faceZoneCount(mesh.faceZones().size(), label(0));
1330 
1331  forAll(ppEdgeFaceZones[ppEdgei], ppEdgeFaceZonei)
1332  {
1333  const label turnedZonei =
1334  ppEdgeFaceZones[ppEdgei][ppEdgeFaceZonei];
1335 
1336  faceZoneCount[mag(turnedZonei) - 1] += sign(turnedZonei);
1337  }
1338 
1339  DynamicList<label> edgeFaceZones(ppEdgeFaceZones[ppEdgei].size());
1340 
1341  forAll(faceZoneCount, faceZonei)
1342  {
1343  if (faceZoneCount[faceZonei] != 0)
1344  {
1345  edgeFaceZones.append
1346  (
1347  sign(faceZoneCount[faceZonei])*(faceZonei + 1)
1348  );
1349  }
1350  }
1351 
1352  ppEdgeFaceZones[ppEdgei].transfer(edgeFaceZones);
1353  }
1354  }
1355 
1356 
1357  // Mark off which edges have been extruded
1358  boolList doneEdge(pp.nEdges(), false);
1359 
1360 
1361  // Create faces. Per face walk connected edges and find string of edges
1362  // between the same two faces and extrude string into a single face.
1363  forAll(pp, patchFacei)
1364  {
1365  const labelList& fEdges = pp.faceEdges()[patchFacei];
1366 
1367  forAll(fEdges, fp)
1368  {
1369  // Get string of edges that needs to be extruded as a single face.
1370  // Returned as indices in fEdges.
1371  const labelPair startAndEndFes
1372  (
1373  getEdgeString
1374  (
1375  pp,
1376  globalEdgeFaces,
1377  doneEdge,
1378  patchFacei,
1379  globalFaces.toGlobal(pp.addressing()[patchFacei])
1380  )
1381  );
1382  const label startFe = startAndEndFes[0];
1383  const label endFe = startAndEndFes[1];
1384 
1385  if (startFe != -1)
1386  {
1387  const face& f = pp.localFaces()[patchFacei];
1388 
1389  // Create a string of vertices to extrude into a face
1390  labelList stringedVerts
1391  (
1392  (endFe + f.size() - startFe + 2) % f.size()
1393  );
1394  label fp = startFe;
1395  for (label i = 0; i < stringedVerts.size() - 1; i ++)
1396  {
1397  stringedVerts[i] = f[fp];
1398  doneEdge[fEdges[fp]] = true;
1399  fp = f.fcIndex(fp);
1400  }
1401  stringedVerts.last() = f[fp];
1402 
1403  // Now stringedVerts contains the vertices in order of face f.
1404  // This is consistent with the order if f becomes the owner cell
1405  // and nbrFacei the neighbour cell. Note that the cells get
1406  // added in order of pp so we can just use face ordering and
1407  // because we loop in incrementing order as well we will
1408  // always have nbrFacei > patchFacei.
1409 
1410  const label startEdgei = fEdges[startFe];
1411 
1412  const label meshEdgei = meshEdges[startEdgei];
1413 
1414  const label numEdgeSideFaces = edgeLayers[startEdgei];
1415 
1416  for (label i = 0; i < numEdgeSideFaces; i++)
1417  {
1418  label vEnd = stringedVerts.last();
1419  label vStart = stringedVerts[0];
1420 
1421  // calculate number of points making up a face
1422  label newFp = 2*stringedVerts.size();
1423 
1424  if (i == 0)
1425  {
1426  // layer 0 gets all the truncation of neighbouring
1427  // faces with more layers.
1428  if (addedPoints_[vEnd].size())
1429  {
1430  newFp +=
1431  addedPoints_[vEnd].size() - numEdgeSideFaces;
1432  }
1433  if (addedPoints_[vStart].size())
1434  {
1435  newFp +=
1436  addedPoints_[vStart].size() - numEdgeSideFaces;
1437  }
1438  }
1439 
1440  face newFace(newFp);
1441 
1442  newFp = 0;
1443 
1444  // For layer 0 get pp points, for all other layers get
1445  // points of layer-1.
1446  if (i == 0)
1447  {
1448  forAll(stringedVerts, stringedI)
1449  {
1450  label v = stringedVerts[stringedI];
1451  addVertex
1452  (
1453  (
1454  addToMesh_
1455  ? meshPoints[v]
1456  : copiedPatchPoints[v]
1457  ),
1458  newFace,
1459  newFp
1460  );
1461  }
1462  }
1463  else
1464  {
1465  forAll(stringedVerts, stringedI)
1466  {
1467  label v = stringedVerts[stringedI];
1468  if (addedPoints_[v].size())
1469  {
1470  label offset =
1471  addedPoints_[v].size() - numEdgeSideFaces;
1472  addVertex
1473  (
1474  addedPoints_[v][i+offset-1],
1475  newFace,
1476  newFp
1477  );
1478  }
1479  else
1480  {
1481  addVertex
1482  (
1483  (
1484  addToMesh_
1485  ? meshPoints[v]
1486  : copiedPatchPoints[v]
1487  ),
1488  newFace,
1489  newFp
1490  );
1491  }
1492  }
1493  }
1494 
1495  // add points between stringed vertices (end)
1496  if (numEdgeSideFaces < addedPoints_[vEnd].size())
1497  {
1498  if (i == 0 && addedPoints_[vEnd].size())
1499  {
1500  label offset =
1501  addedPoints_[vEnd].size() - numEdgeSideFaces;
1502  for (label ioff = 0; ioff < offset; ioff++)
1503  {
1504  addVertex
1505  (
1506  addedPoints_[vEnd][ioff],
1507  newFace,
1508  newFp
1509  );
1510  }
1511  }
1512  }
1513 
1514  forAllReverse(stringedVerts, stringedI)
1515  {
1516  label v = stringedVerts[stringedI];
1517  if (addedPoints_[v].size())
1518  {
1519  label offset =
1520  addedPoints_[v].size() - numEdgeSideFaces;
1521  addVertex
1522  (
1523  addedPoints_[v][i+offset],
1524  newFace,
1525  newFp
1526  );
1527  }
1528  else
1529  {
1530  addVertex
1531  (
1532  (
1533  addToMesh_
1534  ? meshPoints[v]
1535  : copiedPatchPoints[v]
1536  ),
1537  newFace,
1538  newFp
1539  );
1540  }
1541  }
1542 
1543  // add points between stringed vertices (start)
1544  if (numEdgeSideFaces < addedPoints_[vStart].size())
1545  {
1546  if (i == 0 && addedPoints_[vStart].size())
1547  {
1548  label offset =
1549  addedPoints_[vStart].size() - numEdgeSideFaces;
1550  for (label ioff = offset-1; ioff >= 0; ioff--)
1551  {
1552  addVertex
1553  (
1554  addedPoints_[vStart][ioff],
1555  newFace,
1556  newFp
1557  );
1558  }
1559  }
1560  }
1561 
1562  if (newFp >= 3)
1563  {
1564  // Add face in between faces patchFacei and nbrFacei
1565  // (possibly -1 for external edges)
1566 
1567  newFace.setSize(newFp);
1568 
1569  if (debug)
1570  {
1571  labelHashSet verts(2*newFace.size());
1572  forAll(newFace, fp)
1573  {
1574  if (!verts.insert(newFace[fp]))
1575  {
1577  << "Duplicate vertex in face"
1578  << " to be added." << nl
1579  << "newFace:" << newFace << nl
1580  << "points:"
1582  (
1583  meshMod.points(),
1584  newFace
1585  ) << nl
1586  << "Layer:" << i
1587  << " out of:" << numEdgeSideFaces << nl
1588  << "ExtrudeEdge:" << meshEdgei
1589  << " at:"
1590  << mesh.edges()[meshEdgei].line
1591  (
1592  mesh.points()
1593  ) << nl
1594  << "string:" << stringedVerts
1595  << "stringpoints:"
1597  (
1598  pp.localPoints(),
1599  stringedVerts
1600  ) << nl
1601  << "stringNLayers:"
1603  (
1604  nPointLayers,
1605  stringedVerts
1606  ) << nl
1607  << abort(FatalError);
1608  }
1609  }
1610  }
1611 
1612  const label nbrFacei =
1613  nbrFace(pp.edgeFaces(), startEdgei, patchFacei);
1614 
1615  const labelList& meshFaces =
1616  mesh.edgeFaces(meshEdgei, ef);
1617 
1618  const label addedFacei =
1619  addSideFace
1620  (
1621  mesh,
1622  pp,
1623  addedCells,
1624 
1625  newFace, // vertices of new face
1626  sidePatchID[startEdgei],// -1 or patch for face
1627 
1628  patchFacei,
1629  nbrFacei,
1630  i, // layer
1631  numEdgeSideFaces, // num layers
1632  meshFaces, // edgeFaces
1633  meshMod
1634  );
1635 
1636  // Propagate zones to the side face if they are
1637  // consistent across the entire string of edges
1638  auto getFaceEdgeFaceZones = [&](const label fe)
1639  {
1640  labelHashSet result;
1641 
1642  const label edgei = fEdges[fe];
1643 
1644  forAll(ppEdgeFaceZones[edgei], edgeZonei)
1645  {
1646  result.insert
1647  (
1648  ppEdgeFaceZones[edgei][edgeZonei]
1649  *edge::compare
1650  (
1651  mesh.edges()[meshEdgei],
1652  pp[patchFacei].faceEdge(fe)
1653  )
1654  );
1655  }
1656 
1657  return result;
1658  };
1659 
1660  labelHashSet edgeStringFaceZones =
1661  getFaceEdgeFaceZones(startFe);
1662  for
1663  (
1664  label fe = startFe;
1665  fe != endFe;
1666  fe = f.fcIndex(fe)
1667  )
1668  {
1669  edgeStringFaceZones &=
1670  getFaceEdgeFaceZones(f.fcIndex(fe));
1671  }
1672 
1673  forAllConstIter(labelHashSet, edgeStringFaceZones, iter)
1674  {
1675  faceZonesAddedFaces_[mag(iter.key()) - 1].insert
1676  (
1677  addedFacei,
1678  iter.key() > 0
1679  );
1680  }
1681  }
1682  }
1683  }
1684  }
1685  }
1686 }
1687 
1688 
1690 (
1691  const polyTopoChangeMap& map,
1692  const labelList& faceMap, // new to old patch faces
1693  const labelList& pointMap // new to old patch points
1694 )
1695 {
1696  {
1697  labelListList newAddedPoints(pointMap.size());
1698 
1699  forAll(newAddedPoints, newPointi)
1700  {
1701  label oldPointi = pointMap[newPointi];
1702 
1703  const labelList& added = addedPoints_[oldPointi];
1704 
1705  labelList& newAdded = newAddedPoints[newPointi];
1706  newAdded.setSize(added.size());
1707  label newI = 0;
1708 
1709  forAll(added, i)
1710  {
1711  label newPointi = map.reversePointMap()[added[i]];
1712 
1713  if (newPointi >= 0)
1714  {
1715  newAdded[newI++] = newPointi;
1716  }
1717  }
1718  newAdded.setSize(newI);
1719  }
1720  addedPoints_.transfer(newAddedPoints);
1721  }
1722 
1723  {
1724  labelListList newLayerFaces(faceMap.size());
1725 
1726  forAll(newLayerFaces, newFacei)
1727  {
1728  label oldFacei = faceMap[newFacei];
1729 
1730  const labelList& added = layerFaces_[oldFacei];
1731 
1732  labelList& newAdded = newLayerFaces[newFacei];
1733  newAdded.setSize(added.size());
1734  label newI = 0;
1735 
1736  forAll(added, i)
1737  {
1738  label newFacei = map.reverseFaceMap()[added[i]];
1739 
1740  if (newFacei >= 0)
1741  {
1742  newAdded[newI++] = newFacei;
1743  }
1744  }
1745  newAdded.setSize(newI);
1746  }
1747  layerFaces_.transfer(newLayerFaces);
1748  }
1749 }
1750 
1751 
1753 {
1754  // Add the pointZones to the merged mesh
1755  mesh.pointZones().insert(pointZonesAddedPoints_);
1756 
1757  // Add the faceZones to the merged mesh
1758  mesh.faceZones().insert(faceZonesAddedFaces_);
1759 
1760  // Add the cellZones to the merged mesh
1761  mesh.cellZones().insert(cellZonesAddedCells_);
1762 }
1763 
1764 
1765 // ************************************************************************* //
label k
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:445
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
label nEdges() const
Return number of edges in patch.
label nPoints() const
Return number of points supporting patch faces.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const Field< PointType > & points() const
Return reference to global points.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const labelListList & edgeFaces() const
Return edge-face addressing.
const labelListList & faceEdges() const
Return face-edge addressing.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the list.
const List< label > & addressing() const
Return the list addressing.
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
T & last()
Return the last element of the list.
Definition: UListI.H:128
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const Form max
Definition: VectorSpace.H:120
void insert(const List< labelHashSet > &zonesIndices)
Insert given indices into zones.
Definition: ZoneList.C:378
labelList whichZones(const label objectIndex) const
Given a global object index, return the list of zones it is in.
Definition: ZoneList.C:159
label localIndex(const label globalIndex) const
Map storing the local index for every global index. Used to find.
Definition: Zone.C:221
Adds layers of cells to outside of polyPatch. Can optionally create stand-alone extruded mesh (addToM...
static void calcSidePatch(const polyMesh &, const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const indirectPrimitivePatch &pp, labelList &sidePatchID, label &nPatches, Map< label > &nbrProcToPatch, Map< label > &patchToNbrProc)
Boundary edges get extruded into boundary faces. Determine patch.
static labelListList addedCells(const polyMesh &, const labelListList &layerFaces)
Helper: get added cells per patch face.
void updateZones(polyMesh &mesh)
Update the mesh zones.
void topoChange(const polyTopoChangeMap &, const labelList &faceMap, const labelList &pointMap)
Update any locally stored mesh information. Gets additional.
void setRefinement(const polyMesh &mesh, const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const scalarField &expansionRatio, const indirectPrimitivePatch &pp, const labelList &sidePatchID, const labelList &exposedPatchID, const labelList &nFaceLayers, const labelList &nPointLayers, const vectorField &firstLayerDisp, const labelList &faceCellZones, polyTopoChange &meshMod)
Play commands into polyTopoChange to create layers on top.
addPatchCellLayer(const polyMesh &, const bool addToMesh=true)
Construct from mesh.
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Per patch edge the pp faces (in global indices) using it. Uses.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
static int compare(const edge &, const edge &)
Compare edges.
Definition: edgeI.H:36
void insert(const List< Map< bool >> &zonesOrientedIndices)
Insert given oriented indices.
Definition: faceZoneList.C:45
Named list of face indices representing a sub-set of the mesh faces.
Definition: faceZone.H:66
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:262
bool oriented() const
Return true if the faceZone is oriented, i.e. the flipMap is set.
Definition: faceZone.H:256
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
int edgeDirection(const edge &) const
Return the edge direction on the face.
Definition: face.C:424
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label whichProcID(const label i) const
Which processor does global come from? Binary search.
Definition: globalIndexI.H:123
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:95
const labelList & processorPatches() const
Return list of processor patch labels.
Foam::polyBoundaryMesh.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:437
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelList & reversePointMap() const
Reverse point map.
const labelList & reverseFaceMap() const
Reverse face map.
Direct mesh changes based on v1.3 polyTopoChange syntax.
const DynamicList< point > & points() const
Points. Shrunk after constructing mesh (or calling of compact())
label addCell(const label masterCellID)
Add cell and return new cell index.
label addPoint(const point &, const label masterPointID, const bool inCell)
Add point and return new point index.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
label addFace(const face &f, const label own, const label nei, const label masterFaceID, const bool flipFaceFlux, const label patchID)
Add face to cells and return new face index.
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
label nPoints() const
const labelListList & edgeFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
label nEdges() const
label nFaces() const
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
const doubleScalar e
Definition: doubleScalar.H:106
List< label > labelList
A List of labels.
Definition: labelList.H:56
dimensionedScalar sign(const dimensionedScalar &ds)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
edge meshEdge(const PrimitivePatch< FaceList, PointField > &p, const label edgei)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
void offset(label &lst, const label o)
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
Foam::indirectPrimitivePatch.
static const char nl
Definition: Ostream.H:267
Type gMax(const FieldField< Field, Type > &f)
uint8_t direction
Definition: direction.H:45
label newPointi
Definition: readKivaGrid.H:495
labelList f(nPoints)
label nPatches
Definition: readKivaGrid.H:396
List operator to append one list onto another.
Definition: ListOps.H:304