createShellMesh.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) 2019-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 "createShellMesh.H"
30 #include "polyTopoChange.H"
31 #include "meshTools.H"
32 #include "mapPolyMesh.H"
33 #include "polyAddPoint.H"
34 #include "polyAddFace.H"
35 #include "polyModifyFace.H"
36 #include "polyAddCell.H"
37 #include "labelPair.H"
38 #include "indirectPrimitivePatch.H"
39 #include "mapDistribute.H"
40 #include "globalMeshData.H"
41 #include "PatchTools.H"
42 #include "globalIndex.H"
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 defineTypeNameAndDebug(createShellMesh, 0);
49 
50 template<>
52 {
53 public:
54  void operator()(labelPair& x, const labelPair& y) const
55  {
56  x[0] = min(x[0], y[0]);
57  x[1] = min(x[1], y[1]);
58  }
59 };
60 }
61 
62 
63 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64 
65 // Synchronise edges
66 void Foam::createShellMesh::syncEdges
67 (
68  const globalMeshData& globalData,
69 
70  const labelList& patchEdges,
71  const labelList& coupledEdges,
72  const bitSet& sameEdgeOrientation,
73  const bool syncNonCollocated,
74 
75  bitSet& isChangedEdge,
76  DynamicList<label>& changedEdges,
77  labelPairList& allEdgeData
78 )
79 {
80  const mapDistribute& map = globalData.globalEdgeSlavesMap();
81  const bitSet& cppOrientation = globalData.globalEdgeOrientation();
82 
83  // Convert patch-edge data into cpp-edge data
84  labelPairList cppEdgeData
85  (
86  map.constructSize(),
88  );
89 
90  forAll(patchEdges, i)
91  {
92  label patchEdgeI = patchEdges[i];
93  label coupledEdgeI = coupledEdges[i];
94 
95  if (isChangedEdge[patchEdgeI])
96  {
97  const labelPair& data = allEdgeData[patchEdgeI];
98 
99  // Patch-edge data needs to be converted into coupled-edge data
100  // (optionally flipped) and consistent in orientation with
101  // other coupled edge (optionally flipped)
102  if (sameEdgeOrientation[i] == cppOrientation[coupledEdgeI])
103  {
104  cppEdgeData[coupledEdgeI] = data;
105  }
106  else
107  {
108  cppEdgeData[coupledEdgeI] = labelPair(data[1], data[0]);
109  }
110  }
111  }
112 
113  // Synchronise
114  globalData.syncData
115  (
116  cppEdgeData,
117  globalData.globalEdgeSlaves(),
118  (
119  syncNonCollocated
120  ? globalData.globalEdgeTransformedSlaves() // transformed elems
121  : labelListList(globalData.globalEdgeSlaves().size()) //no transformed
122  ),
123  map,
124  minEqOp<labelPair>()
125  );
126 
127  // Back from cpp-edge to patch-edge data
128  forAll(patchEdges, i)
129  {
130  label patchEdgeI = patchEdges[i];
131  label coupledEdgeI = coupledEdges[i];
132 
133  if (cppEdgeData[coupledEdgeI] != labelPair(labelMax, labelMax))
134  {
135  const labelPair& data = cppEdgeData[coupledEdgeI];
136 
137  if (sameEdgeOrientation[i] == cppOrientation[coupledEdgeI])
138  {
139  allEdgeData[patchEdgeI] = data;
140  }
141  else
142  {
143  allEdgeData[patchEdgeI] = labelPair(data[1], data[0]);
144  }
145 
146  if (!isChangedEdge[patchEdgeI])
147  {
148  changedEdges.append(patchEdgeI);
149  isChangedEdge.set(patchEdgeI);
150  }
151  }
152  }
153 }
154 
155 
157 (
158  const globalMeshData& globalData,
159  const primitiveFacePatch& patch,
160  const bitSet& nonManifoldEdge,
161  const bool syncNonCollocated,
162 
163  faceList& pointGlobalRegions,
164  faceList& pointLocalRegions,
165  labelList& localToGlobalRegion
166 )
167 {
168  const indirectPrimitivePatch& cpp = globalData.coupledPatch();
169 
170  // Calculate correspondence between patch and globalData.coupledPatch.
171  labelList patchEdges;
172  labelList coupledEdges;
173  bitSet sameEdgeOrientation;
175  (
176  cpp,
177  patch,
178 
179  coupledEdges,
180  patchEdges,
181  sameEdgeOrientation
182  );
183 
184 
185  // Initial unique regions
186  // ~~~~~~~~~~~~~~~~~~~~~~
187  // These get merged later on across connected edges.
188 
189  // 1. Count
190  label nMaxRegions = 0;
191  forAll(patch.localFaces(), facei)
192  {
193  const face& f = patch.localFaces()[facei];
194  nMaxRegions += f.size();
195  }
196 
197  const globalIndex globalRegions(nMaxRegions);
198 
199  // 2. Assign unique regions
200  label nRegions = 0;
201 
202  pointGlobalRegions.setSize(patch.size());
203  forAll(pointGlobalRegions, facei)
204  {
205  const face& f = patch.localFaces()[facei];
206  labelList& pRegions = pointGlobalRegions[facei];
207  pRegions.setSize(f.size());
208  forAll(pRegions, fp)
209  {
210  pRegions[fp] = globalRegions.toGlobal(nRegions++);
211  }
212  }
213 
214 
215  DynamicList<label> changedEdges(patch.nEdges());
216  labelPairList allEdgeData(patch.nEdges(), labelPair(labelMax, labelMax));
217  bitSet isChangedEdge(patch.nEdges());
218 
219 
220  // Fill initial seed
221  // ~~~~~~~~~~~~~~~~~
222 
223  forAll(patch.edgeFaces(), edgeI)
224  {
225  if (!nonManifoldEdge[edgeI])
226  {
227  // Take over value from one face only.
228  const edge& e = patch.edges()[edgeI];
229  label facei = patch.edgeFaces()[edgeI][0];
230  const face& f = patch.localFaces()[facei];
231 
232  label fp0 = f.find(e[0]);
233  label fp1 = f.find(e[1]);
234  allEdgeData[edgeI] = labelPair
235  (
236  pointGlobalRegions[facei][fp0],
237  pointGlobalRegions[facei][fp1]
238  );
239  if (!isChangedEdge[edgeI])
240  {
241  changedEdges.append(edgeI);
242  isChangedEdge.set(edgeI);
243  }
244  }
245  }
246 
247 
248  syncEdges
249  (
250  globalData,
251 
252  patchEdges,
253  coupledEdges,
254  sameEdgeOrientation,
255  syncNonCollocated,
256 
257  isChangedEdge,
258  changedEdges,
259  allEdgeData
260  );
261 
262 
263  // Edge-Face-Edge walk across patch
264  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
265  // Across edge minimum regions win
266 
267  while (true)
268  {
269  // From edge to face
270  // ~~~~~~~~~~~~~~~~~
271 
272  DynamicList<label> changedFaces(patch.size());
273  bitSet isChangedFace(patch.size());
274 
275  forAll(changedEdges, changedI)
276  {
277  label edgeI = changedEdges[changedI];
278  const labelPair& edgeData = allEdgeData[edgeI];
279 
280  const edge& e = patch.edges()[edgeI];
281  const labelList& eFaces = patch.edgeFaces()[edgeI];
282 
283  forAll(eFaces, i)
284  {
285  label facei = eFaces[i];
286  const face& f = patch.localFaces()[facei];
287 
288  // Combine edgeData with face data
289  label fp0 = f.find(e[0]);
290  if (pointGlobalRegions[facei][fp0] > edgeData[0])
291  {
292  pointGlobalRegions[facei][fp0] = edgeData[0];
293  if (!isChangedFace[facei])
294  {
295  isChangedFace.set(facei);
296  changedFaces.append(facei);
297  }
298  }
299 
300  label fp1 = f.find(e[1]);
301  if (pointGlobalRegions[facei][fp1] > edgeData[1])
302  {
303  pointGlobalRegions[facei][fp1] = edgeData[1];
304  if (!isChangedFace[facei])
305  {
306  isChangedFace.set(facei);
307  changedFaces.append(facei);
308  }
309  }
310  }
311  }
312 
313 
314  label nChangedFaces = returnReduce(changedFaces.size(), sumOp<label>());
315  if (nChangedFaces == 0)
316  {
317  break;
318  }
319 
320 
321  // From face to edge
322  // ~~~~~~~~~~~~~~~~~
323 
324  isChangedEdge = false;
325  changedEdges.clear();
326 
327  forAll(changedFaces, i)
328  {
329  label facei = changedFaces[i];
330  const face& f = patch.localFaces()[facei];
331  const labelList& fEdges = patch.faceEdges()[facei];
332 
333  forAll(fEdges, fp)
334  {
335  label edgeI = fEdges[fp];
336 
337  if (!nonManifoldEdge[edgeI])
338  {
339  const edge& e = patch.edges()[edgeI];
340  label fp0 = f.find(e[0]);
341  label region0 = pointGlobalRegions[facei][fp0];
342  label fp1 = f.find(e[1]);
343  label region1 = pointGlobalRegions[facei][fp1];
344 
345  if
346  (
347  (allEdgeData[edgeI][0] > region0)
348  || (allEdgeData[edgeI][1] > region1)
349  )
350  {
351  allEdgeData[edgeI] = labelPair(region0, region1);
352  if (!isChangedEdge[edgeI])
353  {
354  changedEdges.append(edgeI);
355  isChangedEdge.set(edgeI);
356  }
357  }
358  }
359  }
360  }
361 
362  syncEdges
363  (
364  globalData,
365 
366  patchEdges,
367  coupledEdges,
368  sameEdgeOrientation,
369  syncNonCollocated,
370 
371  isChangedEdge,
372  changedEdges,
373  allEdgeData
374  );
375 
376 
377  label nChangedEdges = returnReduce(changedEdges.size(), sumOp<label>());
378  if (nChangedEdges == 0)
379  {
380  break;
381  }
382  }
383 
384 
385 
386  // Assign local regions
387  // ~~~~~~~~~~~~~~~~~~~~
388 
389  // Calculate addressing from global region back to local region
390  pointLocalRegions.setSize(patch.size());
391  Map<label> globalToLocalRegion(globalRegions.localSize()/4);
392  DynamicList<label> dynLocalToGlobalRegion(globalToLocalRegion.size());
393  forAll(patch.localFaces(), facei)
394  {
395  const face& f = patch.localFaces()[facei];
396  face& pRegions = pointLocalRegions[facei];
397  pRegions.setSize(f.size());
398  forAll(f, fp)
399  {
400  label globalRegionI = pointGlobalRegions[facei][fp];
401 
402  const auto fnd = globalToLocalRegion.cfind(globalRegionI);
403 
404  if (fnd.found())
405  {
406  // Already encountered this global region. Assign same local one
407  pRegions[fp] = fnd();
408  }
409  else
410  {
411  // Region not yet seen. Create new one
412  label localRegionI = globalToLocalRegion.size();
413  pRegions[fp] = localRegionI;
414  globalToLocalRegion.insert(globalRegionI, localRegionI);
415  dynLocalToGlobalRegion.append(globalRegionI);
416  }
417  }
418  }
419  localToGlobalRegion.transfer(dynLocalToGlobalRegion);
420 }
421 
422 
423 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
424 
425 Foam::createShellMesh::createShellMesh
426 (
427  const primitiveFacePatch& patch,
428  const faceList& pointRegions,
429  const labelList& regionPoints
430 )
431 :
432  patch_(patch),
433  pointRegions_(pointRegions),
434  regionPoints_(regionPoints)
435 {
436  if (pointRegions_.size() != patch_.size())
437  {
439  << "nFaces:" << patch_.size()
440  << " pointRegions:" << pointRegions.size()
441  << exit(FatalError);
442  }
443 }
444 
445 
446 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
447 
449 (
450  const pointField& firstLayerDisp,
451  const scalar expansionRatio,
452  const label nLayers,
453  const labelList& topPatchID,
454  const labelList& bottomPatchID,
455  const labelListList& extrudeEdgePatches,
456  polyTopoChange& meshMod
457 )
458 {
459  if (firstLayerDisp.size() != regionPoints_.size())
460  {
462  << "nRegions:" << regionPoints_.size()
463  << " firstLayerDisp:" << firstLayerDisp.size()
464  << exit(FatalError);
465  }
466 
467  if
468  (
469  topPatchID.size() != patch_.size()
470  && bottomPatchID.size() != patch_.size()
471  )
472  {
474  << "nFaces:" << patch_.size()
475  << " topPatchID:" << topPatchID.size()
476  << " bottomPatchID:" << bottomPatchID.size()
477  << exit(FatalError);
478  }
479 
480  if (extrudeEdgePatches.size() != patch_.nEdges())
481  {
483  << "nEdges:" << patch_.nEdges()
484  << " extrudeEdgePatches:" << extrudeEdgePatches.size()
485  << exit(FatalError);
486  }
487 
488 
489 
490  // From cell to patch (trivial)
491  DynamicList<label> cellToFaceMap(nLayers*patch_.size());
492  // From face to patch+turning index
493  DynamicList<label> faceToFaceMap
494  (
495  (nLayers+1)*(patch_.size()+patch_.nEdges())
496  );
497  // From face to patch edge index
498  DynamicList<label> faceToEdgeMap(nLayers*(patch_.nEdges()+patch_.nEdges()));
499  // From point to patch point index
500  DynamicList<label> pointToPointMap((nLayers+1)*patch_.nPoints());
501 
502 
503  // Introduce new cell for every face
504  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
505 
506  labelList addedCells(nLayers*patch_.size());
507  forAll(patch_, facei)
508  {
509  for (label layerI = 0; layerI < nLayers; layerI++)
510  {
511  addedCells[nLayers*facei+layerI] = meshMod.addCell
512  (
513  -1, // masterPointID
514  -1, // masterEdgeID
515  -1, // masterFaceID
516  cellToFaceMap.size(), // masterCellID
517  -1 // zoneID
518  );
519  cellToFaceMap.append(facei);
520  }
521  }
522 
523 
524  // Introduce original points
525  // ~~~~~~~~~~~~~~~~~~~~~~~~~
526 
527  // Original point numbers in local point ordering so no need to store.
528  forAll(patch_.localPoints(), pointi)
529  {
530  //label addedPointi =
531  meshMod.addPoint
532  (
533  patch_.localPoints()[pointi], // point
534  pointToPointMap.size(), // masterPointID
535  -1, // zoneID
536  true // inCell
537  );
538  pointToPointMap.append(pointi);
539 
540  //Pout<< "Added bottom point " << addedPointi
541  // << " at " << patch_.localPoints()[pointi]
542  // << " from point " << pointi
543  // << endl;
544  }
545 
546 
547  // Introduce new points (one for every region)
548  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
549 
550  labelList addedPoints(nLayers*regionPoints_.size());
551  forAll(regionPoints_, regionI)
552  {
553  label pointi = regionPoints_[regionI];
554 
555  point pt = patch_.localPoints()[pointi];
556  point disp = firstLayerDisp[regionI];
557  for (label layerI = 0; layerI < nLayers; layerI++)
558  {
559  pt += disp;
560 
561  addedPoints[nLayers*regionI+layerI] = meshMod.addPoint
562  (
563  pt, // point
564  pointToPointMap.size(), // masterPointID - used only addressing
565  -1, // zoneID
566  true // inCell
567  );
568  pointToPointMap.append(pointi);
569 
570  disp *= expansionRatio;
571  }
572  }
573 
574 
575  // Add face on bottom side
576  forAll(patch_.localFaces(), facei)
577  {
578  meshMod.addFace
579  (
580  patch_.localFaces()[facei].reverseFace(),// vertices
581  addedCells[nLayers*facei], // own
582  -1, // nei
583  -1, // masterPointID
584  -1, // masterEdgeID
585  faceToFaceMap.size(), // masterFaceID : current facei
586  true, // flipFaceFlux
587  bottomPatchID[facei], // patchID
588  -1, // zoneID
589  false // zoneFlip
590  );
591  faceToFaceMap.append(-facei-1); // points to flipped original face
592  faceToEdgeMap.append(-1);
593 
594  //const face newF(patch_.localFaces()[facei].reverseFace());
595  //Pout<< "Added bottom face "
596  // << newF
597  // << " coords:" << UIndirectList<point>(meshMod.points(), newF)
598  // << " own " << addedCells[facei]
599  // << " patch:" << bottomPatchID[facei]
600  // << " at " << patch_.faceCentres()[facei]
601  // << endl;
602  }
603 
604  // Add inbetween faces and face on top
605  forAll(patch_.localFaces(), facei)
606  {
607  // Get face in original ordering
608  const face& f = patch_.localFaces()[facei];
609 
610  face newF(f.size());
611 
612  for (label layerI = 0; layerI < nLayers; layerI++)
613  {
614  // Pick up point based on region and layer
615  forAll(f, fp)
616  {
617  label region = pointRegions_[facei][fp];
618  newF[fp] = addedPoints[region*nLayers+layerI];
619  }
620 
621  label own = addedCells[facei*nLayers+layerI];
622  label nei;
623  label patchi;
624  if (layerI == nLayers-1)
625  {
626  nei = -1;
627  patchi = topPatchID[facei];
628  }
629  else
630  {
631  nei = addedCells[facei*nLayers+layerI+1];
632  patchi = -1;
633  }
634 
635  meshMod.addFace
636  (
637  newF, // vertices
638  own, // own
639  nei, // nei
640  -1, // masterPointID
641  -1, // masterEdgeID
642  faceToFaceMap.size(), // masterFaceID : current facei
643  false, // flipFaceFlux
644  patchi, // patchID
645  -1, // zoneID
646  false // zoneFlip
647  );
648  faceToFaceMap.append(facei+1); // unflipped
649  faceToEdgeMap.append(-1);
650 
651  //Pout<< "Added inbetween face " << newF
652  // << " coords:" << UIndirectList<point>(meshMod.points(), newF)
653  // << " at layer " << layerI
654  // << " own " << own
655  // << " nei " << nei
656  // << " at " << patch_.faceCentres()[facei]
657  // << endl;
658  }
659  }
660 
661 
662  // Add side faces
663  // ~~~~~~~~~~~~~~
664  // Note that we loop over edges multiple times so for edges with
665  // two cyclic faces they get added in two passes (for correct ordering)
666 
667  // Pass1. Internal edges and first face of other edges
668  forAll(extrudeEdgePatches, edgeI)
669  {
670  const labelList& eFaces = patch_.edgeFaces()[edgeI];
671  const labelList& ePatches = extrudeEdgePatches[edgeI];
672 
673  if (ePatches.size() == 0)
674  {
675  // Internal face
676  if (eFaces.size() != 2)
677  {
679  << "edge:" << edgeI
680  << " not internal but does not have side-patches defined."
681  << exit(FatalError);
682  }
683  }
684  else
685  {
686  if (eFaces.size() != ePatches.size())
687  {
689  << "external/feature edge:" << edgeI
690  << " has " << eFaces.size() << " connected extruded faces "
691  << " but only " << ePatches.size()
692  << " boundary faces defined." << exit(FatalError);
693  }
694  }
695 
696 
697 
698  // Make face pointing in to eFaces[0] so out of new master face
699  const face& f = patch_.localFaces()[eFaces[0]];
700  const edge& e = patch_.edges()[edgeI];
701 
702  label fp0 = f.find(e[0]);
703  label fp1 = f.fcIndex(fp0);
704 
705  if (f[fp1] != e[1])
706  {
707  fp1 = fp0;
708  fp0 = f.rcIndex(fp1);
709  }
710 
711  face newF(4);
712 
713  for (label layerI = 0; layerI < nLayers; layerI++)
714  {
715  label region0 = pointRegions_[eFaces[0]][fp0];
716  label region1 = pointRegions_[eFaces[0]][fp1];
717 
718  // Pick up points with correct normal
719  if (layerI == 0)
720  {
721  newF[0] = f[fp0];
722  newF[1] = f[fp1];
723  newF[2] = addedPoints[nLayers*region1+layerI];
724  newF[3] = addedPoints[nLayers*region0+layerI];
725  }
726  else
727  {
728  newF[0] = addedPoints[nLayers*region0+layerI-1];
729  newF[1] = addedPoints[nLayers*region1+layerI-1];
730  newF[2] = addedPoints[nLayers*region1+layerI];
731  newF[3] = addedPoints[nLayers*region0+layerI];
732  }
733 
734  // Optionally rotate so e[0] is always 0th vertex. Note that
735  // this normally is automatically done by coupled face ordering
736  // but with NOORDERING we have to do it ourselves.
737  if (f[fp0] != e[0])
738  {
739  // rotate one back to get newF[1] (originating from e[0])
740  // into newF[0]
741  label v0 = newF[0];
742  for (label i = 0; i < newF.size()-1; i++)
743  {
744  newF[i] = newF[newF.fcIndex(i)];
745  }
746  newF.last() = v0;
747  }
748 
749 
750  label minCelli = addedCells[nLayers*eFaces[0]+layerI];
751  label maxCelli;
752  label patchi;
753  if (ePatches.size() == 0)
754  {
755  maxCelli = addedCells[nLayers*eFaces[1]+layerI];
756  if (minCelli > maxCelli)
757  {
758  // Swap
759  std::swap(minCelli, maxCelli);
760  newF.flip();
761  }
762  patchi = -1;
763  }
764  else
765  {
766  maxCelli = -1;
767  patchi = ePatches[0];
768  }
769 
770  //{
771  // Pout<< "Adding from face:" << patch_.faceCentres()[eFaces[0]]
772  // << " from edge:"
773  // << patch_.localPoints()[f[fp0]]
774  // << patch_.localPoints()[f[fp1]]
775  // << " at layer:" << layerI
776  // << " with new points:" << newF
777  // << " locations:"
778  // << UIndirectList<point>(meshMod.points(), newF)
779  // << " own:" << minCelli
780  // << " nei:" << maxCelli
781  // << endl;
782  //}
783 
784 
785  // newF already outwards pointing.
786  meshMod.addFace
787  (
788  newF, // vertices
789  minCelli, // own
790  maxCelli, // nei
791  -1, // masterPointID
792  -1, // masterEdgeID
793  faceToFaceMap.size(), // masterFaceID
794  false, // flipFaceFlux
795  patchi, // patchID
796  -1, // zoneID
797  false // zoneFlip
798  );
799  faceToFaceMap.append(0);
800  faceToEdgeMap.append(edgeI);
801  }
802  }
803 
804  // Pass2. Other faces of boundary edges
805  forAll(extrudeEdgePatches, edgeI)
806  {
807  const labelList& eFaces = patch_.edgeFaces()[edgeI];
808  const labelList& ePatches = extrudeEdgePatches[edgeI];
809 
810  if (ePatches.size() >= 2)
811  {
812  for (label i = 1; i < ePatches.size(); i++)
813  {
814  // Extrude eFaces[i]
815  label minFacei = eFaces[i];
816 
817  // Make face pointing in to eFaces[0] so out of new master face
818  const face& f = patch_.localFaces()[minFacei];
819 
820  const edge& e = patch_.edges()[edgeI];
821  label fp0 = f.find(e[0]);
822  label fp1 = f.fcIndex(fp0);
823 
824  if (f[fp1] != e[1])
825  {
826  fp1 = fp0;
827  fp0 = f.rcIndex(fp1);
828  }
829 
830  face newF(4);
831  for (label layerI = 0; layerI < nLayers; layerI++)
832  {
833  label region0 = pointRegions_[minFacei][fp0];
834  label region1 = pointRegions_[minFacei][fp1];
835 
836  if (layerI == 0)
837  {
838  newF[0] = f[fp0];
839  newF[1] = f[fp1];
840  newF[2] = addedPoints[nLayers*region1+layerI];
841  newF[3] = addedPoints[nLayers*region0+layerI];
842  }
843  else
844  {
845  newF[0] = addedPoints[nLayers*region0+layerI-1];
846  newF[1] = addedPoints[nLayers*region1+layerI-1];
847  newF[2] = addedPoints[nLayers*region1+layerI];
848  newF[3] = addedPoints[nLayers*region0+layerI];
849  }
850 
851 
852  // Optionally rotate so e[0] is always 0th vertex. Note that
853  // this normally is automatically done by coupled face
854  // ordering but with NOORDERING we have to do it ourselves.
855  if (f[fp0] != e[0])
856  {
857  // rotate one back to get newF[1] (originating
858  // from e[0]) into newF[0].
859  label v0 = newF[0];
860  for (label i = 0; i < newF.size()-1; i++)
861  {
862  newF[i] = newF[newF.fcIndex(i)];
863  }
864  newF.last() = v0;
865  }
867  //{
868  // Pout<< "Adding from MULTI face:"
869  // << patch_.faceCentres()[minFacei]
870  // << " from edge:"
871  // << patch_.localPoints()[f[fp0]]
872  // << patch_.localPoints()[f[fp1]]
873  // << " at layer:" << layerI
874  // << " with new points:" << newF
875  // << " locations:"
876  // << UIndirectList<point>(meshMod.points(), newF)
877  // << endl;
878  //}
879 
880  // newF already outwards pointing.
881  meshMod.addFace
882  (
883  newF, // vertices
884  addedCells[nLayers*minFacei+layerI], // own
885  -1, // nei
886  -1, // masterPointID
887  -1, // masterEdgeID
888  faceToFaceMap.size(), // masterFaceID
889  false, // flipFaceFlux
890  ePatches[i], // patchID
891  -1, // zoneID
892  false // zoneFlip
893  );
894  faceToFaceMap.append(0);
895  faceToEdgeMap.append(edgeI);
896  }
897  }
898  }
899  }
900 
901 
902  cellToFaceMap_.transfer(cellToFaceMap);
903  faceToFaceMap_.transfer(faceToFaceMap);
904  faceToEdgeMap_.transfer(faceToEdgeMap);
905  pointToPointMap_.transfer(pointToPointMap);
906 }
907 
908 
910 {
911  inplaceReorder(map.reverseCellMap(), cellToFaceMap_);
912  inplaceReorder(map.reverseFaceMap(), faceToFaceMap_);
913  inplaceReorder(map.reverseFaceMap(), faceToEdgeMap_);
914  inplaceReorder(map.reversePointMap(), pointToPointMap_);
915 }
916 
917 
918 // ************************************************************************* //
Foam::polyTopoChange::addCell
label addCell(const label masterPointID, const label masterEdgeID, const label masterFaceID, const label masterCellID, const label zoneID)
Add cell. Return new cell label.
Definition: polyTopoChange.C:2895
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
meshTools.H
Foam::labelMax
constexpr label labelMax
Definition: label.H:61
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:63
PatchTools.H
Foam::DynamicList< label >
polyAddFace.H
globalMeshData.H
mapPolyMesh.H
globalIndex.H
Foam::labelPairList
List< labelPair > labelPairList
List of labelPairs.
Definition: labelPair.H:64
Foam::polyTopoChange::addFace
label addFace(const face &f, const label own, const label nei, const label masterPointID, const label masterEdgeID, const label masterFaceID, const bool flipFaceFlux, const label patchID, const label zoneID, const bool zoneFlip)
Add face to cells. Return new face label.
Definition: polyTopoChange.C:2747
polyTopoChange.H
Foam::polyTopoChange::addPoint
label addPoint(const point &pt, const label masterPointID, const label zoneID, const bool inCell)
Add point. Return new point label.
Definition: polyTopoChange.C:2610
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:99
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::Map< label >
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:187
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::sumOp
Definition: ops.H:213
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::createShellMesh::updateMesh
void updateMesh(const mapPolyMesh &)
Update any locally stored mesh information.
Definition: createShellMesh.C:909
createShellMesh.H
Foam::Field< vector >
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
polyAddPoint.H
indirectPrimitivePatch.H
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::FatalError
error FatalError
Foam::face::flip
void flip()
Flip the face in-place.
Definition: face.C:500
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:107
Foam::createShellMesh::setRefinement
void setRefinement(const pointField &firstLayerThickness, const scalar expansionRatio, const label nLayers, const labelList &topPatchID, const labelList &bottomPatchID, const labelListList &extrudeEdgePatches, polyTopoChange &meshMod)
Play commands into polyTopoChange to create layer mesh.
Definition: createShellMesh.C:449
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
polyModifyFace.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:501
polyAddCell.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::mapPolyMesh::reverseCellMap
const labelList & reverseCellMap() const
Reverse cell map.
Definition: mapPolyMesh.H:532
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pair< label >
mapDistribute.H
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mapPolyMesh::reversePointMap
const labelList & reversePointMap() const
Reverse point map.
Definition: mapPolyMesh.H:469
Foam::Vector< scalar >
Foam::createShellMesh::calcPointRegions
static void calcPointRegions(const globalMeshData &globalData, const primitiveFacePatch &patch, const bitSet &nonManifoldEdge, const bool syncNonCollocated, faceList &pointGlobalRegions, faceList &pointLocalRegions, labelList &localToGlobalRegion)
Helper: calculate point regions. The point region is the.
Definition: createShellMesh.C:157
Foam::List< face >
Foam::globalMeshData::coupledPatch
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Definition: globalMeshData.C:2046
Foam::minEqOp< labelPair >::operator()
void operator()(labelPair &x, const labelPair &y) const
Definition: createShellMesh.C:54
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::inplaceReorder
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
Definition: ListOpsTemplates.C:124
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::VectorSpace::size
static constexpr direction size() noexcept
The number of elements in the VectorSpace = Ncmpts.
Definition: VectorSpace.H:176
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:240
y
scalar y
Definition: LISASMDCalcMethod1.H:14
labelPair.H
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79