addPatchCellLayer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-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 "addPatchCellLayer.H"
30 #include "polyMesh.H"
31 #include "polyTopoChange.H"
32 #include "meshTools.H"
33 #include "mapPolyMesh.H"
34 #include "syncTools.H"
35 #include "polyAddPoint.H"
36 #include "polyAddFace.H"
37 #include "polyModifyFace.H"
38 #include "polyAddCell.H"
39 #include "globalIndex.H"
40 #include "PatchTools.H"
41 #include "dummyTransform.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
48 
49  // Reduction class to get minimum value over face.
50  class minEqOpFace
51  {
52  public:
53 
54  void operator()(face& x, const face& y) const
55  {
56  if (x.size())
57  {
58  if (y.size())
59  {
60  label j = 0;
61  forAll(x, i)
62  {
63  x[i] = min(x[i], y[j]);
64 
65  j = y.rcIndex(j);
66  }
67  }
68  }
69  else if (y.size())
70  {
71  x.setSize(y.size());
72  label j = 0;
73  forAll(x, i)
74  {
75  x[i] = y[j];
76  j = y.rcIndex(j);
77  }
78  }
79  }
80  };
81 
82 }
83 
84 
85 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
86 
87 Foam::label Foam::addPatchCellLayer::nbrFace
88 (
89  const labelListList& edgeFaces,
90  const label edgei,
91  const label facei
92 )
93 {
94  const labelList& eFaces = edgeFaces[edgei];
95 
96  if (eFaces.size() == 2)
97  {
98  return (eFaces[0] != facei ? eFaces[0] : eFaces[1]);
99  }
100  else
101  {
102  return -1;
103  }
104 }
105 
106 
107 void Foam::addPatchCellLayer::addVertex
108 (
109  const label pointi,
110  face& f,
111  label& fp
112 )
113 {
114  if (fp == 0)
115  {
116  f[fp++] = pointi;
117  }
118  else
119  {
120  if (f[fp-1] != pointi && f[0] != pointi)
121  {
122  f[fp++] = pointi;
123  }
124  }
125 }
126 
127 
128 // Is edge to the same neighbour? (and needs extrusion and has not been
129 // dealt with already)
130 bool Foam::addPatchCellLayer::sameEdgeNeighbour
131 (
132  const indirectPrimitivePatch& pp,
133  const labelListList& globalEdgeFaces,
134  const boolList& doneEdge,
135  const label thisGlobalFacei,
136  const label nbrGlobalFacei,
137  const label edgei
138 ) const
139 {
140  const edge& e = pp.edges()[edgei];
141 
142  return
143  !doneEdge[edgei] // not yet handled
144  && (
145  addedPoints_[e[0]].size() // is extruded
146  || addedPoints_[e[1]].size()
147  )
148  && (
149  nbrFace(globalEdgeFaces, edgei, thisGlobalFacei)
150  == nbrGlobalFacei // is to same neighbour
151  );
152 }
153 
154 
155 // Collect consecutive string of edges that connects the same two
156 // (possibly coupled) faces. Returns -1 if no unvisited edge can be found.
157 // Otherwise returns start and end index in face.
158 Foam::labelPair Foam::addPatchCellLayer::getEdgeString
159 (
160  const indirectPrimitivePatch& pp,
161  const labelListList& globalEdgeFaces,
162  const boolList& doneEdge,
163  const label patchFacei,
164  const label globalFacei
165 ) const
166 {
167  const labelList& fEdges = pp.faceEdges()[patchFacei];
168 
169  label startFp = -1;
170  label endFp = -1;
171 
172  // Get edge that hasn't been done yet but needs extrusion
173  forAll(fEdges, fp)
174  {
175  label edgei = fEdges[fp];
176  const edge& e = pp.edges()[edgei];
177 
178  if
179  (
180  !doneEdge[edgei]
181  && ( addedPoints_[e[0]].size() || addedPoints_[e[1]].size() )
182  )
183  {
184  startFp = fp;
185  break;
186  }
187  }
188 
189  if (startFp != -1)
190  {
191  // We found an edge that needs extruding but hasn't been done yet.
192  // Now find the face on the other side
193  label nbrGlobalFacei = nbrFace
194  (
195  globalEdgeFaces,
196  fEdges[startFp],
197  globalFacei
198  );
199 
200  if (nbrGlobalFacei == -1)
201  {
202  // Proper boundary edge. Only extrude single edge.
203  endFp = startFp;
204  }
205  else
206  {
207  // Search back for edge
208  // - which hasn't been handled yet
209  // - with same neighbour
210  // - that needs extrusion
211 
212  const label initFp = startFp;
213  while (true)
214  {
215  label prevFp = fEdges.rcIndex(startFp);
216 
217  if (prevFp == initFp)
218  {
219  const edge& e = pp.edges()[fEdges[initFp]];
220  const face& localF = pp.localFaces()[patchFacei];
221 
223  << "On face:" << patchFacei
224  << " fc:" << pp.faceCentres()[patchFacei]
225  << " vertices:" << localF
226  << " points:"
227  << UIndirectList<point>(pp.points(), pp[patchFacei])
228  << " edges:" << fEdges
229  << " All edges of face seem to have same neighbour "
230  << nbrGlobalFacei
231  << " starting walking from edge " << e
232  << exit(FatalError);
233  }
234 
235  if
236  (
237  !sameEdgeNeighbour
238  (
239  pp,
240  globalEdgeFaces,
241  doneEdge,
242  globalFacei,
243  nbrGlobalFacei,
244  fEdges[prevFp]
245  )
246  )
247  {
248  break;
249  }
250  startFp = prevFp;
251  }
252 
253  // Search forward for end of string
254  endFp = startFp;
255  while (true)
256  {
257  label nextFp = fEdges.fcIndex(endFp);
258 
259  if
260  (
261  !sameEdgeNeighbour
262  (
263  pp,
264  globalEdgeFaces,
265  doneEdge,
266  globalFacei,
267  nbrGlobalFacei,
268  fEdges[nextFp]
269  )
270  )
271  {
272  break;
273  }
274  endFp = nextFp;
275  }
276  }
277  }
278 
279  return labelPair(startFp, endFp);
280 }
281 
282 
283 Foam::label Foam::addPatchCellLayer::addSideFace
284 (
285  const indirectPrimitivePatch& pp,
286  const labelListList& addedCells, // per pp face the new extruded cell
287  const face& newFace,
288  const label newPatchID,
289  const label zonei,
290  const bool edgeFlip,
291  const label inflateFacei,
292 
293  const label ownFacei, // pp face that provides owner
294  const label nbrFacei,
295  const label meshEdgei, // corresponding mesh edge
296  const label layeri, // layer
297  const label numEdgeFaces, // number of layers for edge
298  const labelList& meshFaces, // precalculated edgeFaces
299  polyTopoChange& meshMod
300 ) const
301 {
302  // Adds a side face i.e. extrudes a patch edge.
303 
304  label addedFacei = -1;
305 
306 
307  // Is patch edge external edge of indirectPrimitivePatch?
308  if (nbrFacei == -1)
309  {
310  // External edge so external face.
311 
312  // Determine if different number of layer on owner and neighbour side
313  // (relevant only for coupled faces). See section for internal edge
314  // below.
315 
316  label layerOwn;
317 
318  if (addedCells[ownFacei].size() < numEdgeFaces)
319  {
320  label offset = numEdgeFaces - addedCells[ownFacei].size();
321  if (layeri <= offset)
322  {
323  layerOwn = 0;
324  }
325  else
326  {
327  layerOwn = layeri - offset;
328  }
329  }
330  else
331  {
332  layerOwn = layeri;
333  }
334 
335 
336  //Pout<< "Added boundary face:" << newFace
337  // << " atfc:" << newFace.centre(meshMod.points())
338  // << " n:" << newFace.unitNormal(meshMod.points())
339  // << " own:" << addedCells[ownFacei][layerOwn]
340  // << " patch:" << newPatchID
341  // << endl;
342 
343  addedFacei = meshMod.setAction
344  (
345  polyAddFace
346  (
347  newFace, // face
348  addedCells[ownFacei][layerOwn], // owner
349  -1, // neighbour
350  -1, // master point
351  -1, // master edge
352  inflateFacei, // master face
353  false, // flux flip
354  newPatchID, // patch for face
355  zonei, // zone for face
356  edgeFlip // face zone flip
357  )
358  );
359  }
360  else
361  {
362  // When adding side faces we need to modify neighbour and owners
363  // in region where layer mesh is stopped. Determine which side
364  // has max number of faces and make sure layers match closest to
365  // original pp if there are different number of layers.
366 
367  label layerNbr;
368  label layerOwn;
369 
370  if (addedCells[ownFacei].size() > addedCells[nbrFacei].size())
371  {
372  label offset =
373  addedCells[ownFacei].size() - addedCells[nbrFacei].size();
374 
375  layerOwn = layeri;
376 
377  if (layeri <= offset)
378  {
379  layerNbr = 0;
380  }
381  else
382  {
383  layerNbr = layeri - offset;
384  }
385  }
386  else if (addedCells[nbrFacei].size() > addedCells[ownFacei].size())
387  {
388  label offset =
389  addedCells[nbrFacei].size() - addedCells[ownFacei].size();
390 
391  layerNbr = layeri;
392 
393  if (layeri <= offset)
394  {
395  layerOwn = 0;
396  }
397  else
398  {
399  layerOwn = layeri - offset;
400  }
401  }
402  else
403  {
404  // Same number of layers on both sides.
405  layerNbr = layeri;
406  layerOwn = layeri;
407  }
408 
409 
410  // Check mesh internal faces using edge to initialise
411  label inflateEdgei = -1;
412  if (addToMesh_)
413  {
414  forAll(meshFaces, i)
415  {
416  if (mesh_.isInternalFace(meshFaces[i]))
417  {
418  // meshEdge uses internal faces so ok to inflate from it
419  inflateEdgei = meshEdgei;
420  break;
421  }
422  }
423  }
424 
425 
426  addedFacei = meshMod.setAction
427  (
428  polyAddFace
429  (
430  newFace, // face
431  addedCells[ownFacei][layerOwn], // owner
432  addedCells[nbrFacei][layerNbr], // neighbour
433  -1, // master point
434  inflateEdgei, // master edge
435  -1, // master face
436  false, // flux flip
437  -1, // patch for face
438  zonei, // zone for face
439  edgeFlip // face zone flip
440  )
441  );
442 
443  //Pout<< "Added internal face:" << newFace
444  // << " atfc:" << newFace.centre(meshMod.points())
445  // << " n:" << newFace.unitNormal(meshMod.points())
446  // << " own:" << addedCells[ownFacei][layerOwn]
447  // << " nei:" << addedCells[nbrFacei][layerNbr]
448  // << endl;
449  }
450 
451  return addedFacei;
452 }
453 
454 
455 void Foam::addPatchCellLayer::setFaceProps
456 (
457  const polyMesh& mesh,
458  const label facei,
459 
460  label& patchi,
461  label& zonei,
462  bool& zoneFlip
463 )
464 {
465  patchi = mesh.boundaryMesh().whichPatch(facei);
466  zonei = mesh.faceZones().whichZone(facei);
467  if (zonei != -1)
468  {
469  label index = mesh.faceZones()[zonei].whichFace(facei);
470  zoneFlip = mesh.faceZones()[zonei].flipMap()[index];
471  }
472 }
473 
474 
475 void Foam::addPatchCellLayer::setFaceProps
476 (
477  const polyMesh& mesh,
478  const indirectPrimitivePatch& pp,
479  const label ppEdgeI,
480  const label faceI,
481 
482  label& patchI,
483  label& zoneI,
484  bool& zoneFlip,
485  label& inflateFaceI
486 )
487 {
488  setFaceProps
489  (
490  mesh,
491  faceI,
492 
493  patchI,
494  zoneI,
495  zoneFlip
496  );
497 
498  if (patchI != -1 || zoneI != -1)
499  {
500  inflateFaceI = faceI;
501  }
502 
503  if (zoneI != -1)
504  {
505  // Correct flip for patch edge ordering
506  const edge& ppEdge = pp.edges()[ppEdgeI];
507  const edge patchEdge
508  (
509  pp.meshPoints()[ppEdge[0]],
510  pp.meshPoints()[ppEdge[1]]
511  );
512 
513  const face& f = mesh.faces()[faceI];
514  bool found = false;
515  forAll(f, fp)
516  {
517  const edge e(f[fp], f.nextLabel(fp));
518  int stat = edge::compare(e, patchEdge);
519  if (stat == 1)
520  {
521  found = true;
522  break;
523  }
524  else if (stat == -1)
525  {
526  found = true;
527  zoneFlip = !zoneFlip;
528  break;
529  }
530  }
531 
532  if (!found)
533  {
534  //FatalErrorInFunction
536  << "Problem: cannot find patch edge " << ppEdgeI
537  << " with mesh vertices " << patchEdge
538  << " at " << patchEdge.line(mesh.points())
539  << " in face " << faceI << " with mesh vertices "
540  << f
541  << " at " << pointField(mesh.points(), f)
542  << endl
543  << "Continuing with potentially incorrect faceZone orientation"
544  //<< exit(FatalError);
545  << endl;
546  }
547  }
548 }
549 
550 
551 void Foam::addPatchCellLayer::findZoneFace
552 (
553  const bool useInternalFaces,
554  const bool useBoundaryFaces,
555 
556  const polyMesh& mesh,
557  const indirectPrimitivePatch& pp,
558  const label ppEdgeI,
559  const labelUIndList& excludeFaces,
560  const labelList& meshFaces,
561 
562  label& inflateFaceI,
563  label& patchI,
564  label& zoneI,
565  bool& zoneFlip
566 )
567 {
568  inflateFaceI = -1;
569  patchI = -1;
570  zoneI = -1;
571  zoneFlip = false;
572 
573  forAll(meshFaces, k)
574  {
575  label faceI = meshFaces[k];
576 
577  if
578  (
579  !excludeFaces.found(faceI)
580  && (
581  (mesh.isInternalFace(faceI) && useInternalFaces)
582  || (!mesh.isInternalFace(faceI) && useBoundaryFaces)
583  )
584  )
585  {
586  setFaceProps
587  (
588  mesh,
589  pp,
590  ppEdgeI,
591  faceI,
592 
593  patchI,
594  zoneI,
595  zoneFlip,
596  inflateFaceI
597  );
598 
599  if (zoneI != -1 || patchI != -1)
600  {
601  break;
602  }
603  }
604  }
605 }
606 
607 
608 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
609 
610 // Construct from mesh
611 Foam::addPatchCellLayer::addPatchCellLayer
612 (
613  const polyMesh& mesh,
614  const bool addToMesh
615 )
616 :
617  mesh_(mesh),
618  addToMesh_(addToMesh),
619  addedPoints_(0),
620  layerFaces_(0)
621 {}
622 
623 
624 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
625 
627 (
628  const polyMesh& mesh,
629  const labelListList& layerFaces
630 )
631 {
632  labelListList layerCells(layerFaces.size());
633 
634  forAll(layerFaces, patchFacei)
635  {
636  const labelList& faceLabels = layerFaces[patchFacei];
637 
638  if (faceLabels.size())
639  {
640  labelList& added = layerCells[patchFacei];
641  added.setSize(faceLabels.size()-1);
642 
643  for (label i = 0; i < faceLabels.size()-1; i++)
644  {
645  added[i] = mesh.faceNeighbour()[faceLabels[i]];
646  }
647  }
648  }
649  return layerCells;
650 }
651 
652 
654 {
655  return addedCells(mesh_, layerFaces_);
656 }
657 
658 
660 (
661  const polyMesh& mesh,
662  const globalIndex& globalFaces,
664 )
665 {
666  // Precalculate mesh edges for pp.edges.
667  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
668 
669  // From mesh edge to global face labels. Non-empty sublists only for
670  // pp edges.
671  labelListList globalEdgeFaces(mesh.nEdges());
672 
673  const labelListList& edgeFaces = pp.edgeFaces();
674 
675  forAll(edgeFaces, edgeI)
676  {
677  label meshEdgeI = meshEdges[edgeI];
678 
679  const labelList& eFaces = edgeFaces[edgeI];
680 
681  // Store face and processor as unique tag.
682  labelList& globalEFaces = globalEdgeFaces[meshEdgeI];
683  globalEFaces.setSize(eFaces.size());
684  forAll(eFaces, i)
685  {
686  globalEFaces[i] = globalFaces.toGlobal(pp.addressing()[eFaces[i]]);
687  }
688  }
689 
690  // Synchronise across coupled edges.
692  (
693  mesh,
694  globalEdgeFaces,
695  ListOps::uniqueEqOp<label>(),
696  labelList() // null value
697  );
698 
699  // Extract pp part
700  return labelListList(UIndirectList<labelList>(globalEdgeFaces, meshEdges));
701 }
702 
703 
704 void Foam::addPatchCellLayer::markPatchEdges
705 (
706  const polyMesh& mesh,
707  const indirectPrimitivePatch& pp,
708  const labelListList& edgeGlobalFaces,
709  const labelList& meshEdges,
710 
711  bitSet& isPatchEdge,
712  bitSet& isPatchBoundaryEdge
713 )
714 {
715  // Mark (mesh) edges:
716  // - anywhere on extrusion
717  // - where the extrusion ends
718 
719  isPatchEdge.setSize(mesh.nEdges());
720  isPatchEdge = false;
721  isPatchEdge.set(meshEdges);
722  // Synchronise across coupled edges
724  (
725  mesh,
726  isPatchEdge,
727  orEqOp<unsigned int>(),
728  false // initial value
729  );
730 
731  isPatchBoundaryEdge.setSize(mesh.nEdges());
732  isPatchBoundaryEdge = false;
733  forAll(edgeGlobalFaces, edgei)
734  {
735  // Test that edge has single global extruded face.
736  // Mark on processor that holds the face (since edgeGlobalFaces
737  // only gets filled from pp faces so if there is only one this
738  // is it)
739  if (edgeGlobalFaces[edgei].size() == 1)
740  {
741  isPatchBoundaryEdge.set(meshEdges[edgei]);
742  }
743  }
744  // Synchronise across coupled edges
746  (
747  mesh,
748  isPatchBoundaryEdge,
749  orEqOp<unsigned int>(),
750  false // initial value
751  );
752 }
753 
754 
755 void Foam::addPatchCellLayer::globalEdgeInfo
756 (
757  const bool zoneFromAnyFace,
758 
759  const polyMesh& mesh,
760  const globalIndex& globalFaces,
761  const labelListList& edgeGlobalFaces,
762  const indirectPrimitivePatch& pp,
763  const labelList& meshEdges,
764 
765  labelList& patchEdgeToFace, // face (in globalFaces index)
766  labelList& patchEdgeToPatch, // patch on face (or -1 for internal faces)
767  labelList& patchEdgeToZone, // zone on face
768  bitSet& patchEdgeToFlip // flip orientation on face
769 )
770 {
771  // For every edge on the outside of the patch return a potential patch/
772  // faceZone to extrude into.
773 
774  // Mark (mesh) edges on pp.
775  bitSet isExtrudeEdge;
776  bitSet isBoundaryEdge;
777  markPatchEdges
778  (
779  mesh,
780  pp,
781  edgeGlobalFaces,
782  meshEdges,
783 
784  isExtrudeEdge,
785  isBoundaryEdge
786  );
787 
788  // Build map of pp edges (in mesh point indexing). Note that this
789  // is now also on processors that do not have pp (but do have the edge)
790  EdgeMap<label> isBoundaryEdgeSet(pp.nEdges());
791  for (const label edgei : isBoundaryEdge)
792  {
793  isBoundaryEdgeSet.insert(mesh.edges()[edgei], edgei);
794  }
795  EdgeMap<label> isExtrudeEdgeSet(pp.nEdges());
796  for (const label edgei : isExtrudeEdge)
797  {
798  isExtrudeEdgeSet.insert(mesh.edges()[edgei], edgei);
799  }
800 
801 
802  const faceZoneMesh& fzs = mesh.faceZones();
803 
804  // Extract zone info into mesh face indexing for ease of addressing
805  labelList faceToZone(mesh.nFaces(), -1);
806  bitSet faceToFlip(mesh.nFaces());
807  for (const faceZone& fz: fzs)
808  {
809  const labelList& addressing = fz;
810  UIndirectList<label>(faceToZone, addressing) = fz.index();
811 
812  const boolList& fm = fz.flipMap();
813  forAll(addressing, i)
814  {
815  faceToFlip[addressing[i]] = fm[i];
816  }
817  }
818 
819 
820  // Storage (over all mesh edges)
821  // - face that data originates from (in globalFaces indexing)
822  labelList meshEdgeToFace(mesh.nEdges(), -1);
823  // - patch (for boundary faces)
824  labelList meshEdgeToPatch(mesh.nEdges(), -1);
825  // - faceZone
826  labelList meshEdgeToZone(mesh.nEdges(), -1);
827  // - faceZone orientation
828  bitSet meshEdgeToFlip(mesh.nEdges());
829 
830  //if (useInternalFaces)
831  {
832  const bitSet isInternalOrCoupled
833  (
835  );
836 
837  // Loop over edges of the face to find any faceZone.
838  // Edges kept as point pair so we don't construct mesh.faceEdges etc.
839 
840  for (const label facei : isInternalOrCoupled)
841  {
842  const face& f = mesh.faces()[facei];
843 
844  label prevPointi = f.last();
845  for (const label pointi : f)
846  {
847  const edge e(prevPointi, pointi);
848 
849  // Check if edge is internal to extrusion. Take over faceZone
850  // etc from internal face.
851  const auto eFnd = isExtrudeEdgeSet.cfind(e);
852  if (eFnd)
853  {
854  const label edgei = eFnd();
855 
856  if (faceToZone[facei] != -1)
857  {
858  // Found a zoned internal face. Use.
859  meshEdgeToFace[edgei] = globalFaces.toGlobal(facei);
860  meshEdgeToZone[edgei] = faceToZone[facei];
861  const edge& meshE = mesh.edges()[edgei];
862  const int d = edge::compare(e, meshE);
863  if (d == 1)
864  {
865  meshEdgeToFlip[edgei] = faceToFlip[facei];
866  }
867  else if (d == -1)
868  {
869  meshEdgeToFlip[edgei] = !faceToFlip[facei];
870  }
871  else
872  {
873  FatalErrorInFunction << "big problem"
874  << exit(FatalError);
875  }
876  }
877  }
878 
879  prevPointi = pointi;
880  }
881  }
882  }
883 
884 
885  //if (useBoundaryFaces)
886  {
887  // Loop over all patches and find 'best' one (non-coupled,
888  // non-extrusion, non-constraint(?)). Note that logic is slightly
889  // different from internal faces loop above - first patch face
890  // is being used instead of first zoned face for internal faces
891 
892  const polyBoundaryMesh& patches = mesh.boundaryMesh();
893 
894  bitSet isPpFace(mesh.nFaces());
895  isPpFace.set(pp.addressing());
896  // Note: no need to sync ppFace since does not include processor patches
897 
898  for (const polyPatch& pp : patches)
899  {
900  if (!pp.coupled())
901  {
902  // TBD. Check for constraint? This is usually a geometric
903  // constraint and not a topological one so should
904  // be handled in the extrusion vector calculation instead?
905 
906  forAll(pp, i)
907  {
908  const label facei = pp.start()+i;
909 
910  if (!isPpFace[facei])
911  {
912  const face& f = pp[i];
913 
914  label prevPointi = f.last();
915  for (const label pointi : f)
916  {
917  const edge e(prevPointi, pointi);
918  const auto eFnd =
919  (
920  zoneFromAnyFace
921  ? isExtrudeEdgeSet.cfind(e)
922  : isBoundaryEdgeSet.cfind(e)
923  );
924  if (eFnd)
925  {
926  const label edgei = eFnd();
927  if (meshEdgeToFace[edgei] == -1)
928  {
929  // Found unassigned face. Use its
930  // information.
931  // Note that we use the lowest numbered
932  // patch face.
933 
934  meshEdgeToFace[edgei] =
935  globalFaces.toGlobal(facei);
936  }
937 
938  // Override any patch info. Note that
939  // meshEdgeToFace might be an internal face.
940  if (meshEdgeToPatch[edgei] == -1)
941  {
942  meshEdgeToPatch[edgei] = pp.index();
943  }
944 
945  // Override any zone info
946  if (meshEdgeToZone[edgei] == -1)
947  {
948  meshEdgeToZone[edgei] =
949  faceToZone[facei];
950  const edge& meshE = mesh.edges()[edgei];
951  const int d = edge::compare(e, meshE);
952  if (d == 1)
953  {
954  meshEdgeToFlip[edgei] =
955  faceToFlip[facei];
956  }
957  else if (d == -1)
958  {
959  meshEdgeToFlip[edgei] =
960  !faceToFlip[facei];
961  }
962  else
963  {
965  << "big problem"
966  << exit(FatalError);
967  }
968  }
969  }
970 
971  prevPointi = pointi;
972  }
973  }
974  }
975  }
976  }
977  }
978 
979 
980  // Synchronise across coupled edges. Max patch/face/faceZone wins
982  (
983  mesh,
984  meshEdgeToFace,
985  maxEqOp<label>(),
986  label(-1)
987  );
989  (
990  mesh,
991  meshEdgeToPatch,
992  maxEqOp<label>(),
993  label(-1)
994  );
996  (
997  mesh,
998  meshEdgeToZone,
999  maxEqOp<label>(),
1000  label(-1)
1001  );
1002 // // Note: flipMap not yet done. Needs edge orientation. This is handled
1003 // // later on.
1004 // if (Pstream::parRun())
1005 // {
1006 // const globalMeshData& gd = mesh.globalData();
1007 // const indirectPrimitivePatch& cpp = gd.coupledPatch();
1008 //
1009 // labelList patchEdges;
1010 // labelList coupledEdges;
1011 // bitSet sameEdgeOrientation;
1012 // PatchTools::matchEdges
1013 // (
1014 // pp,
1015 // cpp,
1016 // patchEdges,
1017 // coupledEdges,
1018 // sameEdgeOrientation
1019 // );
1020 //
1021 // // Convert data on pp edges to data on coupled patch
1022 // labelList cppEdgeZoneID(cpp.nEdges(), -1);
1023 // boolList cppEdgeFlip(cpp.nEdges(), false);
1024 // forAll(coupledEdges, i)
1025 // {
1026 // label cppEdgei = coupledEdges[i];
1027 // label ppEdgei = patchEdges[i];
1028 //
1029 // cppEdgeZoneID[cppEdgei] = edgeZoneID[ppEdgei];
1030 // if (sameEdgeOrientation[i])
1031 // {
1032 // cppEdgeFlip[cppEdgei] = edgeFlip[ppEdgei];
1033 // }
1034 // else
1035 // {
1036 // cppEdgeFlip[cppEdgei] = !edgeFlip[ppEdgei];
1037 // }
1038 // }
1039 //
1040 // // Sync
1041 // const globalIndexAndTransform& git = gd.globalTransforms();
1042 // const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
1043 //
1044 // globalMeshData::syncData
1045 // (
1046 // cppEdgeZoneID,
1047 // gd.globalEdgeSlaves(),
1048 // gd.globalEdgeTransformedSlaves(),
1049 // edgeMap,
1050 // git,
1051 // maxEqOp<label>(),
1052 // dummyTransform()
1053 // );
1054 // globalMeshData::syncData
1055 // (
1056 // cppEdgeFlip,
1057 // gd.globalEdgeSlaves(),
1058 // gd.globalEdgeTransformedSlaves(),
1059 // edgeMap,
1060 // git,
1061 // andEqOp<bool>(),
1062 // dummyTransform()
1063 // );
1064 //
1065 // // Convert data on coupled edges to pp edges
1066 // forAll(coupledEdges, i)
1067 // {
1068 // label cppEdgei = coupledEdges[i];
1069 // label ppEdgei = patchEdges[i];
1070 //
1071 // edgeZoneID[ppEdgei] = cppEdgeZoneID[cppEdgei];
1072 // if (sameEdgeOrientation[i])
1073 // {
1074 // edgeFlip[ppEdgei] = cppEdgeFlip[cppEdgei];
1075 // }
1076 // else
1077 // {
1078 // edgeFlip[ppEdgei] = !cppEdgeFlip[cppEdgei];
1079 // }
1080 // }
1081 // }
1082 
1084  (
1085  mesh,
1086  meshEdgeToFlip,
1087  orEqOp<unsigned int>(),
1088  0
1089  );
1090 
1091  // Extract pp info
1092  patchEdgeToFace = UIndirectList<label>(meshEdgeToFace, meshEdges);
1093  patchEdgeToPatch = UIndirectList<label>(meshEdgeToPatch, meshEdges);
1094  patchEdgeToZone = UIndirectList<label>(meshEdgeToZone, meshEdges);
1095  patchEdgeToFlip.setSize(meshEdges.size());
1096  patchEdgeToFlip = false;
1097  forAll(meshEdges, i)
1098  {
1099  patchEdgeToFlip[i] = meshEdgeToFlip[meshEdges[i]];
1100  }
1101 }
1102 
1103 
1105 (
1106  const bool zoneFromAnyFace,
1107 
1108  const polyMesh& mesh,
1109  const globalIndex& globalFaces,
1110  const labelListList& globalEdgeFaces,
1111  const indirectPrimitivePatch& pp,
1112 
1113  labelList& edgePatchID,
1114  label& nPatches,
1115  Map<label>& nbrProcToPatch,
1116  Map<label>& patchToNbrProc,
1117  labelList& edgeZoneID,
1118  boolList& edgeFlip,
1119  labelList& inflateFaceID
1120 )
1121 {
1122  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1123  const globalMeshData& gd = mesh.globalData();
1124 
1125  // Precalculate mesh edges for pp.edges.
1126  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
1127 
1128  edgePatchID.setSize(pp.nEdges());
1129  edgePatchID = -1;
1130  nPatches = patches.size();
1131  edgeZoneID.setSize(pp.nEdges());
1132  edgeZoneID = -1;
1133  edgeFlip.setSize(pp.nEdges());
1134  edgeFlip = false;
1135  inflateFaceID.setSize(pp.nEdges(), -1);
1136 
1137 
1138  // Determine properties for faces extruded from edges
1139  // - edge inbetween two different processors:
1140  // - extrude as patch face on correct processor
1141  // - edge at end of patch (so edgeFaces.size() == 1):
1142  // - use mesh face that is a boundary face
1143  // - edge internal to patch (so edgeFaces.size() == 2):
1144 
1145 
1146  // Pass1:
1147  // For all edges inbetween two processors: see if matches to existing
1148  // processor patch or create interprocessor-patch if necessary.
1149  // Sets edgePatchID[edgeI] but none of the other quantities
1150 
1151  forAll(globalEdgeFaces, edgei)
1152  {
1153  const labelList& eGlobalFaces = globalEdgeFaces[edgei];
1154  if
1155  (
1156  eGlobalFaces.size() == 2
1157  && pp.edgeFaces()[edgei].size() == 1
1158  )
1159  {
1160  // Locally but not globally a boundary edge. Hence a coupled
1161  // edge. Find the patch to use if on different processors.
1162 
1163  label f0 = eGlobalFaces[0];
1164  label f1 = eGlobalFaces[1];
1165 
1166  label otherProci = -1;
1167  if (globalFaces.isLocal(f0) && !globalFaces.isLocal(f1))
1168  {
1169  otherProci = globalFaces.whichProcID(f1);
1170  }
1171  else if (!globalFaces.isLocal(f0) && globalFaces.isLocal(f1))
1172  {
1173  otherProci = globalFaces.whichProcID(f0);
1174  }
1175 
1176 
1177  if (otherProci != -1)
1178  {
1179  // Use existing processorPolyPatch to otherProci?
1180 
1181  label procPatchi =
1182  gd.topology().procPatchLookup(otherProci);
1183 
1184  if (procPatchi < 0)
1185  {
1186  // No existing processorPolyPatch to otherProci.
1187  // See if already marked for addition
1188  procPatchi = nbrProcToPatch.lookup(otherProci, -1);
1189 
1190  if (procPatchi < 0)
1191  {
1192  // Add new proc-patch, mark for addition.
1193 
1194  procPatchi = nPatches;
1195  ++nPatches;
1196 
1197  nbrProcToPatch.insert(otherProci, procPatchi);
1198  patchToNbrProc.insert(procPatchi, otherProci);
1199  }
1200  }
1201 
1202  edgePatchID[edgei] = procPatchi;
1203  }
1204  }
1205  }
1206 
1207 
1208  // Pass2: determine face properties for all other edges
1209  // ----------------------------------------------------
1210 
1211  // Info for edges of pp
1212  labelList edgeToFace;
1213  labelList edgeToPatch;
1214  labelList edgeToZone;
1215  bitSet edgeToFlip;
1216  globalEdgeInfo
1217  (
1218  zoneFromAnyFace, // internal edge info also from boundary faces
1219 
1220  mesh,
1221  globalFaces,
1222  globalEdgeFaces,
1223  pp,
1224  meshEdges,
1225 
1226  edgeToFace, // face (in globalFaces index)
1227  edgeToPatch, // patch on face (or -1 for internal faces)
1228  edgeToZone, // zone on face
1229  edgeToFlip // flip orientation on face
1230  );
1231 
1232  const labelListList& edgeFaces = pp.edgeFaces();
1233 
1234  DynamicList<label> dynMeshEdgeFaces;
1235 
1236  forAll(edgeFaces, edgei)
1237  {
1238  if (edgePatchID[edgei] == -1)
1239  {
1240  if (edgeFaces[edgei].size() == 2)
1241  {
1242  // Internal edge. Look at any face (internal or boundary) to
1243  // determine extrusion properties. First one that has zone
1244  // info wins
1245  if (globalFaces.isLocal(edgeToFace[edgei]))
1246  {
1247  inflateFaceID[edgei] =
1248  globalFaces.toLocal(edgeToFace[edgei]);
1249  }
1250  edgeZoneID[edgei] = edgeToZone[edgei];
1251  edgeFlip[edgei] = edgeToFlip[edgei];
1252  }
1253  else
1254  {
1255  // Proper, uncoupled patch edge. Boundary to get info from
1256  // might be on a different processor!
1257 
1258  if (globalFaces.isLocal(edgeToFace[edgei]))
1259  {
1260  inflateFaceID[edgei] =
1261  globalFaces.toLocal(edgeToFace[edgei]);
1262  }
1263  edgePatchID[edgei] = edgeToPatch[edgei];
1264  edgeZoneID[edgei] = edgeToZone[edgei];
1265  edgeFlip[edgei] = edgeToFlip[edgei];
1266  }
1267  }
1268  }
1269 
1270 
1271 
1272  // Now hopefully every boundary edge has a edge patch. Check
1273  if (debug)
1274  {
1275  forAll(edgeFaces, edgei)
1276  {
1277  if (edgeFaces[edgei].size() == 1 && edgePatchID[edgei] == -1)
1278  {
1279  const edge& e = pp.edges()[edgei];
1281  << "Have no sidePatchID for edge " << edgei << " points "
1282  << pp.points()[pp.meshPoints()[e[0]]]
1283  << pp.points()[pp.meshPoints()[e[1]]]
1284  << endl;
1285  }
1286  }
1287  }
1288 
1289 
1290 
1291  // Pass3: for any faces set in pass1 see if we can find a processor face
1292  // to inherit from (we only have a patch, not a patch face)
1293  forAll(edgeFaces, edgei)
1294  {
1295  if
1296  (
1297  edgeFaces[edgei].size() == 1
1298  && globalEdgeFaces[edgei].size() == 2
1299  && edgePatchID[edgei] != -1
1300  && inflateFaceID[edgei] == -1
1301  )
1302  {
1303  // 1. Do we have a local boundary face to inflate from
1304 
1305  label myFaceI = pp.addressing()[edgeFaces[edgei][0]];
1306 
1307  // Pick up any boundary face on this edge and use its properties
1308  label meshEdgei = meshEdges[edgei];
1309  const labelList& meshFaces = mesh.edgeFaces
1310  (
1311  meshEdgei,
1312  dynMeshEdgeFaces
1313  );
1314 
1315  forAll(meshFaces, k)
1316  {
1317  label facei = meshFaces[k];
1318 
1319  if (facei != myFaceI && !mesh.isInternalFace(facei))
1320  {
1321  if (patches.whichPatch(facei) == edgePatchID[edgei])
1322  {
1323  setFaceProps
1324  (
1325  mesh,
1326  pp,
1327  edgei,
1328  facei,
1329 
1330  edgePatchID[edgei],
1331  edgeZoneID[edgei],
1332  edgeFlip[edgei],
1333  inflateFaceID[edgei]
1334  );
1335  break;
1336  }
1337  }
1338  }
1339  }
1340  }
1341 
1342 
1343 
1344  // Sync all data:
1345  // - edgePatchID : might be local in case of processor patch. Do not
1346  // sync for now
1347  // - inflateFaceID: local. Do not sync
1348  // - edgeZoneID : global. sync.
1349  // - edgeFlip : global. sync.
1350 
1351  if (Pstream::parRun())
1352  {
1353  const globalMeshData& gd = mesh.globalData();
1354  const indirectPrimitivePatch& cpp = gd.coupledPatch();
1355 
1356  labelList patchEdges;
1357  labelList coupledEdges;
1358  bitSet sameEdgeOrientation;
1360  (
1361  pp,
1362  cpp,
1363  patchEdges,
1364  coupledEdges,
1365  sameEdgeOrientation
1366  );
1367 
1368  // Convert data on pp edges to data on coupled patch
1369  labelList cppEdgeZoneID(cpp.nEdges(), -1);
1370  boolList cppEdgeFlip(cpp.nEdges(), false);
1371  forAll(coupledEdges, i)
1372  {
1373  label cppEdgei = coupledEdges[i];
1374  label ppEdgei = patchEdges[i];
1375 
1376  cppEdgeZoneID[cppEdgei] = edgeZoneID[ppEdgei];
1377  if (sameEdgeOrientation[i])
1378  {
1379  cppEdgeFlip[cppEdgei] = edgeFlip[ppEdgei];
1380  }
1381  else
1382  {
1383  cppEdgeFlip[cppEdgei] = !edgeFlip[ppEdgei];
1384  }
1385  }
1386 
1387  // Sync
1388  const globalIndexAndTransform& git = gd.globalTransforms();
1389  const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
1390 
1392  (
1393  cppEdgeZoneID,
1394  gd.globalEdgeSlaves(),
1395  gd.globalEdgeTransformedSlaves(),
1396  edgeMap,
1397  git,
1398  maxEqOp<label>(),
1399  dummyTransform()
1400  );
1402  (
1403  cppEdgeFlip,
1404  gd.globalEdgeSlaves(),
1405  gd.globalEdgeTransformedSlaves(),
1406  edgeMap,
1407  git,
1408  andEqOp<bool>(),
1409  dummyTransform()
1410  );
1411 
1412  // Convert data on coupled edges to pp edges
1413  forAll(coupledEdges, i)
1414  {
1415  label cppEdgei = coupledEdges[i];
1416  label ppEdgei = patchEdges[i];
1417 
1418  edgeZoneID[ppEdgei] = cppEdgeZoneID[cppEdgei];
1419  if (sameEdgeOrientation[i])
1420  {
1421  edgeFlip[ppEdgei] = cppEdgeFlip[cppEdgei];
1422  }
1423  else
1424  {
1425  edgeFlip[ppEdgei] = !cppEdgeFlip[cppEdgei];
1426  }
1427  }
1428  }
1429 }
1430 
1431 
1433 (
1434  const globalIndex& globalFaces,
1435  const labelListList& globalEdgeFaces,
1436  const scalarField& expansionRatio,
1437  const indirectPrimitivePatch& pp,
1438  const bitSet& ppFlip,
1439 
1440  const labelList& edgePatchID,
1441  const labelList& edgeZoneID,
1442  const boolList& edgeFlip,
1443  const labelList& inflateFaceID,
1444 
1445  const labelList& exposedPatchID,
1446  const labelList& nFaceLayers,
1447  const labelList& nPointLayers,
1448  const vectorField& firstLayerDisp,
1449  polyTopoChange& meshMod
1450 )
1451 {
1452  if (debug)
1453  {
1454  Pout<< "addPatchCellLayer::setRefinement : Adding up to "
1455  << gMax(nPointLayers)
1456  << " layers of cells to indirectPrimitivePatch with "
1457  << pp.nPoints() << " points" << endl;
1458  }
1459 
1460  if
1461  (
1462  pp.nPoints() != firstLayerDisp.size()
1463  || pp.nPoints() != nPointLayers.size()
1464  || pp.size() != nFaceLayers.size()
1465  || pp.size() != ppFlip.size()
1466  )
1467  {
1469  << "Size of new points is not same as number of points used by"
1470  << " the face subset" << endl
1471  << " patch.nPoints:" << pp.nPoints()
1472  << " displacement:" << firstLayerDisp.size()
1473  << " nPointLayers:" << nPointLayers.size() << nl
1474  << " patch.nFaces:" << pp.size()
1475  << " flip map:" << ppFlip.size()
1476  << " nFaceLayers:" << nFaceLayers.size()
1477  << abort(FatalError);
1478  }
1479  if (!addToMesh_)
1480  {
1481  // flip map should be false
1482  if (ppFlip.count())
1483  {
1485  << "In generating stand-alone mesh the flip map should be empty"
1486  << ". Instead it is " << ppFlip.count()
1487  << abort(FatalError);
1488  }
1489  }
1490  else
1491  {
1492  // Maybe check for adding to neighbour of boundary faces? How about
1493  // coupled faces where the faceZone flipMap is negated
1494 
1495  // For all boundary faces:
1496  // -1 : not extruded
1497  // 0 : extruded from owner outwards (flip = false)
1498  // 1 : extrude from neighbour outwards
1499  labelList stateAndFlip(mesh_.nBoundaryFaces(), 0);
1500  forAll(pp.addressing(), patchFacei)
1501  {
1502  if (nFaceLayers[patchFacei] > 0)
1503  {
1504  const label meshFacei = pp.addressing()[patchFacei];
1505  const label bFacei = meshFacei-mesh_.nInternalFaces();
1506  if (bFacei >= 0)
1507  {
1508  stateAndFlip[bFacei] = label(ppFlip[patchFacei]);
1509  }
1510  }
1511  }
1512  // Make sure uncoupled patches do not trigger the error below
1513  for (const auto& patch : mesh_.boundaryMesh())
1514  {
1515  if (!patch.coupled())
1516  {
1517  forAll(patch, i)
1518  {
1519  label& state = stateAndFlip[patch.offset()+i];
1520  state = (state == 0 ? 1 : 0);
1521  }
1522  }
1523  }
1524  syncTools::swapBoundaryFaceList(mesh_, stateAndFlip);
1525 
1526  forAll(pp.addressing(), patchFacei)
1527  {
1528  if (nFaceLayers[patchFacei] > 0)
1529  {
1530  const label meshFacei = pp.addressing()[patchFacei];
1531  const label bFacei = meshFacei-mesh_.nInternalFaces();
1532  if (bFacei >= 0)
1533  {
1534  if (stateAndFlip[bFacei] == -1)
1535  {
1537  << "At extruded face:" << meshFacei
1538  << " at:" << mesh_.faceCentres()[meshFacei]
1539  << " locally have nLayers:"
1540  << nFaceLayers[patchFacei]
1541  << " but remotely have none" << exit(FatalError);
1542  }
1543  else if (stateAndFlip[bFacei] == label(ppFlip[patchFacei]))
1544  {
1546  << "At extruded face:" << meshFacei
1547  << " at:" << mesh_.faceCentres()[meshFacei]
1548  << " locally have flip:" << ppFlip[patchFacei]
1549  << " which is not the opposite of coupled version "
1550  << stateAndFlip[bFacei]
1551  << exit(FatalError);
1552  }
1553  }
1554  }
1555  }
1556  }
1557 
1558 
1559  forAll(nPointLayers, i)
1560  {
1561  if (nPointLayers[i] < 0)
1562  {
1564  << "Illegal number of layers " << nPointLayers[i]
1565  << " at patch point " << i << abort(FatalError);
1566  }
1567  }
1568  forAll(nFaceLayers, i)
1569  {
1570  if (nFaceLayers[i] < 0)
1571  {
1573  << "Illegal number of layers " << nFaceLayers[i]
1574  << " at patch face " << i << abort(FatalError);
1575  }
1576  }
1577 
1578  forAll(globalEdgeFaces, edgei)
1579  {
1580  if (globalEdgeFaces[edgei].size() > 2)
1581  {
1582  const edge& e = pp.edges()[edgei];
1583 
1584  if (nPointLayers[e[0]] > 0 || nPointLayers[e[1]] > 0)
1585  {
1587  << "Trying to extrude edge "
1588  << e.line(pp.localPoints())
1589  << " which is non-manifold (has "
1590  << globalEdgeFaces[edgei].size()
1591  << " local or coupled faces using it)"
1592  << " of which "
1593  << pp.edgeFaces()[edgei].size()
1594  << " local"
1595  << abort(FatalError);
1596  }
1597  }
1598  }
1599 
1600 
1601  const labelList& meshPoints = pp.meshPoints();
1602 
1603 
1604  // Determine which points are on which side of the extrusion
1605  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1606 
1607  bitSet isBlockedFace(mesh_.nFaces());
1608  forAll(nFaceLayers, patchFacei)
1609  {
1610  if (nFaceLayers[patchFacei] > 0)
1611  {
1612  isBlockedFace.set(pp.addressing()[patchFacei]);
1613  }
1614  }
1615 
1616  // Some storage for edge-face-addressing.
1617  DynamicList<label> ef;
1618 
1619  // Precalculate mesh edges for pp.edges.
1620  const labelList meshEdges(pp.meshEdges(mesh_.edges(), mesh_.pointEdges()));
1621 
1622  if (debug)
1623  {
1624  // Check synchronisation
1625  // ~~~~~~~~~~~~~~~~~~~~~
1626 
1627  {
1628  labelList n(mesh_.nPoints(), Zero);
1629  labelUIndList(n, meshPoints) = nPointLayers;
1630  syncTools::syncPointList(mesh_, n, maxEqOp<label>(), label(0));
1631 
1632  // Non-synced
1633  forAll(meshPoints, i)
1634  {
1635  label meshPointi = meshPoints[i];
1636 
1637  if (n[meshPointi] != nPointLayers[i])
1638  {
1640  << "At mesh point:" << meshPointi
1641  << " coordinate:" << mesh_.points()[meshPointi]
1642  << " specified nLayers:" << nPointLayers[i] << endl
1643  << "On coupled point a different nLayers:"
1644  << n[meshPointi] << " was specified."
1645  << abort(FatalError);
1646  }
1647  }
1648 
1649 
1650  // Check that nPointLayers equals the max layers of connected faces
1651  // (or 0). Anything else makes no sense.
1652  labelList nFromFace(mesh_.nPoints(), Zero);
1653  forAll(nFaceLayers, i)
1654  {
1655  const face& f = pp[i];
1656 
1657  forAll(f, fp)
1658  {
1659  label pointi = f[fp];
1660 
1661  nFromFace[pointi] = max(nFromFace[pointi], nFaceLayers[i]);
1662  }
1663  }
1665  (
1666  mesh_,
1667  nFromFace,
1668  maxEqOp<label>(),
1669  label(0)
1670  );
1671 
1672  forAll(nPointLayers, i)
1673  {
1674  label meshPointi = meshPoints[i];
1675 
1676  if
1677  (
1678  nPointLayers[i] > 0
1679  && nPointLayers[i] != nFromFace[meshPointi]
1680  )
1681  {
1683  << "At mesh point:" << meshPointi
1684  << " coordinate:" << mesh_.points()[meshPointi]
1685  << " specified nLayers:" << nPointLayers[i] << endl
1686  << "but the max nLayers of surrounding faces is:"
1687  << nFromFace[meshPointi]
1688  << abort(FatalError);
1689  }
1690  }
1691  }
1692 
1693  {
1694  pointField d(mesh_.nPoints(), vector::max);
1695  UIndirectList<point>(d, meshPoints) = firstLayerDisp;
1697  (
1698  mesh_,
1699  d,
1700  minEqOp<vector>(),
1701  vector::max
1702  );
1703 
1704  forAll(meshPoints, i)
1705  {
1706  label meshPointi = meshPoints[i];
1707 
1708  if (mag(d[meshPointi] - firstLayerDisp[i]) > SMALL)
1709  {
1711  << "At mesh point:" << meshPointi
1712  << " coordinate:" << mesh_.points()[meshPointi]
1713  << " specified displacement:" << firstLayerDisp[i]
1714  << endl
1715  << "On coupled point a different displacement:"
1716  << d[meshPointi] << " was specified."
1717  << abort(FatalError);
1718  }
1719  }
1720  }
1721 
1722  // Check that edges of pp (so ones that become boundary faces)
1723  // connect to only one boundary face. Guarantees uniqueness of
1724  // patch that they go into so if this is a coupled patch both
1725  // sides decide the same.
1726  // ~~~~~~~~~~~~~~~~~~~~~~
1727 
1728  for (label edgei = pp.nInternalEdges(); edgei < pp.nEdges(); edgei++)
1729  {
1730  const edge& e = pp.edges()[edgei];
1731 
1732  if (nPointLayers[e[0]] > 0 || nPointLayers[e[1]] > 0)
1733  {
1734  // Edge is to become a face
1735 
1736  const labelList& eFaces = pp.edgeFaces()[edgei];
1737 
1738  // First check: pp should be single connected.
1739  if (eFaces.size() != 1)
1740  {
1742  << "boundary-edge-to-be-extruded:"
1743  << pp.points()[meshPoints[e[0]]]
1744  << pp.points()[meshPoints[e[1]]]
1745  << " has more than two faces using it:" << eFaces
1746  << abort(FatalError);
1747  }
1748 
1749  //label myFacei = pp.addressing()[eFaces[0]];
1750  //
1751  //label meshEdgei = meshEdges[edgei];
1752  //
1754  //const labelList& meshFaces = mesh_.edgeFaces(meshEdgei, ef);
1755  //
1757  //const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1758  //
1759  //label bFacei = -1;
1760  //
1761  //forAll(meshFaces, i)
1762  //{
1763  // label facei = meshFaces[i];
1764  //
1765  // if (facei != myFacei)
1766  // {
1767  // if (!mesh_.isInternalFace(facei))
1768  // {
1769  // if (bFacei == -1)
1770  // {
1771  // bFacei = facei;
1772  // }
1773  // else
1774  // {
1775  // //FatalErrorInFunction
1776  // WarningInFunction
1777  // << "boundary-edge-to-be-extruded:"
1778  // << pp.points()[meshPoints[e[0]]]
1779  // << pp.points()[meshPoints[e[1]]]
1780  // << " has more than one boundary faces"
1781  // << " using it:"
1782  // << bFacei << " fc:"
1783  // << mesh_.faceCentres()[bFacei]
1784  // << " patch:" << patches.whichPatch(bFacei)
1785  // << " and " << facei << " fc:"
1786  // << mesh_.faceCentres()[facei]
1787  // << " patch:" << patches.whichPatch(facei)
1788  // << endl;
1789  // //abort(FatalError);
1790  // }
1791  // }
1792  // }
1793  //}
1794  }
1795  }
1796  }
1797 
1798 
1799  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1800 
1801  // Precalculated patchID for each patch face
1802  labelList patchID(pp.size());
1803 
1804  forAll(pp, patchFacei)
1805  {
1806  label meshFacei = pp.addressing()[patchFacei];
1807 
1808  patchID[patchFacei] = patches.whichPatch(meshFacei);
1809  }
1810 
1811 
1812  // From master point (in patch point label) to added points (in mesh point
1813  // label)
1814  addedPoints_.setSize(pp.nPoints());
1815 
1816  // Mark points that do not get extruded by setting size of addedPoints_ to 0
1817  label nTruncated = 0;
1818 
1819  forAll(nPointLayers, patchPointi)
1820  {
1821  if (nPointLayers[patchPointi] > 0)
1822  {
1823  addedPoints_[patchPointi].setSize(nPointLayers[patchPointi]);
1824  }
1825  else
1826  {
1827  nTruncated++;
1828  }
1829  }
1830 
1831  if (debug)
1832  {
1833  Pout<< "Not adding points at " << nTruncated << " out of "
1834  << pp.nPoints() << " points" << endl;
1835  }
1836 
1837 
1838  // Store per face whether it uses the duplicated point or the original one
1839  // Also mark any affected cells. We could transport the duplicated point
1840  // itself but since it is a processor-local index only we only transport
1841  // a boolean.
1842 
1843  // Per face, per index in face either labelMax or a valid index. Note:
1844  // most faces are not affected in which case the face will be zero size
1845  // and only have a nullptr and a size.
1846  faceList baseFaces(mesh_.nFaces());
1847  bitSet isAffectedCell(mesh_.nCells());
1848  {
1849  const faceList& localFaces = pp.localFaces();
1850  forAll(localFaces, patchFacei)
1851  {
1852  const face& f = localFaces[patchFacei];
1853  forAll(f, fp)
1854  {
1855  const label patchPointi = f[fp];
1856  if (nPointLayers[patchPointi] > 0)
1857  {
1858  const label meshFacei = pp.addressing()[patchFacei];
1859  face& baseF = baseFaces[meshFacei];
1860  // Initialise to labelMax if not yet sized
1861  baseF.setSize(f.size(), labelMax);
1862  baseF[fp] = pp.meshPoints()[patchPointi];
1863 
1864  if (ppFlip[patchFacei])
1865  {
1866  // Neighbour stays. Affected points on the owner side.
1867  const label celli = mesh_.faceOwner()[meshFacei];
1868  isAffectedCell.set(celli);
1869  }
1870  else if (mesh_.isInternalFace(meshFacei))
1871  {
1872  // Owner unaffected. Unaffected points on neighbour side
1873  const label celli = mesh_.faceNeighbour()[meshFacei];
1874  isAffectedCell.set(celli);
1875  }
1876  }
1877  }
1878  }
1879  }
1880 
1881  // Transport affected side across faces. Could do across edges: say we have
1882  // a loose cell edge-(but not face-)connected to face-to-be-extruded do
1883  // we want it to move with the extrusion or stay connected to the original?
1884  // For now just keep it connected to the original.
1885  {
1886  // Work space
1887  Map<label> minPointValue;
1888  faceList oldBoundaryFaces(mesh_.nBoundaryFaces());
1889 
1890  while (true)
1891  {
1892  bitSet newIsAffectedCell(mesh_.nCells());
1893 
1894  label nChanged = 0;
1895  for (const label celli : isAffectedCell)
1896  {
1897  const cell& cFaces = mesh_.cells()[celli];
1898 
1899  // 1. Determine marked base points. Inside a single cell all
1900  // faces use the same 'instance' of a point.
1901  minPointValue.clear();
1902  for (const label facei : cFaces)
1903  {
1904  const face& baseF = baseFaces[facei];
1905  const face& f = mesh_.faces()[facei];
1906 
1907  if (baseF.size())
1908  {
1909  forAll(f, fp)
1910  {
1911  if (baseF[fp] != labelMax)
1912  {
1913  // Could check here for inconsistent patchPoint
1914  // e.g. cell using both sides of a
1915  // face-to-be-extruded. Is not possible!
1916  minPointValue.insert(f[fp], baseF[fp]);
1917  }
1918  }
1919  }
1920  }
1921 
1922  //Pout<< "For cell:" << celli
1923  // << " at:" << mesh_.cellCentres()[celli]
1924  // << " have minPointValue:" << minPointValue
1925  // << endl;
1926 
1927  // 2. Transport marked points on all cell points
1928  for (const label facei : cFaces)
1929  {
1930  const face& f = mesh_.faces()[facei];
1931  face& baseF = baseFaces[facei];
1932 
1933  const label oldNChanged = nChanged;
1934  forAll(f, fp)
1935  {
1936  const auto fnd = minPointValue.find(f[fp]);
1937  if (fnd.good())
1938  {
1939  baseF.setSize(f.size(), labelMax);
1940  if (baseF[fp] == labelMax)
1941  {
1942  baseF[fp] = fnd();
1943  nChanged++;
1944 
1945  //Pout<< "For cell:" << celli
1946  // << " at:" << mesh_.cellCentres()[celli]
1947  // << " on face:" << facei
1948  // << " points:"
1949  // << UIndirectList<point>(mesh_.points(), f)
1950  // << " now have baseFace:" << baseF
1951  // << endl;
1952  }
1953  }
1954  }
1955 
1956  if (!isBlockedFace(facei) && nChanged > oldNChanged)
1957  {
1958  // Mark neighbouring cells
1959  const label own = mesh_.faceOwner()[facei];
1960  if (!isAffectedCell[own])
1961  {
1962  newIsAffectedCell.set(own);
1963  }
1964  if (mesh_.isInternalFace(facei))
1965  {
1966  const label nei = mesh_.faceNeighbour()[facei];
1967  if (!isAffectedCell[nei])
1968  {
1969  newIsAffectedCell.set(nei);
1970  }
1971  }
1972  }
1973  }
1974  }
1975 
1976  if (debug)
1977  {
1978  Pout<< "isAffectedCell:" << isAffectedCell.count() << endl;
1979  Pout<< "newIsAffectedCell:" << newIsAffectedCell.count()
1980  << endl;
1981  Pout<< "nChanged:" << nChanged << endl;
1982  }
1983 
1984  if (!returnReduceOr(nChanged))
1985  {
1986  break;
1987  }
1988 
1989 
1990  // Transport minimum across coupled faces
1991  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1992 
1993  SubList<face> l
1994  (
1995  baseFaces,
1996  mesh_.nBoundaryFaces(),
1997  mesh_.nInternalFaces()
1998  );
1999  oldBoundaryFaces = l;
2001  (
2002  mesh_,
2003  l,
2004  minEqOpFace(),
2005  Foam::dummyTransform() // dummy transformation
2006  );
2007 
2008  forAll(l, bFacei)
2009  {
2010  // Note: avoid special handling of comparing zero-sized faces
2011  // (see face::operator==). Review.
2012  const labelUList& baseVts = l[bFacei];
2013  const labelUList& oldVts = oldBoundaryFaces[bFacei];
2014  if (baseVts != oldVts)
2015  {
2016  const label facei = mesh_.nInternalFaces()+bFacei;
2017  const label own = mesh_.faceOwner()[facei];
2018  if (!isAffectedCell[own])
2019  {
2020  newIsAffectedCell.set(own);
2021  }
2022  }
2023  }
2024 
2025  isAffectedCell = newIsAffectedCell;
2026  }
2027  }
2028 
2029 
2030  //
2031  // Create new points
2032  //
2033 
2034  // If creating new mesh: copy existing patch points
2035  labelList copiedPatchPoints;
2036  if (!addToMesh_)
2037  {
2038  copiedPatchPoints.setSize(firstLayerDisp.size());
2039  forAll(firstLayerDisp, patchPointi)
2040  {
2041  if (addedPoints_[patchPointi].size())
2042  {
2043  label meshPointi = meshPoints[patchPointi];
2044  label zoneI = mesh_.pointZones().whichZone(meshPointi);
2045  copiedPatchPoints[patchPointi] = meshMod.setAction
2046  (
2047  polyAddPoint
2048  (
2049  mesh_.points()[meshPointi], // point
2050  -1, // master point
2051  zoneI, // zone for point
2052  true // supports a cell
2053  )
2054  );
2055  }
2056  }
2057  }
2058 
2059 
2060  // Create points for additional layers
2061  forAll(firstLayerDisp, patchPointi)
2062  {
2063  if (addedPoints_[patchPointi].size())
2064  {
2065  const label meshPointi = meshPoints[patchPointi];
2066  const label zoneI = mesh_.pointZones().whichZone(meshPointi);
2067 
2068  point pt = mesh_.points()[meshPointi];
2069 
2070  vector disp = firstLayerDisp[patchPointi];
2071 
2072  forAll(addedPoints_[patchPointi], i)
2073  {
2074  pt += disp;
2075 
2076  const label addedVertI = meshMod.setAction
2077  (
2078  polyAddPoint
2079  (
2080  pt, // point
2081  (addToMesh_ ? meshPointi : -1), // master point
2082  zoneI, // zone for point
2083  true // supports a cell
2084  )
2085  );
2086 
2087 
2088  //Pout<< "Adding point:" << addedVertI << " at:" << pt
2089  // << " from point:" << meshPointi
2090  // << " at:" << mesh_.points()[meshPointi]
2091  // << endl;
2092 
2093  addedPoints_[patchPointi][i] = addedVertI;
2094 
2095  disp *= expansionRatio[patchPointi];
2096  }
2097  }
2098  }
2099 
2100 
2101  //
2102  // Add cells to all boundaryFaces
2103  //
2104 
2105  labelListList addedCells(pp.size());
2106 
2107  forAll(pp, patchFacei)
2108  {
2109  if (nFaceLayers[patchFacei] > 0)
2110  {
2111  const label meshFacei = pp.addressing()[patchFacei];
2112 
2113  label extrudeCelli = -2;
2114  label extrudeZonei;
2115  if (!addToMesh_)
2116  {
2117  extrudeCelli = -1;
2118  const label ownCelli = mesh_.faceOwner()[meshFacei];
2119  extrudeZonei = mesh_.cellZones().whichZone(ownCelli);
2120  }
2121  else if (!ppFlip[patchFacei])
2122  {
2123  // Normal: extrude from owner face
2124  extrudeCelli = mesh_.faceOwner()[meshFacei];
2125  extrudeZonei = mesh_.cellZones().whichZone(extrudeCelli);
2126  }
2127  else if (mesh_.isInternalFace(meshFacei))
2128  {
2129  // Extrude from neighbour face (if internal). Might be
2130  // that it is a coupled face and the other side is
2131  // extruded
2132  extrudeCelli = mesh_.faceNeighbour()[meshFacei];
2133  extrudeZonei = mesh_.cellZones().whichZone(extrudeCelli);
2134  }
2135 
2136  if (extrudeCelli != -2)
2137  {
2138  addedCells[patchFacei].setSize(nFaceLayers[patchFacei]);
2139 
2140  for (label i = 0; i < nFaceLayers[patchFacei]; i++)
2141  {
2142  // Note: add from cell (owner of patch face) or from face?
2143  // for now add from cell so we can map easily.
2144  addedCells[patchFacei][i] = meshMod.setAction
2145  (
2146  polyAddCell
2147  (
2148  -1, // master point
2149  -1, // master edge
2150  -1, // master face
2151  extrudeCelli, // master cell
2152  extrudeZonei // zone for cell
2153  )
2154  );
2155 
2156  //Pout<< "Added cell:" << addedCells[patchFacei][i]
2157  // << " from master:" << extrudeCelli
2158  // << " at:" << mesh_.cellCentres()[extrudeCelli]
2159  // << endl;
2160  }
2161  }
2162  }
2163  }
2164 
2165 
2166 
2167  // Create faces on top of the original patch faces.
2168  // These faces are created from original patch faces outwards so in order
2169  // of increasing cell number. So orientation should be same as original
2170  // patch face for them to have owner<neighbour.
2171 
2172  layerFaces_.setSize(pp.size());
2173 
2174  forAll(pp.localFaces(), patchFacei)
2175  {
2176  label meshFacei = pp.addressing()[patchFacei];
2177 
2178  if (addedCells[patchFacei].size())
2179  {
2180  layerFaces_[patchFacei].setSize(addedCells[patchFacei].size() + 1);
2181 
2182  // Get duplicated vertices on the patch face.
2183  const face& f = pp.localFaces()[patchFacei];
2184 
2185  face newFace(f.size());
2186 
2187  forAll(addedCells[patchFacei], i)
2188  {
2189  forAll(f, fp)
2190  {
2191  if (addedPoints_[f[fp]].empty())
2192  {
2193  // Keep original point
2194  newFace[fp] =
2195  (
2196  addToMesh_
2197  ? meshPoints[f[fp]]
2198  : copiedPatchPoints[f[fp]]
2199  );
2200  }
2201  else
2202  {
2203  // Get new outside point
2204  label offset =
2205  addedPoints_[f[fp]].size()
2206  - addedCells[patchFacei].size();
2207  newFace[fp] = addedPoints_[f[fp]][i+offset];
2208  }
2209  }
2210  //Pout<< " newFace:" << newFace << endl;
2211  //Pout<< " coords:"
2212  // << UIndirectList<point>(meshMod.points(), newFace)
2213  // << " normal:" << newFace.unitNormal(meshMod.points())
2214  // << endl;
2215 
2216  // Get new neighbour
2217  label own = addedCells[patchFacei][i];
2218  label nei;
2219  label patchi;
2220  label zoneI = -1;
2221  bool flip = false;
2222  bool fluxFlip = false;
2223 
2224  if (i == addedCells[patchFacei].size()-1)
2225  {
2226  // Top layer so is either patch face or connects to
2227  // the other cell
2228  patchi = patchID[patchFacei];
2229  if (patchi == -1)
2230  {
2231  // Internal face
2232  nei =
2233  (
2234  !ppFlip[patchFacei]
2235  ? mesh_.faceNeighbour()[meshFacei]
2236  : mesh_.faceOwner()[meshFacei]
2237  );
2238 
2239  if (ppFlip[patchFacei])
2240  {
2241  newFace = newFace.reverseFace();
2242  }
2243  //Pout<< "** adding top (internal) face:"
2244  // << " at:" << mesh_.faceCentres()[meshFacei]
2245  // << " own:" << own << " nei:" << nei
2246  // << " patchi:" << patchi
2247  // << " newFace:" << newFace
2248  // << endl;
2249  }
2250  else
2251  {
2252  nei = -1;
2253  }
2254  zoneI = mesh_.faceZones().whichZone(meshFacei);
2255  if (zoneI != -1)
2256  {
2257  const faceZone& fz = mesh_.faceZones()[zoneI];
2258  flip = fz.flipMap()[fz.whichFace(meshFacei)];
2259  }
2260  }
2261  else
2262  {
2263  // Internal face between layer i and i+1
2264  nei = addedCells[patchFacei][i+1];
2265  patchi = -1;
2266  }
2267 
2268  if (nei != -1 && nei < own)
2269  {
2270  // Wrongly oriented internal face
2271  newFace = newFace.reverseFace();
2272  std::swap(own, nei);
2273  flip = !flip;
2274  fluxFlip = true;
2275 
2276  //Pout<< "Flipped newFace:"
2277  // << newFace.unitNormal(meshMod.points())
2278  // << " own:" << own
2279  // << " nei:" << nei
2280  // << endl;
2281  }
2282 
2283 
2284 
2285  layerFaces_[patchFacei][i+1] = meshMod.setAction
2286  (
2287  polyAddFace
2288  (
2289  newFace, // face
2290  own, // owner
2291  nei, // neighbour
2292  -1, // master point
2293  -1, // master edge
2294  (addToMesh_ ? meshFacei : -1), // master face
2295  fluxFlip, // flux flip
2296  patchi, // patch for face
2297  zoneI, // zone for face
2298  flip // face zone flip
2299  )
2300  );
2301 
2302  //Pout<< "added layer face:" << layerFaces_[patchFacei][i+1]
2303  // << " verts:" << newFace
2304  // << " newFc:" << newFace.centre(meshMod.points())
2305  // << " originalFc:" << mesh_.faceCentres()[meshFacei]
2306  // << nl
2307  // << " n:" << newFace.unitNormal(meshMod.points())
2308  // << " own:" << own << " nei:" << nei
2309  // << " patchi:" << patchi
2310  // << endl;
2311  }
2312  }
2313  }
2314 
2315  //
2316  // Modify owner faces to have addedCells as neighbour
2317  //
2318 
2319  if (addToMesh_)
2320  {
2321  forAll(pp, patchFacei)
2322  {
2323  if (addedCells[patchFacei].size())
2324  {
2325  label meshFacei = pp.addressing()[patchFacei];
2326 
2327  layerFaces_[patchFacei][0] = meshFacei;
2328  const face& f = pp[patchFacei];
2329 
2330  const label own =
2331  (
2332  !ppFlip[patchFacei]
2333  ? mesh_.faceOwner()[meshFacei]
2334  : mesh_.faceNeighbour()[meshFacei]
2335  );
2336  const label nei = addedCells[patchFacei][0];
2337 
2338  meshMod.setAction
2339  (
2340  polyModifyFace
2341  (
2342  (ppFlip[patchFacei] ? f.reverseFace() : f),// verts
2343  meshFacei, // label of face
2344  own, // owner
2345  nei, // neighbour
2346  ppFlip[patchFacei], // face flip
2347  -1, // patch for face
2348  true, //false, // remove from zone
2349  -1, //zoneI, // zone for face
2350  false // face flip in zone
2351  )
2352  );
2353 
2354  //Pout<< "Modified bottom face " << meshFacei
2355  // << " at:" << mesh_.faceCentres()[meshFacei]
2356  // << " new own:" << own << " new nei:" << nei
2357  // << " verts:" << meshMod.faces()[meshFacei]
2358  // << " n:"
2359  // << meshMod.faces()[meshFacei].unitNormal(meshMod.points())
2360  // << endl;
2361  }
2362  }
2363  }
2364  else
2365  {
2366  // If creating new mesh: reverse original faces and put them
2367  // in the exposed patch ID.
2368  forAll(pp, patchFacei)
2369  {
2370  if (addedCells[patchFacei].size())
2371  {
2372  label meshFacei = pp.addressing()[patchFacei];
2373  label zoneI = mesh_.faceZones().whichZone(meshFacei);
2374  bool zoneFlip = false;
2375  if (zoneI != -1)
2376  {
2377  const faceZone& fz = mesh_.faceZones()[zoneI];
2378  zoneFlip = !fz.flipMap()[fz.whichFace(meshFacei)];
2379  }
2380 
2381  // Reverse and renumber old patch face.
2382  face f(pp.localFaces()[patchFacei].reverseFace());
2383  forAll(f, fp)
2384  {
2385  f[fp] = copiedPatchPoints[f[fp]];
2386  }
2387 
2388  layerFaces_[patchFacei][0] = meshMod.setAction
2389  (
2390  polyAddFace
2391  (
2392  f, // modified face
2393  addedCells[patchFacei][0], // owner
2394  -1, // neighbour
2395  -1, // masterPoint
2396  -1, // masterEdge
2397  -1, // masterFace
2398  true, // face flip
2399  exposedPatchID[patchFacei], // patch for face
2400  zoneI, // zone for face
2401  zoneFlip // face flip in zone
2402  )
2403  );
2404  }
2405  }
2406  }
2407 
2408 
2409 
2410  //
2411  // Create 'side' faces, one per edge that is being extended.
2412  //
2413 
2414  const labelListList& faceEdges = pp.faceEdges();
2415  const faceList& localFaces = pp.localFaces();
2416  const edgeList& edges = pp.edges();
2417 
2418  // Get number of layers per edge. This is 0 if edge is not extruded;
2419  // max of connected faces otherwise.
2420  labelList edgeLayers(pp.nEdges());
2421 
2422  {
2423  // Use list over mesh.nEdges() since syncTools does not yet support
2424  // partial list synchronisation.
2425  labelList meshEdgeLayers(mesh_.nEdges(), -1);
2426 
2427  forAll(meshEdges, edgei)
2428  {
2429  const edge& e = edges[edgei];
2430 
2431  label meshEdgei = meshEdges[edgei];
2432 
2433  if ((nPointLayers[e[0]] == 0) && (nPointLayers[e[1]] == 0))
2434  {
2435  meshEdgeLayers[meshEdgei] = 0;
2436  }
2437  else
2438  {
2439  const labelList& eFaces = pp.edgeFaces()[edgei];
2440 
2441  forAll(eFaces, i)
2442  {
2443  meshEdgeLayers[meshEdgei] = max
2444  (
2445  nFaceLayers[eFaces[i]],
2446  meshEdgeLayers[meshEdgei]
2447  );
2448  }
2449  }
2450  }
2451 
2453  (
2454  mesh_,
2455  meshEdgeLayers,
2456  maxEqOp<label>(),
2457  label(0) // initial value
2458  );
2459 
2460  forAll(meshEdges, edgei)
2461  {
2462  edgeLayers[edgei] = meshEdgeLayers[meshEdges[edgei]];
2463  }
2464  }
2465 
2466 
2467  // Mark off which edges have been extruded
2468  boolList doneEdge(pp.nEdges(), false);
2469 
2470 
2471  // Create faces. Per face walk connected edges and find string of edges
2472  // between the same two faces and extrude string into a single face.
2473  forAll(pp, patchFacei)
2474  {
2475  const labelList& fEdges = faceEdges[patchFacei];
2476 
2477  forAll(fEdges, fp)
2478  {
2479  // Get string of edges that needs to be extruded as a single face.
2480  // Returned as indices in fEdges.
2481  labelPair indexPair
2482  (
2483  getEdgeString
2484  (
2485  pp,
2486  globalEdgeFaces,
2487  doneEdge,
2488  patchFacei,
2489  globalFaces.toGlobal(pp.addressing()[patchFacei])
2490  )
2491  );
2492 
2493  //Pout<< "Found unextruded edges in edges:" << fEdges
2494  // << " start:" << indexPair[0]
2495  // << " end:" << indexPair[1]
2496  // << endl;
2497 
2498  const label startFp = indexPair[0];
2499  const label endFp = indexPair[1];
2500 
2501  if (startFp != -1)
2502  {
2503  // Extrude edges from indexPair[0] up to indexPair[1]
2504  // (note indexPair = indices of edges. There is one more vertex
2505  // than edges)
2506  const face& f = localFaces[patchFacei];
2507 
2508  labelList stringedVerts;
2509  if (endFp >= startFp)
2510  {
2511  stringedVerts.setSize(endFp-startFp+2);
2512  }
2513  else
2514  {
2515  stringedVerts.setSize(endFp+f.size()-startFp+2);
2516  }
2517 
2518  label fp = startFp;
2519 
2520  for (label i = 0; i < stringedVerts.size()-1; i++)
2521  {
2522  stringedVerts[i] = f[fp];
2523  doneEdge[fEdges[fp]] = true;
2524  fp = f.fcIndex(fp);
2525  }
2526  stringedVerts.last() = f[fp];
2527 
2528 
2529  // Now stringedVerts contains the vertices in order of face f.
2530  // This is consistent with the order if f becomes the owner cell
2531  // and nbrFacei the neighbour cell. Note that the cells get
2532  // added in order of pp so we can just use face ordering and
2533  // because we loop in incrementing order as well we will
2534  // always have nbrFacei > patchFacei.
2535 
2536  label startEdgei = fEdges[startFp];
2537 
2538  label meshEdgei = meshEdges[startEdgei];
2539 
2540  label numEdgeSideFaces = edgeLayers[startEdgei];
2541 
2542  for (label i = 0; i < numEdgeSideFaces; i++)
2543  {
2544  label vEnd = stringedVerts.last();
2545  label vStart = stringedVerts[0];
2546 
2547  // calculate number of points making up a face
2548  label newFp = 2*stringedVerts.size();
2549 
2550  if (i == 0)
2551  {
2552  // layer 0 gets all the truncation of neighbouring
2553  // faces with more layers.
2554  if (addedPoints_[vEnd].size())
2555  {
2556  newFp +=
2557  addedPoints_[vEnd].size() - numEdgeSideFaces;
2558  }
2559  if (addedPoints_[vStart].size())
2560  {
2561  newFp +=
2562  addedPoints_[vStart].size() - numEdgeSideFaces;
2563  }
2564  }
2565 
2566  face newFace(newFp);
2567 
2568  newFp = 0;
2569 
2570  // For layer 0 get pp points, for all other layers get
2571  // points of layer-1.
2572  if (i == 0)
2573  {
2574  forAll(stringedVerts, stringedI)
2575  {
2576  label v = stringedVerts[stringedI];
2577  addVertex
2578  (
2579  (
2580  addToMesh_
2581  ? meshPoints[v]
2582  : copiedPatchPoints[v]
2583  ),
2584  newFace,
2585  newFp
2586  );
2587  }
2588  }
2589  else
2590  {
2591  forAll(stringedVerts, stringedI)
2592  {
2593  label v = stringedVerts[stringedI];
2594  if (addedPoints_[v].size())
2595  {
2596  label offset =
2597  addedPoints_[v].size() - numEdgeSideFaces;
2598  addVertex
2599  (
2600  addedPoints_[v][i+offset-1],
2601  newFace,
2602  newFp
2603  );
2604  }
2605  else
2606  {
2607  addVertex
2608  (
2609  (
2610  addToMesh_
2611  ? meshPoints[v]
2612  : copiedPatchPoints[v]
2613  ),
2614  newFace,
2615  newFp
2616  );
2617  }
2618  }
2619  }
2620 
2621  // add points between stringed vertices (end)
2622  if (numEdgeSideFaces < addedPoints_[vEnd].size())
2623  {
2624  if (i == 0 && addedPoints_[vEnd].size())
2625  {
2626  label offset =
2627  addedPoints_[vEnd].size() - numEdgeSideFaces;
2628  for (label ioff = 0; ioff < offset; ioff++)
2629  {
2630  addVertex
2631  (
2632  addedPoints_[vEnd][ioff],
2633  newFace,
2634  newFp
2635  );
2636  }
2637  }
2638  }
2639 
2640  forAllReverse(stringedVerts, stringedI)
2641  {
2642  label v = stringedVerts[stringedI];
2643  if (addedPoints_[v].size())
2644  {
2645  label offset =
2646  addedPoints_[v].size() - numEdgeSideFaces;
2647  addVertex
2648  (
2649  addedPoints_[v][i+offset],
2650  newFace,
2651  newFp
2652  );
2653  }
2654  else
2655  {
2656  addVertex
2657  (
2658  (
2659  addToMesh_
2660  ? meshPoints[v]
2661  : copiedPatchPoints[v]
2662  ),
2663  newFace,
2664  newFp
2665  );
2666  }
2667  }
2668 
2669 
2670  // add points between stringed vertices (start)
2671  if (numEdgeSideFaces < addedPoints_[vStart].size())
2672  {
2673  if (i == 0 && addedPoints_[vStart].size())
2674  {
2675  label offset =
2676  addedPoints_[vStart].size() - numEdgeSideFaces;
2677  for (label ioff = offset-1; ioff >= 0; ioff--)
2678  {
2679  addVertex
2680  (
2681  addedPoints_[vStart][ioff],
2682  newFace,
2683  newFp
2684  );
2685  }
2686  }
2687  }
2688 
2689  if (newFp >= 3)
2690  {
2691  // Add face inbetween faces patchFacei and nbrFacei
2692  // (possibly -1 for external edges)
2693 
2694  newFace.setSize(newFp);
2695 
2696  // Walked edges as if owner face was extruded. Reverse
2697  // for neighbour face extrusion.
2698  if (ppFlip[patchFacei])
2699  {
2700  newFace = newFace.reverseFace();
2701  }
2702 
2703  if (debug)
2704  {
2705  labelHashSet verts(2*newFace.size());
2706  forAll(newFace, fp)
2707  {
2708  if (!verts.insert(newFace[fp]))
2709  {
2711  << "Duplicate vertex in face"
2712  << " to be added." << nl
2713  << "newFace:" << newFace << nl
2714  << "points:"
2715  << UIndirectList<point>
2716  (
2717  meshMod.points(),
2718  newFace
2719  ) << nl
2720  << "Layer:" << i
2721  << " out of:" << numEdgeSideFaces << nl
2722  << "ExtrudeEdge:" << meshEdgei
2723  << " at:"
2724  << mesh_.edges()[meshEdgei].line
2725  (
2726  mesh_.points()
2727  ) << nl
2728  << "string:" << stringedVerts
2729  << "stringpoints:"
2730  << UIndirectList<point>
2731  (
2732  pp.localPoints(),
2733  stringedVerts
2734  ) << nl
2735  << "stringNLayers:"
2736  << labelUIndList
2737  (
2738  nPointLayers,
2739  stringedVerts
2740  ) << nl
2741  << abort(FatalError);
2742  }
2743  }
2744  }
2745 
2746  label nbrFacei = nbrFace
2747  (
2748  pp.edgeFaces(),
2749  startEdgei,
2750  patchFacei
2751  );
2752 
2753  const labelList& meshFaces = mesh_.edgeFaces
2754  (
2755  meshEdgei,
2756  ef
2757  );
2758 
2759  // Because we walk in order of patch face and in order
2760  // of face edges so face orientation will be opposite
2761  // that of the patch edge
2762  bool zoneFlip = false;
2763  if (edgeZoneID[startEdgei] != -1)
2764  {
2765  zoneFlip = !edgeFlip[startEdgei];
2766  }
2767 
2768  addSideFace
2769  (
2770  pp,
2771  addedCells,
2772 
2773  newFace, // vertices of new face
2774  edgePatchID[startEdgei],// -1 or patch for face
2775  edgeZoneID[startEdgei],
2776  zoneFlip,
2777  inflateFaceID[startEdgei],
2778 
2779  patchFacei,
2780  nbrFacei,
2781  meshEdgei, // (mesh) edge to inflate
2782  i, // layer
2783  numEdgeSideFaces, // num layers
2784  meshFaces, // edgeFaces
2785  meshMod
2786  );
2787  }
2788  }
2789  }
2790  }
2791  }
2792 
2793 
2794  // Adjust side faces if they're on the side where points were duplicated
2795  // (i.e. adding to internal faces). Like duplicatePoints::setRefinement.
2796  if (addToMesh_)
2797  {
2798  face newFace;
2799 
2800  forAll(baseFaces, facei)
2801  {
2802  const face& f = mesh_.faces()[facei];
2803  const face& baseF = baseFaces[facei];
2804 
2805  if (isBlockedFace(facei) || baseF.empty())
2806  {
2807  // Either part of patch or no duplicated points on face
2808  continue;
2809  }
2810 
2811  // Start off from original face
2812  newFace = f;
2813  forAll(f, fp)
2814  {
2815  const label meshPointi = f[fp];
2816  if (baseF[fp] != labelMax)
2817  {
2818  // Duplicated point
2819  const label patchPointi = pp.meshPointMap()[meshPointi];
2820  const label addedPointi = addedPoints_[patchPointi].last();
2821 
2822  //Pout<< " For point:" << meshPointi
2823  // << " at:" << mesh_.points()[meshPointi]
2824  // << " at:" << pp.localPoints()[patchPointi]
2825  // << " using addedpoint:" << addedPointi
2826  // << " at:" << meshMod.points()[addedPointi]
2827  // << endl;
2828 
2829  newFace[fp] = addedPointi;
2830  }
2831  }
2832 
2833  //Pout<< "for face:" << facei << nl
2834  // << " old:" << f
2835  // << " n:" << newFace.unitNormal(meshMod.points())
2836  // << " coords:" << UIndirectList<point>(mesh_.points(), f)
2837  // << nl
2838  // << " new:" << newFace
2839  // << " n:" << newFace.unitNormal(meshMod.points())
2840  // << " coords:"
2841  // << UIndirectList<point>(meshMod.points(), newFace)
2842  // << endl;
2843 
2844  // Get current zone info
2845  label zoneID = mesh_.faceZones().whichZone(facei);
2846  bool zoneFlip = false;
2847  if (zoneID >= 0)
2848  {
2849  const faceZone& fZone = mesh_.faceZones()[zoneID];
2850  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
2851  }
2852 
2853 
2854  if (mesh_.isInternalFace(facei))
2855  {
2856  meshMod.modifyFace
2857  (
2858  newFace, // modified face
2859  facei, // label of modified face
2860  mesh_.faceOwner()[facei], // owner
2861  mesh_.faceNeighbour()[facei], // neighbour
2862  false, // face flip
2863  -1, // patch for face
2864  zoneID, // zone for face
2865  zoneFlip // face flip in zone
2866  );
2867  }
2868  else
2869  {
2870  meshMod.modifyFace
2871  (
2872  newFace, // modified face
2873  facei, // label of modified face
2874  mesh_.faceOwner()[facei], // owner
2875  -1, // neighbour
2876  false, // face flip
2877  patches.whichPatch(facei), // patch for face
2878  zoneID, // zone for face
2879  zoneFlip // face flip in zone
2880  );
2881  }
2882  }
2883  }
2884 }
2885 
2886 
2888 (
2889  const mapPolyMesh& morphMap,
2890  const labelList& faceMap, // new to old patch faces
2891  const labelList& pointMap // new to old patch points
2892 )
2893 {
2894  {
2895  labelListList newAddedPoints(pointMap.size());
2896 
2897  forAll(newAddedPoints, newPointi)
2898  {
2899  label oldPointi = pointMap[newPointi];
2900 
2901  const labelList& added = addedPoints_[oldPointi];
2902 
2903  labelList& newAdded = newAddedPoints[newPointi];
2904  newAdded.setSize(added.size());
2905  label newI = 0;
2906 
2907  forAll(added, i)
2908  {
2909  label newPointi = morphMap.reversePointMap()[added[i]];
2910 
2911  if (newPointi >= 0)
2912  {
2913  newAdded[newI++] = newPointi;
2914  }
2915  }
2916  newAdded.setSize(newI);
2917  }
2918  addedPoints_.transfer(newAddedPoints);
2919  }
2920 
2921  {
2922  labelListList newLayerFaces(faceMap.size());
2923 
2924  forAll(newLayerFaces, newFacei)
2925  {
2926  label oldFacei = faceMap[newFacei];
2927 
2928  const labelList& added = layerFaces_[oldFacei];
2929 
2930  labelList& newAdded = newLayerFaces[newFacei];
2931  newAdded.setSize(added.size());
2932  label newI = 0;
2933 
2934  forAll(added, i)
2935  {
2936  label newFacei = morphMap.reverseFaceMap()[added[i]];
2937 
2938  if (newFacei >= 0)
2939  {
2940  newAdded[newI++] = newFacei;
2941  }
2942  }
2943  newAdded.setSize(newI);
2944  }
2945  layerFaces_.transfer(newLayerFaces);
2946  }
2947 }
2948 
2949 
2950 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
label nPoints() const
Number of points supporting patch faces.
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=UPstream::parRun())
Synchronize values on boundary faces only.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const Field< point_type > & localPoints() const
Return pointField of points in patch.
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
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
const labelListList & pointEdges() const
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:489
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static void matchEdges(const PrimitivePatch< FaceList1, PointField1 > &p1, const PrimitivePatch< FaceList2, PointField2 > &p2, labelList &p1EdgeLabels, labelList &p2EdgeLabels, bitSet &sameOrientation)
Find corresponding edges on patches sharing the same points.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void calcExtrudeInfo(const bool zoneFromAnyFace, const polyMesh &, const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const indirectPrimitivePatch &pp, labelList &edgePatchID, label &nPatches, Map< label > &nbrProcToPatch, Map< label > &patchToNbrProc, labelList &edgeZoneID, boolList &edgeFlip, labelList &inflateFaceID)
Determine extrude information per patch edge:
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
label nInternalEdges() const
Number of internal edges.
label k
Boltzmann constant.
label nFaces() const noexcept
Number of mesh faces.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
labelList faceLabels(nFaceLabels)
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Get per face whether it is internal or coupled.
Definition: syncTools.C:169
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static void syncData(List< Type > &elems, const labelListList &slaves, const labelListList &transformedSlaves, const mapDistribute &slavesMap, const globalIndexAndTransform &, const CombineOp &cop, const TransformOp &top)
Helper: synchronise data with transforms.
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
Dummy transform to be used with syncTools.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
scalar y
A list of faces which address into the list of points.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
const Map< label > & meshPointMap() const
Mesh point map.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
Adds layers of cells to outside of polyPatch. Can optionally create stand-alone extruded mesh (addToM...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const Field< point_type > & faceCentres() const
Return face centres for patch.
const labelListList & edgeFaces() const
Return edge-face addressing.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Vector< scalar > vector
Definition: vector.H:57
const Field< point_type > & points() const noexcept
Return reference to global points.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
void setRefinement(const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const scalarField &expansionRatio, const indirectPrimitivePatch &pp, const bitSet &flip, const labelList &sidePatchID, const labelList &sideZoneID, const boolList &sideFlip, const labelList &inflateFaceID, const labelList &exposedPatchID, const labelList &nFaceLayers, const labelList &nPointLayers, const vectorField &firstLayerDisp, polyTopoChange &meshMod)
Play commands into polyTopoChange to create layers on top.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
label nEdges() const
Number of mesh edges.
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
label nEdges() const
Number of edges in patch.
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
Type gMax(const FieldField< Field, Type > &f)
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:867
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:346
void updateMesh(const mapPolyMesh &, const labelList &faceMap, const labelList &pointMap)
Update any locally stored mesh information. Gets additional.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
vector point
Point is a vector.
Definition: point.H:37
label toGlobal(const label proci, const label i) const
From local to global on proci.
Definition: globalIndexI.H:308
#define WarningInFunction
Report a warning using Foam::Warning.
label newPointi
Definition: readKivaGrid.H:496
void operator()(face &x, const face &y) const
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Return labels of patch edges in the global edge list using cell addressing.
const polyBoundaryMesh & patches
const std::string patch
OpenFOAM patch number as a std::string.
constexpr label labelMax
Definition: label.H:55
const labelListList & faceEdges() const
Return face-edge addressing.
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Per patch edge the pp faces (in global indices) using it.
label n
labelListList addedCells() const
Added cells given current mesh & layerfaces.
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:437
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
List< label > labelList
A List of labels.
Definition: List.H:62
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< bool > boolList
A List of bools.
Definition: List.H:60
bool found
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const labelListList & edgeFaces() const
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:485
static int compare(const edge &a, const edge &b)
Compare edges.
Definition: edgeI.H:26
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127