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-2017,2020-2021 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 {
47  defineTypeNameAndDebug(addPatchCellLayer, 0);
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 Foam::label Foam::addPatchCellLayer::nbrFace
54 (
55  const labelListList& edgeFaces,
56  const label edgei,
57  const label facei
58 )
59 {
60  const labelList& eFaces = edgeFaces[edgei];
61 
62  if (eFaces.size() == 2)
63  {
64  return (eFaces[0] != facei ? eFaces[0] : eFaces[1]);
65  }
66  else
67  {
68  return -1;
69  }
70 }
71 
72 
73 void Foam::addPatchCellLayer::addVertex
74 (
75  const label pointi,
76  face& f,
77  label& fp
78 )
79 {
80  if (fp == 0)
81  {
82  f[fp++] = pointi;
83  }
84  else
85  {
86  if (f[fp-1] != pointi && f[0] != pointi)
87  {
88  f[fp++] = pointi;
89  }
90  }
91 }
92 
93 
94 // Is edge to the same neighbour? (and needs extrusion and has not been
95 // dealt with already)
96 bool Foam::addPatchCellLayer::sameEdgeNeighbour
97 (
98  const indirectPrimitivePatch& pp,
99  const labelListList& globalEdgeFaces,
100  const boolList& doneEdge,
101  const label thisGlobalFacei,
102  const label nbrGlobalFacei,
103  const label edgei
104 ) const
105 {
106  const edge& e = pp.edges()[edgei];
107 
108  return
109  !doneEdge[edgei] // not yet handled
110  && (
111  addedPoints_[e[0]].size() // is extruded
112  || addedPoints_[e[1]].size()
113  )
114  && (
115  nbrFace(globalEdgeFaces, edgei, thisGlobalFacei)
116  == nbrGlobalFacei // is to same neighbour
117  );
118 }
119 
120 
121 // Collect consecutive string of edges that connects the same two
122 // (possibly coupled) faces. Returns -1 if no unvisited edge can be found.
123 // Otherwise returns start and end index in face.
124 Foam::labelPair Foam::addPatchCellLayer::getEdgeString
125 (
126  const indirectPrimitivePatch& pp,
127  const labelListList& globalEdgeFaces,
128  const boolList& doneEdge,
129  const label patchFacei,
130  const label globalFacei
131 ) const
132 {
133  const labelList& fEdges = pp.faceEdges()[patchFacei];
134 
135  label startFp = -1;
136  label endFp = -1;
137 
138  // Get edge that hasn't been done yet but needs extrusion
139  forAll(fEdges, fp)
140  {
141  label edgei = fEdges[fp];
142  const edge& e = pp.edges()[edgei];
143 
144  if
145  (
146  !doneEdge[edgei]
147  && ( addedPoints_[e[0]].size() || addedPoints_[e[1]].size() )
148  )
149  {
150  startFp = fp;
151  break;
152  }
153  }
154 
155  if (startFp != -1)
156  {
157  // We found an edge that needs extruding but hasn't been done yet.
158  // Now find the face on the other side
159  label nbrGlobalFacei = nbrFace
160  (
161  globalEdgeFaces,
162  fEdges[startFp],
163  globalFacei
164  );
165 
166  if (nbrGlobalFacei == -1)
167  {
168  // Proper boundary edge. Only extrude single edge.
169  endFp = startFp;
170  }
171  else
172  {
173  // Search back for edge
174  // - which hasn't been handled yet
175  // - with same neighbour
176  // - that needs extrusion
177 
178  const label initFp = startFp;
179  while (true)
180  {
181  label prevFp = fEdges.rcIndex(startFp);
182 
183  if (prevFp == initFp)
184  {
185  const edge& e = pp.edges()[fEdges[initFp]];
186  const face& localF = pp.localFaces()[patchFacei];
187 
189  << "On face:" << patchFacei
190  << " fc:" << pp.faceCentres()[patchFacei]
191  << " vertices:" << localF
192  << " points:"
193  << UIndirectList<point>(pp.points(), pp[patchFacei])
194  << " edges:" << fEdges
195  << " All edges of face seem to have same neighbour "
196  << nbrGlobalFacei
197  << " starting walking from edge " << e
198  << exit(FatalError);
199  }
200 
201  if
202  (
203  !sameEdgeNeighbour
204  (
205  pp,
206  globalEdgeFaces,
207  doneEdge,
208  globalFacei,
209  nbrGlobalFacei,
210  fEdges[prevFp]
211  )
212  )
213  {
214  break;
215  }
216  startFp = prevFp;
217  }
218 
219  // Search forward for end of string
220  endFp = startFp;
221  while (true)
222  {
223  label nextFp = fEdges.fcIndex(endFp);
224 
225  if
226  (
227  !sameEdgeNeighbour
228  (
229  pp,
230  globalEdgeFaces,
231  doneEdge,
232  globalFacei,
233  nbrGlobalFacei,
234  fEdges[nextFp]
235  )
236  )
237  {
238  break;
239  }
240  endFp = nextFp;
241  }
242  }
243  }
244 
245  return labelPair(startFp, endFp);
246 }
247 
248 
249 Foam::label Foam::addPatchCellLayer::addSideFace
250 (
251  const indirectPrimitivePatch& pp,
252  const labelListList& addedCells, // per pp face the new extruded cell
253  const face& newFace,
254  const label newPatchID,
255  const label zonei,
256  const bool edgeFlip,
257  const label inflateFacei,
258 
259  const label ownFacei, // pp face that provides owner
260  const label nbrFacei,
261  const label meshEdgei, // corresponding mesh edge
262  const label layeri, // layer
263  const label numEdgeFaces, // number of layers for edge
264  const labelList& meshFaces, // precalculated edgeFaces
265  polyTopoChange& meshMod
266 ) const
267 {
268  // Adds a side face i.e. extrudes a patch edge.
269 
270  label addedFacei = -1;
271 
272 
273  // Is patch edge external edge of indirectPrimitivePatch?
274  if (nbrFacei == -1)
275  {
276  // External edge so external face.
277 
278  // Determine if different number of layer on owner and neighbour side
279  // (relevant only for coupled faces). See section for internal edge
280  // below.
281 
282  label layerOwn;
283 
284  if (addedCells[ownFacei].size() < numEdgeFaces)
285  {
286  label offset = numEdgeFaces - addedCells[ownFacei].size();
287  if (layeri <= offset)
288  {
289  layerOwn = 0;
290  }
291  else
292  {
293  layerOwn = layeri - offset;
294  }
295  }
296  else
297  {
298  layerOwn = layeri;
299  }
300 
301 
302  //Pout<< "Added boundary face:" << newFace
303  // << " own:" << addedCells[ownFacei][layerOwn]
304  // << " patch:" << newPatchID
305  // << endl;
306 
307  addedFacei = meshMod.setAction
308  (
309  polyAddFace
310  (
311  newFace, // face
312  addedCells[ownFacei][layerOwn], // owner
313  -1, // neighbour
314  -1, // master point
315  -1, // master edge
316  inflateFacei, // master face
317  false, // flux flip
318  newPatchID, // patch for face
319  zonei, // zone for face
320  edgeFlip // face zone flip
321  )
322  );
323  }
324  else
325  {
326  // When adding side faces we need to modify neighbour and owners
327  // in region where layer mesh is stopped. Determine which side
328  // has max number of faces and make sure layers match closest to
329  // original pp if there are different number of layers.
330 
331  label layerNbr;
332  label layerOwn;
333 
334  if (addedCells[ownFacei].size() > addedCells[nbrFacei].size())
335  {
336  label offset =
337  addedCells[ownFacei].size() - addedCells[nbrFacei].size();
338 
339  layerOwn = layeri;
340 
341  if (layeri <= offset)
342  {
343  layerNbr = 0;
344  }
345  else
346  {
347  layerNbr = layeri - offset;
348  }
349  }
350  else if (addedCells[nbrFacei].size() > addedCells[ownFacei].size())
351  {
352  label offset =
353  addedCells[nbrFacei].size() - addedCells[ownFacei].size();
354 
355  layerNbr = layeri;
356 
357  if (layeri <= offset)
358  {
359  layerOwn = 0;
360  }
361  else
362  {
363  layerOwn = layeri - offset;
364  }
365  }
366  else
367  {
368  // Same number of layers on both sides.
369  layerNbr = layeri;
370  layerOwn = layeri;
371  }
372 
373 
374  // Check mesh internal faces using edge to initialise
375  label inflateEdgei = -1;
376  if (addToMesh_)
377  {
378  forAll(meshFaces, i)
379  {
380  if (mesh_.isInternalFace(meshFaces[i]))
381  {
382  // meshEdge uses internal faces so ok to inflate from it
383  inflateEdgei = meshEdgei;
384  break;
385  }
386  }
387  }
388 
389 
390  addedFacei = meshMod.setAction
391  (
392  polyAddFace
393  (
394  newFace, // face
395  addedCells[ownFacei][layerOwn], // owner
396  addedCells[nbrFacei][layerNbr], // neighbour
397  -1, // master point
398  inflateEdgei, // master edge
399  -1, // master face
400  false, // flux flip
401  -1, // patch for face
402  zonei, // zone for face
403  edgeFlip // face zone flip
404  )
405  );
406 
407  //Pout<< "Added internal face:" << newFace
408  // << " own:" << addedCells[ownFacei][layerOwn]
409  // << " nei:" << addedCells[nbrFacei][layerNbr]
410  // << endl;
411  }
412 
413  return addedFacei;
414 }
415 
416 
417 Foam::label Foam::addPatchCellLayer::findProcPatch
418 (
419  const polyMesh& mesh,
420  const label nbrProcID
421 )
422 {
423  const polyBoundaryMesh& patches = mesh.boundaryMesh();
424 
426  {
427  label patchi = mesh.globalData().processorPatches()[i];
428 
429  if
430  (
431  refCast<const processorPolyPatch>(patches[patchi]).neighbProcNo()
432  == nbrProcID
433  )
434  {
435  return patchi;
436  }
437  }
438  return -1;
439 }
440 
441 
442 void Foam::addPatchCellLayer::setFaceProps
443 (
444  const polyMesh& mesh,
445  const label facei,
446 
447  label& patchi,
448  label& zonei,
449  bool& zoneFlip
450 )
451 {
452  patchi = mesh.boundaryMesh().whichPatch(facei);
453  zonei = mesh.faceZones().whichZone(facei);
454  if (zonei != -1)
455  {
456  label index = mesh.faceZones()[zonei].whichFace(facei);
457  zoneFlip = mesh.faceZones()[zonei].flipMap()[index];
458  }
459 }
460 
461 
462 void Foam::addPatchCellLayer::setFaceProps
463 (
464  const polyMesh& mesh,
465  const indirectPrimitivePatch& pp,
466  const label ppEdgeI,
467  const label faceI,
468 
469  label& patchI,
470  label& zoneI,
471  bool& zoneFlip,
472  label& inflateFaceI
473 )
474 {
475  setFaceProps
476  (
477  mesh,
478  faceI,
479 
480  patchI,
481  zoneI,
482  zoneFlip
483  );
484 
485  if (patchI != -1 || zoneI != -1)
486  {
487  inflateFaceI = faceI;
488  }
489 
490  if (zoneI != -1)
491  {
492  // Correct flip for patch edge ordering
493  const edge& ppEdge = pp.edges()[ppEdgeI];
494  const edge patchEdge
495  (
496  pp.meshPoints()[ppEdge[0]],
497  pp.meshPoints()[ppEdge[1]]
498  );
499 
500  const face& f = mesh.faces()[faceI];
501  bool found = false;
502  forAll(f, fp)
503  {
504  const edge e(f[fp], f.nextLabel(fp));
505  int stat = edge::compare(e, patchEdge);
506  if (stat == 1)
507  {
508  found = true;
509  break;
510  }
511  else if (stat == -1)
512  {
513  found = true;
514  zoneFlip = !zoneFlip;
515  break;
516  }
517  }
518 
519  if (!found)
520  {
521  //FatalErrorInFunction
523  << "Problem: cannot find patch edge " << ppEdgeI
524  << " with mesh vertices " << patchEdge
525  << " at " << patchEdge.line(mesh.points())
526  << " in face " << faceI << " with mesh vertices "
527  << f
528  << " at " << pointField(mesh.points(), f)
529  << endl
530  << "Continuing with potentially incorrect faceZone orientation"
531  //<< exit(FatalError);
532  << endl;
533  }
534  }
535 }
536 
537 
538 void Foam::addPatchCellLayer::findZoneFace
539 (
540  const bool useInternalFaces,
541  const bool useBoundaryFaces,
542 
543  const polyMesh& mesh,
544  const indirectPrimitivePatch& pp,
545  const label ppEdgeI,
546  const labelUIndList& excludeFaces,
547  const labelList& meshFaces,
548 
549  label& inflateFaceI,
550  label& patchI,
551  label& zoneI,
552  bool& zoneFlip
553 )
554 {
555  inflateFaceI = -1;
556  patchI = -1;
557  zoneI = -1;
558  zoneFlip = false;
559 
560  forAll(meshFaces, k)
561  {
562  label faceI = meshFaces[k];
563 
564  if
565  (
566  !excludeFaces.found(faceI)
567  && (
568  (mesh.isInternalFace(faceI) && useInternalFaces)
569  || (!mesh.isInternalFace(faceI) && useBoundaryFaces)
570  )
571  )
572  {
573  setFaceProps
574  (
575  mesh,
576  pp,
577  ppEdgeI,
578  faceI,
579 
580  patchI,
581  zoneI,
582  zoneFlip,
583  inflateFaceI
584  );
585 
586  if (zoneI != -1 || patchI != -1)
587  {
588  break;
589  }
590  }
591  }
592 }
593 
594 
595 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
596 
597 // Construct from mesh
598 Foam::addPatchCellLayer::addPatchCellLayer
599 (
600  const polyMesh& mesh,
601  const bool addToMesh
602 )
603 :
604  mesh_(mesh),
605  addToMesh_(addToMesh),
606  addedPoints_(0),
607  layerFaces_(0)
608 {}
609 
610 
611 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
612 
614 (
615  const polyMesh& mesh,
616  const labelListList& layerFaces
617 )
618 {
619  labelListList layerCells(layerFaces.size());
620 
621  forAll(layerFaces, patchFacei)
622  {
623  const labelList& faceLabels = layerFaces[patchFacei];
624 
625  if (faceLabels.size())
626  {
627  labelList& added = layerCells[patchFacei];
628  added.setSize(faceLabels.size()-1);
629 
630  for (label i = 0; i < faceLabels.size()-1; i++)
631  {
632  added[i] = mesh.faceNeighbour()[faceLabels[i]];
633  }
634  }
635  }
636  return layerCells;
637 }
638 
639 
641 {
642  return addedCells(mesh_, layerFaces_);
643 }
644 
645 
647 (
648  const polyMesh& mesh,
649  const globalIndex& globalFaces,
650  const indirectPrimitivePatch& pp
651 )
652 {
653  // Precalculate mesh edges for pp.edges.
654  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
655 
656  // From mesh edge to global face labels. Non-empty sublists only for
657  // pp edges.
658  labelListList globalEdgeFaces(mesh.nEdges());
659 
660  const labelListList& edgeFaces = pp.edgeFaces();
661 
662  forAll(edgeFaces, edgeI)
663  {
664  label meshEdgeI = meshEdges[edgeI];
665 
666  const labelList& eFaces = edgeFaces[edgeI];
667 
668  // Store face and processor as unique tag.
669  labelList& globalEFaces = globalEdgeFaces[meshEdgeI];
670  globalEFaces.setSize(eFaces.size());
671  forAll(eFaces, i)
672  {
673  globalEFaces[i] = globalFaces.toGlobal(pp.addressing()[eFaces[i]]);
674  }
675  }
676 
677  // Synchronise across coupled edges.
679  (
680  mesh,
681  globalEdgeFaces,
683  labelList() // null value
684  );
685 
686  // Extract pp part
687  return labelListList(UIndirectList<labelList>(globalEdgeFaces, meshEdges));
688 }
689 
690 
691 void Foam::addPatchCellLayer::markPatchEdges
692 (
693  const polyMesh& mesh,
694  const indirectPrimitivePatch& pp,
695  const labelListList& edgeGlobalFaces,
696  const labelList& meshEdges,
697 
698  bitSet& isPatchEdge,
699  bitSet& isPatchBoundaryEdge
700 )
701 {
702  // Mark (mesh) edges:
703  // - anywhere on extrusion
704  // - where the extrusion ends
705 
706  isPatchEdge.setSize(mesh.nEdges());
707  isPatchEdge = false;
708  isPatchEdge.set(meshEdges);
709  // Synchronise across coupled edges
711  (
712  mesh,
713  isPatchEdge,
715  false // initial value
716  );
717 
718  isPatchBoundaryEdge.setSize(mesh.nEdges());
719  isPatchBoundaryEdge = false;
720  forAll(edgeGlobalFaces, edgei)
721  {
722  // Test that edge has single global extruded face.
723  // Mark on processor that holds the face (since edgeGlobalFaces
724  // only gets filled from pp faces so if there is only one this
725  // is it)
726  if (edgeGlobalFaces[edgei].size() == 1)
727  {
728  isPatchBoundaryEdge.set(meshEdges[edgei]);
729  }
730  }
731  // Synchronise across coupled edges
733  (
734  mesh,
735  isPatchBoundaryEdge,
736  orEqOp<unsigned int>(),
737  false // initial value
738  );
739 }
740 
741 
742 void Foam::addPatchCellLayer::globalEdgeInfo
743 (
744  const bool zoneFromAnyFace,
745 
746  const polyMesh& mesh,
747  const globalIndex& globalFaces,
748  const labelListList& edgeGlobalFaces,
749  const indirectPrimitivePatch& pp,
750  const labelList& meshEdges,
751 
752  labelList& patchEdgeToFace, // face (in globalFaces index)
753  labelList& patchEdgeToPatch, // patch on face (or -1 for internal faces)
754  labelList& patchEdgeToZone, // zone on face
755  bitSet& patchEdgeToFlip // flip orientation on face
756 )
757 {
758  // For every edge on the outside of the patch return a potential patch/
759  // faceZone to extrude into.
760 
761  // Mark (mesh) edges on pp.
762  bitSet isExtrudeEdge;
763  bitSet isBoundaryEdge;
764  markPatchEdges
765  (
766  mesh,
767  pp,
768  edgeGlobalFaces,
769  meshEdges,
770 
771  isExtrudeEdge,
772  isBoundaryEdge
773  );
774 
775  // Build map of pp edges (in mesh point indexing). Note that this
776  // is now also on processors that do not have pp (but do have the edge)
777  EdgeMap<label> isBoundaryEdgeSet(pp.nEdges());
778  for (const label edgei : isBoundaryEdge)
779  {
780  isBoundaryEdgeSet.insert(mesh.edges()[edgei], edgei);
781  }
782  EdgeMap<label> isExtrudeEdgeSet(pp.nEdges());
783  for (const label edgei : isExtrudeEdge)
784  {
785  isExtrudeEdgeSet.insert(mesh.edges()[edgei], edgei);
786  }
787 
788 
789  const faceZoneMesh& fzs = mesh.faceZones();
790 
791  // Extract zone info into mesh face indexing for ease of addressing
792  labelList faceToZone(mesh.nFaces(), -1);
793  bitSet faceToFlip(mesh.nFaces());
794  for (const faceZone& fz: fzs)
795  {
796  const labelList& addressing = fz;
797  UIndirectList<label>(faceToZone, addressing) = fz.index();
798 
799  const boolList& fm = fz.flipMap();
800  forAll(addressing, i)
801  {
802  faceToFlip[addressing[i]] = fm[i];
803  }
804  }
805 
806 
807  // Storage (over all mesh edges)
808  // - face that data originates from (in globalFaces indexing)
809  labelList meshEdgeToFace(mesh.nEdges(), -1);
810  // - patch (for boundary faces)
811  labelList meshEdgeToPatch(mesh.nEdges(), -1);
812  // - faceZone
813  labelList meshEdgeToZone(mesh.nEdges(), -1);
814  // - faceZone orientation
815  bitSet meshEdgeToFlip(mesh.nEdges());
816 
817  //if (useInternalFaces)
818  {
819  const bitSet isInternalOrCoupled
820  (
822  );
823 
824  // Loop over edges of the face to find any faceZone.
825  // Edges kept as point pair so we don't construct mesh.faceEdges etc.
826 
827  for (const label facei : isInternalOrCoupled)
828  {
829  const face& f = mesh.faces()[facei];
830 
831  label prevPointi = f.last();
832  for (const label pointi : f)
833  {
834  const edge e(prevPointi, pointi);
835 
836  // Check if edge is internal to extrusion. Take over faceZone
837  // etc from internal face.
838  const auto eFnd = isExtrudeEdgeSet.cfind(e);
839  if (eFnd)
840  {
841  const label edgei = eFnd();
842 
843  if (faceToZone[facei] != -1)
844  {
845  // Found a zoned internal face. Use.
846  meshEdgeToFace[edgei] = globalFaces.toGlobal(facei);
847  meshEdgeToZone[edgei] = faceToZone[facei];
848  const edge& meshE = mesh.edges()[edgei];
849  const int d = edge::compare(e, meshE);
850  if (d == 1)
851  {
852  meshEdgeToFlip[edgei] = faceToFlip[facei];
853  }
854  else if (d == -1)
855  {
856  meshEdgeToFlip[edgei] = !faceToFlip[facei];
857  }
858  else
859  {
860  FatalErrorInFunction << "big problem"
861  << exit(FatalError);
862  }
863  }
864  }
865 
866  prevPointi = pointi;
867  }
868  }
869  }
870 
871 
872  //if (useBoundaryFaces)
873  {
874  // Loop over all patches and find 'best' one (non-coupled,
875  // non-extrusion, non-constraint(?)). Note that logic is slightly
876  // different from internal faces loop above - first patch face
877  // is being used instead of first zoned face for internal faces
878 
879  const polyBoundaryMesh& patches = mesh.boundaryMesh();
880 
881  bitSet isPpFace(mesh.nFaces());
882  isPpFace.set(pp.addressing());
883  // Note: no need to sync ppFace since does not include processor patches
884 
885  for (const polyPatch& pp : patches)
886  {
887  if (!pp.coupled())
888  {
889  // TBD. Check for constraint? This is usually a geometric
890  // constraint and not a topological one so should
891  // be handled in the extrusion vector calculation instead?
892 
893  forAll(pp, i)
894  {
895  const label facei = pp.start()+i;
896 
897  if (!isPpFace[facei])
898  {
899  const face& f = pp[i];
900 
901  label prevPointi = f.last();
902  for (const label pointi : f)
903  {
904  const edge e(prevPointi, pointi);
905  const auto eFnd =
906  (
907  zoneFromAnyFace
908  ? isExtrudeEdgeSet.cfind(e)
909  : isBoundaryEdgeSet.cfind(e)
910  );
911  if (eFnd)
912  {
913  const label edgei = eFnd();
914  if (meshEdgeToFace[edgei] == -1)
915  {
916  // Found unassigned face. Use its
917  // information.
918  // Note that we use the lowest numbered
919  // patch face.
920 
921  meshEdgeToFace[edgei] =
922  globalFaces.toGlobal(facei);
923 
924  // Override any patch info
925  if (meshEdgeToPatch[edgei] == -1)
926  {
927  meshEdgeToPatch[edgei] = pp.index();
928  }
929 
930  // Override any zone info
931  if (meshEdgeToZone[edgei] == -1)
932  {
933  meshEdgeToZone[edgei] =
934  faceToZone[facei];
935  const edge& meshE = mesh.edges()[edgei];
936  const int d = edge::compare(e, meshE);
937  if (d == 1)
938  {
939  meshEdgeToFlip[edgei] =
940  faceToFlip[facei];
941  }
942  else if (d == -1)
943  {
944  meshEdgeToFlip[edgei] =
945  !faceToFlip[facei];
946  }
947  else
948  {
950  << "big problem"
951  << exit(FatalError);
952  }
953  }
954  }
955  }
956 
957  prevPointi = pointi;
958  }
959  }
960  }
961  }
962  }
963  }
964 
965  // Synchronise across coupled edges. Max patch/face/faceZone wins
967  (
968  mesh,
969  meshEdgeToFace,
970  maxEqOp<label>(),
971  label(-1)
972  );
974  (
975  mesh,
976  meshEdgeToPatch,
977  maxEqOp<label>(),
978  label(-1)
979  );
981  (
982  mesh,
983  meshEdgeToZone,
984  maxEqOp<label>(),
985  label(-1)
986  );
987 // // Note: flipMap not yet done. Needs edge orientation. This is handled
988 // // later on.
989 // if (Pstream::parRun())
990 // {
991 // const globalMeshData& gd = mesh.globalData();
992 // const indirectPrimitivePatch& cpp = gd.coupledPatch();
993 //
994 // labelList patchEdges;
995 // labelList coupledEdges;
996 // bitSet sameEdgeOrientation;
997 // PatchTools::matchEdges
998 // (
999 // pp,
1000 // cpp,
1001 // patchEdges,
1002 // coupledEdges,
1003 // sameEdgeOrientation
1004 // );
1005 //
1006 // // Convert data on pp edges to data on coupled patch
1007 // labelList cppEdgeZoneID(cpp.nEdges(), -1);
1008 // boolList cppEdgeFlip(cpp.nEdges(), false);
1009 // forAll(coupledEdges, i)
1010 // {
1011 // label cppEdgei = coupledEdges[i];
1012 // label ppEdgei = patchEdges[i];
1013 //
1014 // cppEdgeZoneID[cppEdgei] = edgeZoneID[ppEdgei];
1015 // if (sameEdgeOrientation[i])
1016 // {
1017 // cppEdgeFlip[cppEdgei] = edgeFlip[ppEdgei];
1018 // }
1019 // else
1020 // {
1021 // cppEdgeFlip[cppEdgei] = !edgeFlip[ppEdgei];
1022 // }
1023 // }
1024 //
1025 // // Sync
1026 // const globalIndexAndTransform& git = gd.globalTransforms();
1027 // const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
1028 //
1029 // globalMeshData::syncData
1030 // (
1031 // cppEdgeZoneID,
1032 // gd.globalEdgeSlaves(),
1033 // gd.globalEdgeTransformedSlaves(),
1034 // edgeMap,
1035 // git,
1036 // maxEqOp<label>(),
1037 // dummyTransform()
1038 // );
1039 // globalMeshData::syncData
1040 // (
1041 // cppEdgeFlip,
1042 // gd.globalEdgeSlaves(),
1043 // gd.globalEdgeTransformedSlaves(),
1044 // edgeMap,
1045 // git,
1046 // andEqOp<bool>(),
1047 // dummyTransform()
1048 // );
1049 //
1050 // // Convert data on coupled edges to pp edges
1051 // forAll(coupledEdges, i)
1052 // {
1053 // label cppEdgei = coupledEdges[i];
1054 // label ppEdgei = patchEdges[i];
1055 //
1056 // edgeZoneID[ppEdgei] = cppEdgeZoneID[cppEdgei];
1057 // if (sameEdgeOrientation[i])
1058 // {
1059 // edgeFlip[ppEdgei] = cppEdgeFlip[cppEdgei];
1060 // }
1061 // else
1062 // {
1063 // edgeFlip[ppEdgei] = !cppEdgeFlip[cppEdgei];
1064 // }
1065 // }
1066 // }
1067 
1069  (
1070  mesh,
1071  meshEdgeToFlip,
1072  orEqOp<unsigned int>(),
1073  0
1074  );
1075 
1076  // Extract pp info
1077  patchEdgeToFace = UIndirectList<label>(meshEdgeToFace, meshEdges);
1078  patchEdgeToPatch = UIndirectList<label>(meshEdgeToPatch, meshEdges);
1079  patchEdgeToZone = UIndirectList<label>(meshEdgeToZone, meshEdges);
1080  patchEdgeToFlip.setSize(meshEdges.size());
1081  patchEdgeToFlip = false;
1082  forAll(meshEdges, i)
1083  {
1084  patchEdgeToFlip[i] = meshEdgeToFlip[meshEdges[i]];
1085  }
1086 }
1087 
1088 
1091  const bool zoneFromAnyFace,
1092 
1093  const polyMesh& mesh,
1094  const globalIndex& globalFaces,
1095  const labelListList& globalEdgeFaces,
1096  const indirectPrimitivePatch& pp,
1097 
1098  labelList& edgePatchID,
1099  label& nPatches,
1100  Map<label>& nbrProcToPatch,
1101  Map<label>& patchToNbrProc,
1102  labelList& edgeZoneID,
1103  boolList& edgeFlip,
1104  labelList& inflateFaceID
1105 )
1106 {
1108  const globalMeshData& gd = mesh.globalData();
1109 
1110  // Precalculate mesh edges for pp.edges.
1111  const labelList meshEdges(pp.meshEdges(mesh.edges(), mesh.pointEdges()));
1112 
1113  edgePatchID.setSize(pp.nEdges());
1114  edgePatchID = -1;
1115  nPatches = patches.size();
1116  edgeZoneID.setSize(pp.nEdges());
1117  edgeZoneID = -1;
1118  edgeFlip.setSize(pp.nEdges());
1119  edgeFlip = false;
1120  inflateFaceID.setSize(pp.nEdges(), -1);
1121 
1122 
1123  // Determine properties for faces extruded from edges
1124  // - edge inbetween two different processors:
1125  // - extrude as patch face on correct processor
1126  // - edge at end of patch (so edgeFaces.size() == 1):
1127  // - use mesh face that is a boundary face
1128  // - edge internal to patch (so edgeFaces.size() == 2):
1129 
1130 
1131  // Pass1:
1132  // For all edges inbetween two processors: see if matches to existing
1133  // processor patch or create interprocessor-patch if necessary.
1134  // Sets edgePatchID[edgeI] but none of the other quantities
1135 
1136  forAll(globalEdgeFaces, edgei)
1137  {
1138  const labelList& eGlobalFaces = globalEdgeFaces[edgei];
1139  if
1140  (
1141  eGlobalFaces.size() == 2
1142  && pp.edgeFaces()[edgei].size() == 1
1143  )
1144  {
1145  // Locally but not globally a boundary edge. Hence a coupled
1146  // edge. Find the patch to use if on different processors.
1147 
1148  label f0 = eGlobalFaces[0];
1149  label f1 = eGlobalFaces[1];
1150 
1151  label otherProci = -1;
1152  if (globalFaces.isLocal(f0) && !globalFaces.isLocal(f1))
1153  {
1154  otherProci = globalFaces.whichProcID(f1);
1155  }
1156  else if (!globalFaces.isLocal(f0) && globalFaces.isLocal(f1))
1157  {
1158  otherProci = globalFaces.whichProcID(f0);
1159  }
1160 
1161 
1162  if (otherProci != -1)
1163  {
1164  if (gd[Pstream::myProcNo()].found(otherProci))
1165  {
1166  // There is already a processorPolyPatch to otherProci.
1167  // Use it. Note that we can only index procPatchMap
1168  // if the processor actually is a neighbour processor
1169  // so that is why we first check.
1170  edgePatchID[edgei] = gd.procPatchMap()[otherProci];
1171  }
1172  else
1173  {
1174  // Cannot find a patch to processor. See if already
1175  // marked for addition
1176  if (nbrProcToPatch.found(otherProci))
1177  {
1178  edgePatchID[edgei] = nbrProcToPatch[otherProci];
1179  }
1180  else
1181  {
1182  edgePatchID[edgei] = nPatches;
1183  nbrProcToPatch.insert(otherProci, nPatches);
1184  patchToNbrProc.insert(nPatches, otherProci);
1185  nPatches++;
1186  }
1187  }
1188  }
1189  }
1190  }
1191 
1192 
1193  // Pass2: determine face properties for all other edges
1194  // ----------------------------------------------------
1195 
1196  // Info for edges of pp
1197  labelList edgeToFace;
1198  labelList edgeToPatch;
1199  labelList edgeToZone;
1200  bitSet edgeToFlip;
1201  globalEdgeInfo
1202  (
1203  zoneFromAnyFace, // internal edge info also from boundary faces
1204 
1205  mesh,
1206  globalFaces,
1207  globalEdgeFaces,
1208  pp,
1209  meshEdges,
1210 
1211  edgeToFace, // face (in globalFaces index)
1212  edgeToPatch, // patch on face (or -1 for internal faces)
1213  edgeToZone, // zone on face
1214  edgeToFlip // flip orientation on face
1215  );
1216 
1217  const labelListList& edgeFaces = pp.edgeFaces();
1218 
1219  DynamicList<label> dynMeshEdgeFaces;
1220 
1221  forAll(edgeFaces, edgei)
1222  {
1223  if (edgePatchID[edgei] == -1)
1224  {
1225  if (edgeFaces[edgei].size() == 2)
1226  {
1227  // Internal edge. Look at any face (internal or boundary) to
1228  // determine extrusion properties. First one that has zone
1229  // info wins
1230  if (globalFaces.isLocal(edgeToFace[edgei]))
1231  {
1232  inflateFaceID[edgei] =
1233  globalFaces.toLocal(edgeToFace[edgei]);
1234  }
1235  edgeZoneID[edgei] = edgeToZone[edgei];
1236  edgeFlip[edgei] = edgeToFlip[edgei];
1237  }
1238  else
1239  {
1240  // Proper, uncoupled patch edge. Boundary to get info from
1241  // might be on a different processor!
1242 
1243  if (globalFaces.isLocal(edgeToFace[edgei]))
1244  {
1245  inflateFaceID[edgei] =
1246  globalFaces.toLocal(edgeToFace[edgei]);
1247  }
1248  edgePatchID[edgei] = edgeToPatch[edgei];
1249  edgeZoneID[edgei] = edgeToZone[edgei];
1250  edgeFlip[edgei] = edgeToFlip[edgei];
1251  }
1252  }
1253  }
1254 
1255 
1256 
1257  // Now hopefully every boundary edge has a edge patch. Check
1258  if (debug)
1259  {
1260  forAll(edgeFaces, edgei)
1261  {
1262  if (edgeFaces[edgei].size() == 1 && edgePatchID[edgei] == -1)
1263  {
1264  const edge& e = pp.edges()[edgei];
1266  << "Have no sidePatchID for edge " << edgei << " points "
1267  << pp.points()[pp.meshPoints()[e[0]]]
1268  << pp.points()[pp.meshPoints()[e[1]]]
1269  << endl;
1270  }
1271  }
1272  }
1273 
1274 
1275 
1276  // Pass3: for any faces set in pass1 see if we can find a processor face
1277  // to inherit from (we only have a patch, not a patch face)
1278  forAll(edgeFaces, edgei)
1279  {
1280  if
1281  (
1282  edgeFaces[edgei].size() == 1
1283  && globalEdgeFaces[edgei].size() == 2
1284  && edgePatchID[edgei] != -1
1285  && inflateFaceID[edgei] == -1
1286  )
1287  {
1288  // 1. Do we have a local boundary face to inflate from
1289 
1290  label myFaceI = pp.addressing()[edgeFaces[edgei][0]];
1291 
1292  // Pick up any boundary face on this edge and use its properties
1293  label meshEdgei = meshEdges[edgei];
1294  const labelList& meshFaces = mesh.edgeFaces
1295  (
1296  meshEdgei,
1297  dynMeshEdgeFaces
1298  );
1299 
1300  forAll(meshFaces, k)
1301  {
1302  label facei = meshFaces[k];
1303 
1304  if (facei != myFaceI && !mesh.isInternalFace(facei))
1305  {
1306  if (patches.whichPatch(facei) == edgePatchID[edgei])
1307  {
1308  setFaceProps
1309  (
1310  mesh,
1311  pp,
1312  edgei,
1313  facei,
1314 
1315  edgePatchID[edgei],
1316  edgeZoneID[edgei],
1317  edgeFlip[edgei],
1318  inflateFaceID[edgei]
1319  );
1320  break;
1321  }
1322  }
1323  }
1324  }
1325  }
1326 
1327 
1328 
1329  // Sync all data:
1330  // - edgePatchID : might be local in case of processor patch. Do not
1331  // sync for now
1332  // - inflateFaceID: local. Do not sync
1333  // - edgeZoneID : global. sync.
1334  // - edgeFlip : global. sync.
1335 
1336  if (Pstream::parRun())
1337  {
1338  const globalMeshData& gd = mesh.globalData();
1339  const indirectPrimitivePatch& cpp = gd.coupledPatch();
1340 
1341  labelList patchEdges;
1342  labelList coupledEdges;
1343  bitSet sameEdgeOrientation;
1345  (
1346  pp,
1347  cpp,
1348  patchEdges,
1349  coupledEdges,
1350  sameEdgeOrientation
1351  );
1352 
1353  // Convert data on pp edges to data on coupled patch
1354  labelList cppEdgeZoneID(cpp.nEdges(), -1);
1355  boolList cppEdgeFlip(cpp.nEdges(), false);
1356  forAll(coupledEdges, i)
1357  {
1358  label cppEdgei = coupledEdges[i];
1359  label ppEdgei = patchEdges[i];
1360 
1361  cppEdgeZoneID[cppEdgei] = edgeZoneID[ppEdgei];
1362  if (sameEdgeOrientation[i])
1363  {
1364  cppEdgeFlip[cppEdgei] = edgeFlip[ppEdgei];
1365  }
1366  else
1367  {
1368  cppEdgeFlip[cppEdgei] = !edgeFlip[ppEdgei];
1369  }
1370  }
1371 
1372  // Sync
1373  const globalIndexAndTransform& git = gd.globalTransforms();
1374  const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
1375 
1377  (
1378  cppEdgeZoneID,
1379  gd.globalEdgeSlaves(),
1381  edgeMap,
1382  git,
1383  maxEqOp<label>(),
1384  dummyTransform()
1385  );
1387  (
1388  cppEdgeFlip,
1389  gd.globalEdgeSlaves(),
1391  edgeMap,
1392  git,
1393  andEqOp<bool>(),
1394  dummyTransform()
1395  );
1396 
1397  // Convert data on coupled edges to pp edges
1398  forAll(coupledEdges, i)
1399  {
1400  label cppEdgei = coupledEdges[i];
1401  label ppEdgei = patchEdges[i];
1402 
1403  edgeZoneID[ppEdgei] = cppEdgeZoneID[cppEdgei];
1404  if (sameEdgeOrientation[i])
1405  {
1406  edgeFlip[ppEdgei] = cppEdgeFlip[cppEdgei];
1407  }
1408  else
1409  {
1410  edgeFlip[ppEdgei] = !cppEdgeFlip[cppEdgei];
1411  }
1412  }
1413  }
1414 }
1415 
1416 
1419  const globalIndex& globalFaces,
1420  const labelListList& globalEdgeFaces,
1421  const scalarField& expansionRatio,
1422  const indirectPrimitivePatch& pp,
1423 
1424  const labelList& edgePatchID,
1425  const labelList& edgeZoneID,
1426  const boolList& edgeFlip,
1427  const labelList& inflateFaceID,
1428 
1429  const labelList& exposedPatchID,
1430  const labelList& nFaceLayers,
1431  const labelList& nPointLayers,
1432  const vectorField& firstLayerDisp,
1433  polyTopoChange& meshMod
1434 )
1435 {
1436  if (debug)
1437  {
1438  Pout<< "addPatchCellLayer::setRefinement : Adding up to "
1439  << gMax(nPointLayers)
1440  << " layers of cells to indirectPrimitivePatch with "
1441  << pp.nPoints() << " points" << endl;
1442  }
1443 
1444  if
1445  (
1446  pp.nPoints() != firstLayerDisp.size()
1447  || pp.nPoints() != nPointLayers.size()
1448  || pp.size() != nFaceLayers.size()
1449  )
1450  {
1452  << "Size of new points is not same as number of points used by"
1453  << " the face subset" << endl
1454  << " patch.nPoints:" << pp.nPoints()
1455  << " displacement:" << firstLayerDisp.size()
1456  << " nPointLayers:" << nPointLayers.size() << nl
1457  << " patch.nFaces:" << pp.size()
1458  << " nFaceLayers:" << nFaceLayers.size()
1459  << abort(FatalError);
1460  }
1461 
1462  forAll(nPointLayers, i)
1463  {
1464  if (nPointLayers[i] < 0)
1465  {
1467  << "Illegal number of layers " << nPointLayers[i]
1468  << " at patch point " << i << abort(FatalError);
1469  }
1470  }
1471  forAll(nFaceLayers, i)
1472  {
1473  if (nFaceLayers[i] < 0)
1474  {
1476  << "Illegal number of layers " << nFaceLayers[i]
1477  << " at patch face " << i << abort(FatalError);
1478  }
1479  }
1480 
1481  forAll(globalEdgeFaces, edgei)
1482  {
1483  if (globalEdgeFaces[edgei].size() > 2)
1484  {
1485  const edge& e = pp.edges()[edgei];
1486 
1487  if (nPointLayers[e[0]] > 0 || nPointLayers[e[1]] > 0)
1488  {
1490  << "Trying to extrude edge "
1491  << e.line(pp.localPoints())
1492  << " which is non-manifold (has "
1493  << globalEdgeFaces[edgei].size()
1494  << " local or coupled faces using it)"
1495  << " of which "
1496  << pp.edgeFaces()[edgei].size()
1497  << " local"
1498  << abort(FatalError);
1499  }
1500  }
1501  }
1502 
1503 
1504  const labelList& meshPoints = pp.meshPoints();
1505 
1506  // Some storage for edge-face-addressing.
1507  DynamicList<label> ef;
1508 
1509  // Precalculate mesh edges for pp.edges.
1510  const labelList meshEdges(pp.meshEdges(mesh_.edges(), mesh_.pointEdges()));
1511 
1512  if (debug)
1513  {
1514  // Check synchronisation
1515  // ~~~~~~~~~~~~~~~~~~~~~
1516 
1517  {
1518  labelList n(mesh_.nPoints(), Zero);
1519  labelUIndList(n, meshPoints) = nPointLayers;
1520  syncTools::syncPointList(mesh_, n, maxEqOp<label>(), label(0));
1521 
1522  // Non-synced
1523  forAll(meshPoints, i)
1524  {
1525  label meshPointi = meshPoints[i];
1526 
1527  if (n[meshPointi] != nPointLayers[i])
1528  {
1530  << "At mesh point:" << meshPointi
1531  << " coordinate:" << mesh_.points()[meshPointi]
1532  << " specified nLayers:" << nPointLayers[i] << endl
1533  << "On coupled point a different nLayers:"
1534  << n[meshPointi] << " was specified."
1535  << abort(FatalError);
1536  }
1537  }
1538 
1539 
1540  // Check that nPointLayers equals the max layers of connected faces
1541  // (or 0). Anything else makes no sense.
1542  labelList nFromFace(mesh_.nPoints(), Zero);
1543  forAll(nFaceLayers, i)
1544  {
1545  const face& f = pp[i];
1546 
1547  forAll(f, fp)
1548  {
1549  label pointi = f[fp];
1550 
1551  nFromFace[pointi] = max(nFromFace[pointi], nFaceLayers[i]);
1552  }
1553  }
1555  (
1556  mesh_,
1557  nFromFace,
1558  maxEqOp<label>(),
1559  label(0)
1560  );
1561 
1562  forAll(nPointLayers, i)
1563  {
1564  label meshPointi = meshPoints[i];
1565 
1566  if
1567  (
1568  nPointLayers[i] > 0
1569  && nPointLayers[i] != nFromFace[meshPointi]
1570  )
1571  {
1573  << "At mesh point:" << meshPointi
1574  << " coordinate:" << mesh_.points()[meshPointi]
1575  << " specified nLayers:" << nPointLayers[i] << endl
1576  << "but the max nLayers of surrounding faces is:"
1577  << nFromFace[meshPointi]
1578  << abort(FatalError);
1579  }
1580  }
1581  }
1582 
1583  {
1584  pointField d(mesh_.nPoints(), vector::max);
1585  UIndirectList<point>(d, meshPoints) = firstLayerDisp;
1587  (
1588  mesh_,
1589  d,
1590  minEqOp<vector>(),
1591  vector::max
1592  );
1593 
1594  forAll(meshPoints, i)
1595  {
1596  label meshPointi = meshPoints[i];
1597 
1598  if (mag(d[meshPointi] - firstLayerDisp[i]) > SMALL)
1599  {
1601  << "At mesh point:" << meshPointi
1602  << " coordinate:" << mesh_.points()[meshPointi]
1603  << " specified displacement:" << firstLayerDisp[i]
1604  << endl
1605  << "On coupled point a different displacement:"
1606  << d[meshPointi] << " was specified."
1607  << abort(FatalError);
1608  }
1609  }
1610  }
1611 
1612  // Check that edges of pp (so ones that become boundary faces)
1613  // connect to only one boundary face. Guarantees uniqueness of
1614  // patch that they go into so if this is a coupled patch both
1615  // sides decide the same.
1616  // ~~~~~~~~~~~~~~~~~~~~~~
1617 
1618  for (label edgei = pp.nInternalEdges(); edgei < pp.nEdges(); edgei++)
1619  {
1620  const edge& e = pp.edges()[edgei];
1621 
1622  if (nPointLayers[e[0]] > 0 || nPointLayers[e[1]] > 0)
1623  {
1624  // Edge is to become a face
1625 
1626  const labelList& eFaces = pp.edgeFaces()[edgei];
1627 
1628  // First check: pp should be single connected.
1629  if (eFaces.size() != 1)
1630  {
1632  << "boundary-edge-to-be-extruded:"
1633  << pp.points()[meshPoints[e[0]]]
1634  << pp.points()[meshPoints[e[1]]]
1635  << " has more than two faces using it:" << eFaces
1636  << abort(FatalError);
1637  }
1638 
1639  label myFacei = pp.addressing()[eFaces[0]];
1640 
1641  label meshEdgei = meshEdges[edgei];
1642 
1643  // Mesh faces using edge
1644  const labelList& meshFaces = mesh_.edgeFaces(meshEdgei, ef);
1645 
1646  // Check that there is only one patchface using edge.
1647  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1648 
1649  label bFacei = -1;
1650 
1651  forAll(meshFaces, i)
1652  {
1653  label facei = meshFaces[i];
1654 
1655  if (facei != myFacei)
1656  {
1657  if (!mesh_.isInternalFace(facei))
1658  {
1659  if (bFacei == -1)
1660  {
1661  bFacei = facei;
1662  }
1663  else
1664  {
1666  << "boundary-edge-to-be-extruded:"
1667  << pp.points()[meshPoints[e[0]]]
1668  << pp.points()[meshPoints[e[1]]]
1669  << " has more than two boundary faces"
1670  << " using it:"
1671  << bFacei << " fc:"
1672  << mesh_.faceCentres()[bFacei]
1673  << " patch:" << patches.whichPatch(bFacei)
1674  << " and " << facei << " fc:"
1675  << mesh_.faceCentres()[facei]
1676  << " patch:" << patches.whichPatch(facei)
1677  << abort(FatalError);
1678  }
1679  }
1680  }
1681  }
1682  }
1683  }
1684  }
1685 
1686 
1687  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
1688 
1689  // Precalculated patchID for each patch face
1690  labelList patchID(pp.size());
1691 
1692  forAll(pp, patchFacei)
1693  {
1694  label meshFacei = pp.addressing()[patchFacei];
1695 
1696  patchID[patchFacei] = patches.whichPatch(meshFacei);
1697  }
1698 
1699 
1700  // From master point (in patch point label) to added points (in mesh point
1701  // label)
1702  addedPoints_.setSize(pp.nPoints());
1703 
1704  // Mark points that do not get extruded by setting size of addedPoints_ to 0
1705  label nTruncated = 0;
1706 
1707  forAll(nPointLayers, patchPointi)
1708  {
1709  if (nPointLayers[patchPointi] > 0)
1710  {
1711  addedPoints_[patchPointi].setSize(nPointLayers[patchPointi]);
1712  }
1713  else
1714  {
1715  nTruncated++;
1716  }
1717  }
1718 
1719  if (debug)
1720  {
1721  Pout<< "Not adding points at " << nTruncated << " out of "
1722  << pp.nPoints() << " points" << endl;
1723  }
1724 
1725 
1726  //
1727  // Create new points
1728  //
1729 
1730  // If creating new mesh: copy existing patch points
1731  labelList copiedPatchPoints;
1732  if (!addToMesh_)
1733  {
1734  copiedPatchPoints.setSize(firstLayerDisp.size());
1735  forAll(firstLayerDisp, patchPointi)
1736  {
1737  label meshPointi = meshPoints[patchPointi];
1738  label zoneI = mesh_.pointZones().whichZone(meshPointi);
1739  copiedPatchPoints[patchPointi] = meshMod.setAction
1740  (
1741  polyAddPoint
1742  (
1743  mesh_.points()[meshPointi], // point
1744  -1, // master point
1745  zoneI, // zone for point
1746  true // supports a cell
1747  )
1748  );
1749  }
1750  }
1751 
1752 
1753  // Create points for additional layers
1754  forAll(firstLayerDisp, patchPointi)
1755  {
1756  if (addedPoints_[patchPointi].size())
1757  {
1758  label meshPointi = meshPoints[patchPointi];
1759 
1760  label zoneI = mesh_.pointZones().whichZone(meshPointi);
1761 
1762  point pt = mesh_.points()[meshPointi];
1763 
1764  vector disp = firstLayerDisp[patchPointi];
1765 
1766  forAll(addedPoints_[patchPointi], i)
1767  {
1768  pt += disp;
1769 
1770  label addedVertI = meshMod.setAction
1771  (
1772  polyAddPoint
1773  (
1774  pt, // point
1775  (addToMesh_ ? meshPointi : -1), // master point
1776  zoneI, // zone for point
1777  true // supports a cell
1778  )
1779  );
1780 
1781  addedPoints_[patchPointi][i] = addedVertI;
1782 
1783  disp *= expansionRatio[patchPointi];
1784  }
1785  }
1786  }
1787 
1788 
1789  //
1790  // Add cells to all boundaryFaces
1791  //
1792 
1793  labelListList addedCells(pp.size());
1794 
1795  forAll(pp, patchFacei)
1796  {
1797  if (nFaceLayers[patchFacei] > 0)
1798  {
1799  addedCells[patchFacei].setSize(nFaceLayers[patchFacei]);
1800 
1801  label meshFacei = pp.addressing()[patchFacei];
1802 
1803  label ownZoneI = mesh_.cellZones().whichZone
1804  (
1805  mesh_.faceOwner()[meshFacei]
1806  );
1807 
1808  for (label i = 0; i < nFaceLayers[patchFacei]; i++)
1809  {
1810  // Note: add from cell (owner of patch face) or from face?
1811  // for now add from cell so we can map easily.
1812  addedCells[patchFacei][i] = meshMod.setAction
1813  (
1814  polyAddCell
1815  (
1816  -1, // master point
1817  -1, // master edge
1818  -1, // master face
1819  (addToMesh_ ? mesh_.faceOwner()[meshFacei] : -1),
1820  //master
1821  ownZoneI // zone for cell
1822  )
1823  );
1824  }
1825  }
1826  }
1827 
1828 
1829 
1830  // Create faces on top of the original patch faces.
1831  // These faces are created from original patch faces outwards so in order
1832  // of increasing cell number. So orientation should be same as original
1833  // patch face for them to have owner<neighbour.
1834 
1835  layerFaces_.setSize(pp.size());
1836 
1837  forAll(pp.localFaces(), patchFacei)
1838  {
1839  label meshFacei = pp.addressing()[patchFacei];
1840 
1841  if (addedCells[patchFacei].size())
1842  {
1843  layerFaces_[patchFacei].setSize(addedCells[patchFacei].size() + 1);
1844 
1845  // Get duplicated vertices on the patch face.
1846  const face& f = pp.localFaces()[patchFacei];
1847 
1848  face newFace(f.size());
1849 
1850  forAll(addedCells[patchFacei], i)
1851  {
1852  forAll(f, fp)
1853  {
1854  if (addedPoints_[f[fp]].empty())
1855  {
1856  // Keep original point
1857  newFace[fp] =
1858  (
1859  addToMesh_
1860  ? meshPoints[f[fp]]
1861  : copiedPatchPoints[f[fp]]
1862  );
1863  }
1864  else
1865  {
1866  // Get new outside point
1867  label offset =
1868  addedPoints_[f[fp]].size()
1869  - addedCells[patchFacei].size();
1870  newFace[fp] = addedPoints_[f[fp]][i+offset];
1871  }
1872  }
1873 
1874 
1875  // Get new neighbour
1876  label nei;
1877  label patchi;
1878  label zoneI = -1;
1879  bool flip = false;
1880 
1881 
1882  if (i == addedCells[patchFacei].size()-1)
1883  {
1884  // Top layer so is patch face.
1885  nei = -1;
1886  patchi = patchID[patchFacei];
1887  zoneI = mesh_.faceZones().whichZone(meshFacei);
1888  if (zoneI != -1)
1889  {
1890  const faceZone& fz = mesh_.faceZones()[zoneI];
1891  flip = fz.flipMap()[fz.whichFace(meshFacei)];
1892  }
1893  }
1894  else
1895  {
1896  // Internal face between layer i and i+1
1897  nei = addedCells[patchFacei][i+1];
1898  patchi = -1;
1899  }
1900 
1901 
1902  layerFaces_[patchFacei][i+1] = meshMod.setAction
1903  (
1904  polyAddFace
1905  (
1906  newFace, // face
1907  addedCells[patchFacei][i], // owner
1908  nei, // neighbour
1909  -1, // master point
1910  -1, // master edge
1911  (addToMesh_ ? meshFacei : -1), // master face
1912  false, // flux flip
1913  patchi, // patch for face
1914  zoneI, // zone for face
1915  flip // face zone flip
1916  )
1917  );
1918  }
1919  }
1920  }
1921 
1922  //
1923  // Modify old patch faces to be on the inside
1924  //
1925 
1926  if (addToMesh_)
1927  {
1928  forAll(pp, patchFacei)
1929  {
1930  if (addedCells[patchFacei].size())
1931  {
1932  label meshFacei = pp.addressing()[patchFacei];
1933 
1934  layerFaces_[patchFacei][0] = meshFacei;
1935 
1936  meshMod.setAction
1937  (
1939  (
1940  pp[patchFacei], // modified face
1941  meshFacei, // label of face
1942  mesh_.faceOwner()[meshFacei], // owner
1943  addedCells[patchFacei][0], // neighbour
1944  false, // face flip
1945  -1, // patch for face
1946  true, //false, // remove from zone
1947  -1, //zoneI, // zone for face
1948  false // face flip in zone
1949  )
1950  );
1951  }
1952  }
1953  }
1954  else
1955  {
1956  // If creating new mesh: reverse original faces and put them
1957  // in the exposed patch ID.
1958  forAll(pp, patchFacei)
1959  {
1960  if (nFaceLayers[patchFacei] > 0)
1961  {
1962  label meshFacei = pp.addressing()[patchFacei];
1963  label zoneI = mesh_.faceZones().whichZone(meshFacei);
1964  bool zoneFlip = false;
1965  if (zoneI != -1)
1966  {
1967  const faceZone& fz = mesh_.faceZones()[zoneI];
1968  zoneFlip = !fz.flipMap()[fz.whichFace(meshFacei)];
1969  }
1970 
1971  // Reverse and renumber old patch face.
1972  face f(pp.localFaces()[patchFacei].reverseFace());
1973  forAll(f, fp)
1974  {
1975  f[fp] = copiedPatchPoints[f[fp]];
1976  }
1977 
1978  layerFaces_[patchFacei][0] = meshMod.setAction
1979  (
1980  polyAddFace
1981  (
1982  f, // modified face
1983  addedCells[patchFacei][0], // owner
1984  -1, // neighbour
1985  -1, // masterPoint
1986  -1, // masterEdge
1987  -1, // masterFace
1988  true, // face flip
1989  exposedPatchID[patchFacei], // patch for face
1990  zoneI, // zone for face
1991  zoneFlip // face flip in zone
1992  )
1993  );
1994  }
1995  }
1996  }
1997 
1998 
1999 
2000  //
2001  // Create 'side' faces, one per edge that is being extended.
2002  //
2003 
2004  const labelListList& faceEdges = pp.faceEdges();
2005  const faceList& localFaces = pp.localFaces();
2006  const edgeList& edges = pp.edges();
2007 
2008  // Get number of layers per edge. This is 0 if edge is not extruded;
2009  // max of connected faces otherwise.
2010  labelList edgeLayers(pp.nEdges());
2011 
2012  {
2013  // Use list over mesh.nEdges() since syncTools does not yet support
2014  // partial list synchronisation.
2015  labelList meshEdgeLayers(mesh_.nEdges(), -1);
2016 
2017  forAll(meshEdges, edgei)
2018  {
2019  const edge& e = edges[edgei];
2020 
2021  label meshEdgei = meshEdges[edgei];
2022 
2023  if ((nPointLayers[e[0]] == 0) && (nPointLayers[e[1]] == 0))
2024  {
2025  meshEdgeLayers[meshEdgei] = 0;
2026  }
2027  else
2028  {
2029  const labelList& eFaces = pp.edgeFaces()[edgei];
2030 
2031  forAll(eFaces, i)
2032  {
2033  meshEdgeLayers[meshEdgei] = max
2034  (
2035  nFaceLayers[eFaces[i]],
2036  meshEdgeLayers[meshEdgei]
2037  );
2038  }
2039  }
2040  }
2041 
2043  (
2044  mesh_,
2045  meshEdgeLayers,
2046  maxEqOp<label>(),
2047  label(0) // initial value
2048  );
2049 
2050  forAll(meshEdges, edgei)
2051  {
2052  edgeLayers[edgei] = meshEdgeLayers[meshEdges[edgei]];
2053  }
2054  }
2055 
2056 
2057  // Mark off which edges have been extruded
2058  boolList doneEdge(pp.nEdges(), false);
2059 
2060 
2061  // Create faces. Per face walk connected edges and find string of edges
2062  // between the same two faces and extrude string into a single face.
2063  forAll(pp, patchFacei)
2064  {
2065  const labelList& fEdges = faceEdges[patchFacei];
2066 
2067  forAll(fEdges, fp)
2068  {
2069  // Get string of edges that needs to be extruded as a single face.
2070  // Returned as indices in fEdges.
2071  labelPair indexPair
2072  (
2073  getEdgeString
2074  (
2075  pp,
2076  globalEdgeFaces,
2077  doneEdge,
2078  patchFacei,
2079  globalFaces.toGlobal(pp.addressing()[patchFacei])
2080  )
2081  );
2082 
2083  //Pout<< "Found unextruded edges in edges:" << fEdges
2084  // << " start:" << indexPair[0]
2085  // << " end:" << indexPair[1]
2086  // << endl;
2087 
2088  const label startFp = indexPair[0];
2089  const label endFp = indexPair[1];
2090 
2091  if (startFp != -1)
2092  {
2093  // Extrude edges from indexPair[0] up to indexPair[1]
2094  // (note indexPair = indices of edges. There is one more vertex
2095  // than edges)
2096  const face& f = localFaces[patchFacei];
2097 
2098  labelList stringedVerts;
2099  if (endFp >= startFp)
2100  {
2101  stringedVerts.setSize(endFp-startFp+2);
2102  }
2103  else
2104  {
2105  stringedVerts.setSize(endFp+f.size()-startFp+2);
2106  }
2107 
2108  label fp = startFp;
2109 
2110  for (label i = 0; i < stringedVerts.size()-1; i++)
2111  {
2112  stringedVerts[i] = f[fp];
2113  doneEdge[fEdges[fp]] = true;
2114  fp = f.fcIndex(fp);
2115  }
2116  stringedVerts.last() = f[fp];
2117 
2118 
2119  // Now stringedVerts contains the vertices in order of face f.
2120  // This is consistent with the order if f becomes the owner cell
2121  // and nbrFacei the neighbour cell. Note that the cells get
2122  // added in order of pp so we can just use face ordering and
2123  // because we loop in incrementing order as well we will
2124  // always have nbrFacei > patchFacei.
2125 
2126  label startEdgei = fEdges[startFp];
2127 
2128  label meshEdgei = meshEdges[startEdgei];
2129 
2130  label numEdgeSideFaces = edgeLayers[startEdgei];
2131 
2132  for (label i = 0; i < numEdgeSideFaces; i++)
2133  {
2134  label vEnd = stringedVerts.last();
2135  label vStart = stringedVerts[0];
2136 
2137  // calculate number of points making up a face
2138  label newFp = 2*stringedVerts.size();
2139 
2140  if (i == 0)
2141  {
2142  // layer 0 gets all the truncation of neighbouring
2143  // faces with more layers.
2144  if (addedPoints_[vEnd].size())
2145  {
2146  newFp +=
2147  addedPoints_[vEnd].size() - numEdgeSideFaces;
2148  }
2149  if (addedPoints_[vStart].size())
2150  {
2151  newFp +=
2152  addedPoints_[vStart].size() - numEdgeSideFaces;
2153  }
2154  }
2155 
2156  face newFace(newFp);
2157 
2158  newFp = 0;
2159 
2160  // For layer 0 get pp points, for all other layers get
2161  // points of layer-1.
2162  if (i == 0)
2163  {
2164  forAll(stringedVerts, stringedI)
2165  {
2166  label v = stringedVerts[stringedI];
2167  addVertex
2168  (
2169  (
2170  addToMesh_
2171  ? meshPoints[v]
2172  : copiedPatchPoints[v]
2173  ),
2174  newFace,
2175  newFp
2176  );
2177  }
2178  }
2179  else
2180  {
2181  forAll(stringedVerts, stringedI)
2182  {
2183  label v = stringedVerts[stringedI];
2184  if (addedPoints_[v].size())
2185  {
2186  label offset =
2187  addedPoints_[v].size() - numEdgeSideFaces;
2188  addVertex
2189  (
2190  addedPoints_[v][i+offset-1],
2191  newFace,
2192  newFp
2193  );
2194  }
2195  else
2196  {
2197  addVertex
2198  (
2199  (
2200  addToMesh_
2201  ? meshPoints[v]
2202  : copiedPatchPoints[v]
2203  ),
2204  newFace,
2205  newFp
2206  );
2207  }
2208  }
2209  }
2210 
2211  // add points between stringed vertices (end)
2212  if (numEdgeSideFaces < addedPoints_[vEnd].size())
2213  {
2214  if (i == 0 && addedPoints_[vEnd].size())
2215  {
2216  label offset =
2217  addedPoints_[vEnd].size() - numEdgeSideFaces;
2218  for (label ioff = 0; ioff < offset; ioff++)
2219  {
2220  addVertex
2221  (
2222  addedPoints_[vEnd][ioff],
2223  newFace,
2224  newFp
2225  );
2226  }
2227  }
2228  }
2229 
2230  forAllReverse(stringedVerts, stringedI)
2231  {
2232  label v = stringedVerts[stringedI];
2233  if (addedPoints_[v].size())
2234  {
2235  label offset =
2236  addedPoints_[v].size() - numEdgeSideFaces;
2237  addVertex
2238  (
2239  addedPoints_[v][i+offset],
2240  newFace,
2241  newFp
2242  );
2243  }
2244  else
2245  {
2246  addVertex
2247  (
2248  (
2249  addToMesh_
2250  ? meshPoints[v]
2251  : copiedPatchPoints[v]
2252  ),
2253  newFace,
2254  newFp
2255  );
2256  }
2257  }
2258 
2259 
2260  // add points between stringed vertices (start)
2261  if (numEdgeSideFaces < addedPoints_[vStart].size())
2262  {
2263  if (i == 0 && addedPoints_[vStart].size())
2264  {
2265  label offset =
2266  addedPoints_[vStart].size() - numEdgeSideFaces;
2267  for (label ioff = offset-1; ioff >= 0; ioff--)
2268  {
2269  addVertex
2270  (
2271  addedPoints_[vStart][ioff],
2272  newFace,
2273  newFp
2274  );
2275  }
2276  }
2277  }
2278 
2279  if (newFp >= 3)
2280  {
2281  // Add face inbetween faces patchFacei and nbrFacei
2282  // (possibly -1 for external edges)
2283 
2284  newFace.setSize(newFp);
2285 
2286  if (debug)
2287  {
2288  labelHashSet verts(2*newFace.size());
2289  forAll(newFace, fp)
2290  {
2291  if (!verts.insert(newFace[fp]))
2292  {
2294  << "Duplicate vertex in face"
2295  << " to be added." << nl
2296  << "newFace:" << newFace << nl
2297  << "points:"
2299  (
2300  meshMod.points(),
2301  newFace
2302  ) << nl
2303  << "Layer:" << i
2304  << " out of:" << numEdgeSideFaces << nl
2305  << "ExtrudeEdge:" << meshEdgei
2306  << " at:"
2307  << mesh_.edges()[meshEdgei].line
2308  (
2309  mesh_.points()
2310  ) << nl
2311  << "string:" << stringedVerts
2312  << "stringpoints:"
2314  (
2315  pp.localPoints(),
2316  stringedVerts
2317  ) << nl
2318  << "stringNLayers:"
2319  << labelUIndList
2320  (
2321  nPointLayers,
2322  stringedVerts
2323  ) << nl
2324  << abort(FatalError);
2325  }
2326  }
2327  }
2328 
2329  label nbrFacei = nbrFace
2330  (
2331  pp.edgeFaces(),
2332  startEdgei,
2333  patchFacei
2334  );
2335 
2336  const labelList& meshFaces = mesh_.edgeFaces
2337  (
2338  meshEdgei,
2339  ef
2340  );
2341 
2342  // Because we walk in order of patch face and in order
2343  // of face edges so face orientation will be opposite
2344  // that of the patch edge
2345  bool zoneFlip = false;
2346  if (edgeZoneID[startEdgei] != -1)
2347  {
2348  zoneFlip = !edgeFlip[startEdgei];
2349  }
2350 
2351  addSideFace
2352  (
2353  pp,
2354  addedCells,
2355 
2356  newFace, // vertices of new face
2357  edgePatchID[startEdgei],// -1 or patch for face
2358  edgeZoneID[startEdgei],
2359  zoneFlip,
2360  inflateFaceID[startEdgei],
2361 
2362  patchFacei,
2363  nbrFacei,
2364  meshEdgei, // (mesh) edge to inflate
2365  i, // layer
2366  numEdgeSideFaces, // num layers
2367  meshFaces, // edgeFaces
2368  meshMod
2369  );
2370  }
2371  }
2372  }
2373  }
2374  }
2375 }
2376 
2377 
2380  const mapPolyMesh& morphMap,
2381  const labelList& faceMap, // new to old patch faces
2382  const labelList& pointMap // new to old patch points
2383 )
2384 {
2385  {
2386  labelListList newAddedPoints(pointMap.size());
2387 
2388  forAll(newAddedPoints, newPointi)
2389  {
2390  label oldPointi = pointMap[newPointi];
2391 
2392  const labelList& added = addedPoints_[oldPointi];
2393 
2394  labelList& newAdded = newAddedPoints[newPointi];
2395  newAdded.setSize(added.size());
2396  label newI = 0;
2397 
2398  forAll(added, i)
2399  {
2400  label newPointi = morphMap.reversePointMap()[added[i]];
2401 
2402  if (newPointi >= 0)
2403  {
2404  newAdded[newI++] = newPointi;
2405  }
2406  }
2407  newAdded.setSize(newI);
2408  }
2409  addedPoints_.transfer(newAddedPoints);
2410  }
2411 
2412  {
2413  labelListList newLayerFaces(faceMap.size());
2414 
2415  forAll(newLayerFaces, newFacei)
2416  {
2417  label oldFacei = faceMap[newFacei];
2418 
2419  const labelList& added = layerFaces_[oldFacei];
2420 
2421  labelList& newAdded = newLayerFaces[newFacei];
2422  newAdded.setSize(added.size());
2423  label newI = 0;
2424 
2425  forAll(added, i)
2426  {
2427  label newFacei = morphMap.reverseFaceMap()[added[i]];
2428 
2429  if (newFacei >= 0)
2430  {
2431  newAdded[newI++] = newFacei;
2432  }
2433  }
2434  newAdded.setSize(newI);
2435  }
2436  layerFaces_.transfer(newLayerFaces);
2437  }
2438 }
2439 
2440 
2441 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::PrimitivePatch::points
const Field< point_type > & points() const noexcept
Return reference to global points.
Definition: PrimitivePatch.H:299
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::syncTools::syncEdgeList
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.
Definition: syncToolsTemplates.C:800
addPatchCellLayer.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
meshTools.H
Foam::syncTools::getInternalOrCoupledFaces
static bitSet getInternalOrCoupledFaces(const polyMesh &mesh)
Get per face whether it is internal or coupled.
Definition: syncTools.C:176
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeTopological.C:94
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatch.C:262
Foam::faceZone::flipMap
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:272
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:183
Foam::globalMeshData::globalEdgeSlavesMap
const mapDistribute & globalEdgeSlavesMap() const
Definition: globalMeshData.C:2245
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
PatchTools.H
Foam::addPatchCellLayer::globalEdgeFaces
static labelListList globalEdgeFaces(const polyMesh &, const globalIndex &globalFaces, const indirectPrimitivePatch &pp)
Per patch edge the pp faces (in global indices) using it.
Definition: addPatchCellLayer.C:647
Foam::DynamicList< label >
Foam::polyTopoChange::points
const DynamicList< point > & points() const
Points. Shrunk after constructing mesh (or calling of compact())
Definition: polyTopoChange.H:450
polyAddFace.H
Foam::primitiveMesh::edgeFaces
const labelListList & edgeFaces() const
Definition: primitiveMeshEdgeFaces.C:34
mapPolyMesh.H
globalIndex.H
Foam::PrimitivePatch::nEdges
label nEdges() const
Number of edges in patch.
Definition: PrimitivePatch.H:322
polyTopoChange.H
Foam::andEqOp
Definition: ops.H:85
Foam::PrimitivePatch::meshEdges
labelList meshEdges(const edgeList &allEdges, const labelListList &cellEdges, const labelList &faceCells) const
Definition: PrimitivePatchMeshEdges.C:53
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:99
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::primitiveMesh::nEdges
label nEdges() const
Number of mesh edges.
Definition: primitiveMeshI.H:67
dummyTransform.H
Dummy transform to be used with syncTools.
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:574
Foam::Map< label >
Foam::orEqOp
Definition: ops.H:86
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
Foam::globalIndex::isLocal
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:224
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:50
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::PtrList::set
const T * set(const label i) const
Return const pointer to element (can be nullptr),.
Definition: PtrList.H:138
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
Foam::dummyTransform
Definition: dummyTransform.H:47
Foam::PrimitivePatch::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatch.C:275
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::primitiveMesh::pointEdges
const labelListList & pointEdges() const
Definition: primitiveMeshEdges.C:516
polyMesh.H
Foam::HashSet< label, Hash< label > >
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:721
Foam::globalMeshData::processorPatches
const labelList & processorPatches() const noexcept
Return list of processor patch labels.
Definition: globalMeshData.H:381
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::minEqOp
Definition: ops.H:81
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
Foam::faceZoneMesh
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Definition: faceZoneMeshFwd.H:44
n
label n
Definition: TABSMDCalcMethod2.H:31
nPatches
const label nPatches
Definition: printMeshSummary.H:30
Foam::Field< scalar >
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::globalMeshData::globalTransforms
const globalIndexAndTransform & globalTransforms() const
Global transforms numbering.
Definition: globalMeshData.C:2160
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:486
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:812
Foam::globalMeshData::globalEdgeSlaves
const labelListList & globalEdgeSlaves() const
Definition: globalMeshData.C:2214
polyAddPoint.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::globalMeshData::syncData
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.
Definition: globalMeshDataTemplates.C:37
Foam::PatchTools::matchEdges
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.
Definition: PatchToolsMatch.C:77
Foam::PackedList::setSize
void setSize(const label n, unsigned int val=0u)
Alias for resize()
Definition: PackedList.H:558
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:289
newPointi
label newPointi
Definition: readKivaGrid.H:496
Foam::PrimitivePatch::nPoints
label nPoints() const
Number of points supporting patch faces.
Definition: PrimitivePatch.H:316
Foam::FatalError
error FatalError
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:107
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:317
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:372
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::polyAddCell
Class containing data for cell addition.
Definition: polyAddCell.H:48
Foam::ListOps::uniqueEqOp
List helper to append y unique elements onto the end of x.
Definition: ListOps.H:577
Foam::maxEqOp
Definition: ops.H:80
Foam::addPatchCellLayer::updateMesh
void updateMesh(const mapPolyMesh &, const labelList &faceMap, const labelList &pointMap)
Update any locally stored mesh information. Gets additional.
Definition: addPatchCellLayer.C:2379
polyModifyFace.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:359
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
Foam::indirectPrimitivePatch
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
Definition: indirectPrimitivePatch.H:49
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:214
found
bool found
Definition: TABSMDCalcMethod2.H:32
polyAddCell.H
Foam::globalIndex::toLocal
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:305
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::edge::compare
static int compare(const edge &a, const edge &b)
Compare edges.
Definition: edgeI.H:33
Foam::polyTopoChange::setAction
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
Definition: polyTopoChange.C:2481
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::max
static const Vector< scalar > max
Definition: VectorSpace.H:117
Foam::globalMeshData::globalEdgeTransformedSlaves
const labelListList & globalEdgeTransformedSlaves() const
Definition: globalMeshData.C:2224
f
labelList f(nPoints)
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< labelList >
Foam::addPatchCellLayer::calcExtrudeInfo
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:
Definition: addPatchCellLayer.C:1090
Foam::globalMeshData::coupledPatch
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Definition: globalMeshData.C:2046
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::globalIndex::whichProcID
label whichProcID(const label i) const
Which processor does global come from? Binary search.
Definition: globalIndexI.H:311
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::polyAddPoint
Class containing data for point addition.
Definition: polyAddPoint.H:49
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:309
Foam::polyAddFace
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:53
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:330
Foam::globalIndexAndTransform
Determination and storage of the possible independent transforms introduced by coupledPolyPatches,...
Definition: globalIndexAndTransform.H:64
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::ProcessorTopology::procPatchMap
const labelList & procPatchMap() const
From neighbour processor to index in boundaryMesh. Local information.
Definition: ProcessorTopology.H:94
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:240
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
Foam::addPatchCellLayer::setRefinement
void setRefinement(const globalIndex &globalFaces, const labelListList &globalEdgeFaces, const scalarField &expansionRatio, const indirectPrimitivePatch &pp, 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.
Definition: addPatchCellLayer.C:1418
Foam::labelUIndList
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: UIndirectList.H:58
Foam::addPatchCellLayer::addedCells
labelListList addedCells() const
Added cells given current mesh & layerfaces.
Definition: addPatchCellLayer.C:640