mergePolyMesh.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "mergePolyMesh.H"
30 #include "Time.H"
31 #include "polyTopoChanger.H"
32 #include "mapPolyMesh.H"
33 #include "polyAddPoint.H"
34 #include "polyAddCell.H"
35 #include "polyAddFace.H"
36 #include "processorPolyPatch.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(mergePolyMesh, 1);
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 Foam::label Foam::mergePolyMesh::patchIndex(const polyPatch& p)
49 {
50  // Find the patch name on the list. If the patch is already there
51  // and patch types match, return index
52  const word& pType = p.type();
53  const word& pName = p.name();
54 
55  bool nameFound = false;
56 
57  forAll(patchNames_, patchi)
58  {
59  if (patchNames_[patchi] == pName)
60  {
61  if (patchDicts_[patchi].get<word>("type") == pType)
62  {
63  // Found name and types match
64  return patchi;
65  }
66  else
67  {
68  // Found the name, but type is different
69  nameFound = true;
70  }
71  }
72  }
73 
74  // Patch not found. Append to the list
75  {
76  OStringStream os;
77  p.write(os);
78  patchDicts_.append(dictionary(IStringStream(os.str())()));
79  }
80 
81  if (nameFound)
82  {
83  // Duplicate name is not allowed. Create a composite name from the
84  // patch name and case name
85  const word& caseName = p.boundaryMesh().mesh().time().caseName();
86 
87  patchNames_.append(pName + "_" + caseName);
88 
89  Info<< "label patchIndex(const polyPatch& p) : "
90  << "Patch " << p.index() << " named "
91  << pName << " in mesh " << caseName
92  << " already exists, but patch types "
93  << " do not match.\nCreating a composite name as "
94  << patchNames_.last() << endl;
95  }
96  else
97  {
98  patchNames_.append(pName);
99  }
100 
101  return patchNames_.size() - 1;
102 }
103 
104 
105 Foam::label Foam::mergePolyMesh::zoneIndex
106 (
107  DynamicList<word>& names,
108  const word& curName
109 )
110 {
111  forAll(names, zonei)
112  {
113  if (names[zonei] == curName)
114  {
115  return zonei;
116  }
117  }
118 
119  // Not found. Add new name to the list
120  names.append(curName);
121 
122  return names.size() - 1;
123 }
124 
125 
126 void Foam::mergePolyMesh::sortProcessorPatches()
127 {
128  Info<< "Reordering processor patches last" << endl;
129 
130  // Updates boundaryMesh() and meshMod_ to guarantee processor patches
131  // are last. This could be done inside the merge() but it is far easier
132  // to do separately.
133 
134 
135  // 1. Shuffle the patches in the boundaryMesh
136 
137  const polyBoundaryMesh& oldPatches = boundaryMesh();
138 
139  polyPatchList newPatches(oldPatches.size());
140 
141  labelList oldToSorted(oldPatches.size());
142 
143  label nPatches = 0;
144 
145  forAll(oldPatches, patchi)
146  {
147  const polyPatch& pp = oldPatches[patchi];
148 
149  if (!isA<processorPolyPatch>(pp))
150  {
151  newPatches.set
152  (
153  nPatches,
154  pp.clone
155  (
156  oldPatches,
157  nPatches,
158  0,
159  nInternalFaces()
160  )
161  );
162 
163  oldToSorted[patchi] = nPatches;
164  ++nPatches;
165  }
166  }
167 
168  forAll(oldPatches, patchi)
169  {
170  const polyPatch& pp = oldPatches[patchi];
171 
172  if (isA<processorPolyPatch>(pp))
173  {
174  newPatches.set
175  (
176  nPatches,
177  pp.clone
178  (
179  oldPatches,
180  oldToSorted[patchi],
181  0,
182  nInternalFaces()
183  )
184  );
185 
186  oldToSorted[patchi] = nPatches;
187  ++nPatches;
188  }
189  }
190 
191 
192  removeBoundary();
193  addPatches(newPatches);
194 
195 
196  // Update the polyTopoChange
197  DynamicList<label>& patchID = const_cast<DynamicList<label>&>
198  (
199  meshMod_.region()
200  );
201 
202  forAll(patchID, facei)
203  {
204  label patchi = patchID[facei];
205  if (patchi != -1)
206  {
207  patchID[facei] = oldToSorted[patchID[facei]];
208  }
209  }
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
214 
215 Foam::mergePolyMesh::mergePolyMesh(const IOobject& io)
216 :
217  polyMesh(io),
218  meshMod_(*this),
219  patchNames_(2*boundaryMesh().size()),
220  patchDicts_(2*boundaryMesh().size()),
221  pointZoneNames_(),
222  faceZoneNames_(),
223  cellZoneNames_()
224 {
225  // Insert the original patches into the list
226  wordList curPatchNames = boundaryMesh().names();
227 
228  forAll(boundaryMesh(), patchi)
229  {
230  patchNames_.append(boundaryMesh()[patchi].name());
231 
232  OStringStream os;
233  boundaryMesh()[patchi].write(os);
234  patchDicts_.append(dictionary(IStringStream(os.str())()));
235  }
236 
237  // Insert point, face and cell zones into the list
238 
239  // Point zones
240  wordList curPointZoneNames = pointZones().names();
241  if (curPointZoneNames.size())
242  {
243  pointZoneNames_.setCapacity(2*curPointZoneNames.size());
244  pointZoneNames_.append(curPointZoneNames);
245  }
246 
247  // Face zones
248  wordList curFaceZoneNames = faceZones().names();
249  if (curFaceZoneNames.size())
250  {
251  faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
252  faceZoneNames_.append(curFaceZoneNames);
253  }
254 
255  // Cell zones
256  wordList curCellZoneNames = cellZones().names();
257  if (curCellZoneNames.size())
258  {
259  cellZoneNames_.setCapacity(2*curCellZoneNames.size());
260  cellZoneNames_.append(curCellZoneNames);
261  }
262 }
263 
264 
265 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
266 
267 void Foam::mergePolyMesh::addMesh(const polyMesh& m)
268 {
269  // Add all the points, faces and cells of the new mesh
270 
271  // Add points
272 
273  label zoneID = -1;
274 
275  const pointField& p = m.points();
276  labelList renumberPoints(p.size());
277 
278  const pointZoneMesh& pz = m.pointZones();
279  labelList pointZoneIndices(pz.size());
280 
281  forAll(pz, zoneI)
282  {
283  pointZoneIndices[zoneI] = zoneIndex(pointZoneNames_, pz[zoneI].name());
284  }
285 
286  forAll(p, pointi)
287  {
288  // Grab zone ID. If a point is not in a zone, it will return -1
289  zoneID = pz.whichZone(pointi);
290 
291  if (zoneID >= 0)
292  {
293  // Translate zone ID into the new index
294  zoneID = pointZoneIndices[zoneID];
295  }
296 
297  renumberPoints[pointi] =
298  meshMod_.setAction
299  (
300  polyAddPoint
301  (
302  p[pointi], // Point to add
303  -1, // Master point (straight addition)
304  zoneID, // Zone for point
305  pointi < m.nPoints() // Is in cell?
306  )
307  );
308  }
309 
310  // Add cells
311 
312  const cellList& c = m.cells();
313  labelList renumberCells(c.size());
314 
315  const cellZoneMesh& cz = m.cellZones();
316  labelList cellZoneIndices(cz.size());
317 
318  forAll(cz, zoneI)
319  {
320  cellZoneIndices[zoneI] = zoneIndex(cellZoneNames_, cz[zoneI].name());
321  }
322 
323  forAll(c, celli)
324  {
325  // Grab zone ID. If a cell is not in a zone, it will return -1
326  zoneID = cz.whichZone(celli);
327 
328  if (zoneID >= 0)
329  {
330  // Translate zone ID into the new index
331  zoneID = cellZoneIndices[zoneID];
332  }
333 
334  renumberCells[celli] =
335  meshMod_.setAction
336  (
337  polyAddCell
338  (
339  -1, // Master point
340  -1, // Master edge
341  -1, // Master face
342  -1, // Master cell
343  zoneID // Zone for cell
344  )
345  );
346  }
347 
348  // Add faces
349  const polyBoundaryMesh& bm = m.boundaryMesh();
350 
351  // Gather the patch indices
352  labelList patchIndices(bm.size());
353 
354  forAll(patchIndices, patchi)
355  {
356  patchIndices[patchi] = patchIndex(bm[patchi]);
357  }
358 
359  // Temporary: update number of allowable patches. This should be
360  // determined at the top - before adding anything.
361  meshMod_.setNumPatches(patchNames_.size());
362 
363 
364 
365  const faceZoneMesh& fz = m.faceZones();
366  labelList faceZoneIndices(fz.size());
367 
368  forAll(fz, zoneI)
369  {
370  faceZoneIndices[zoneI] = zoneIndex(faceZoneNames_, fz[zoneI].name());
371  }
372 
373  const faceList& f = m.faces();
374  labelList renumberFaces(f.size());
375 
376  const labelList& own = m.faceOwner();
377  const labelList& nei = m.faceNeighbour();
378 
379  label newOwn, newNei, newPatch, newZone;
380  bool newZoneFlip;
381 
382  forAll(f, facei)
383  {
384  const face& curFace = f[facei];
385 
386  face newFace(curFace.size());
387 
388  forAll(curFace, pointi)
389  {
390  newFace[pointi] = renumberPoints[curFace[pointi]];
391  }
392 
393  if (debug)
394  {
395  // Check that the face is valid
396  if (min(newFace) < 0)
397  {
399  << "Error in point mapping for face " << facei
400  << ". Old face: " << curFace << " New face: " << newFace
401  << abort(FatalError);
402  }
403  }
404 
405  if (facei < m.nInternalFaces() || facei >= m.nFaces())
406  {
407  newPatch = -1;
408  }
409  else
410  {
411  newPatch = patchIndices[bm.whichPatch(facei)];
412  }
413 
414  newOwn = own[facei];
415  if (newOwn > -1) newOwn = renumberCells[newOwn];
416 
417  if (newPatch > -1)
418  {
419  newNei = -1;
420  }
421  else
422  {
423  newNei = nei[facei];
424  newNei = renumberCells[newNei];
425  }
426 
427 
428  newZone = fz.whichZone(facei);
429  newZoneFlip = false;
430 
431  if (newZone >= 0)
432  {
433  newZoneFlip = fz[newZone].flipMap()[fz[newZone].whichFace(facei)];
434 
435  // Grab the new zone
436  newZone = faceZoneIndices[newZone];
437  }
438 
439  renumberFaces[facei] =
440  meshMod_.setAction
441  (
442  polyAddFace
443  (
444  newFace,
445  newOwn,
446  newNei,
447  -1,
448  -1,
449  -1,
450  false,
451  newPatch,
452  newZone,
453  newZoneFlip
454  )
455  );
456  }
457 }
458 
459 
461 {
462  Info<< "patch names: " << patchNames_ << nl
463  << "patch dicts: " << patchDicts_ << nl
464  << "point zone names: " << pointZoneNames_ << nl
465  << "face zone names: " << faceZoneNames_ << nl
466  << "cell zone names: " << cellZoneNames_ << endl;
467 
468  // Add the patches if necessary
469  if (patchNames_.size() != boundaryMesh().size())
470  {
471  Info<< "Copying old patches" << endl;
472 
473  polyPatchList newPatches(patchNames_.size());
474 
475  const polyBoundaryMesh& oldPatches = boundaryMesh();
476 
477  // Note. Re-using counter in two for loops
478  label patchi = 0;
479 
480  for (patchi = 0; patchi < oldPatches.size(); patchi++)
481  {
482  newPatches.set
483  (
484  patchi,
485  oldPatches[patchi].clone(oldPatches)
486  );
487  }
488 
489  Info<< "Adding new patches. " << endl;
490 
491  label endOfLastPatch =
492  patchi == 0
493  ? 0
494  : oldPatches[patchi - 1].start() + oldPatches[patchi - 1].size();
495 
496  for (; patchi < patchNames_.size(); patchi++)
497  {
498  // Add a patch
499  dictionary dict(patchDicts_[patchi]);
500  dict.set("nFaces", 0);
501  dict.set("startFace", endOfLastPatch);
502 
503  newPatches.set
504  (
505  patchi,
506  (
508  (
509  patchNames_[patchi],
510  dict,
511  patchi,
512  oldPatches
513  )
514  )
515  );
516  }
517 
518  removeBoundary();
519  addPatches(newPatches);
520  }
521 
522  // Add the zones if necessary
523  if (pointZoneNames_.size() > pointZones().size())
524  {
525  Info<< "Adding new pointZones." << endl;
526 
527  label zonei = pointZones().size(); // continue from here
528 
529  const label nZones = pointZoneNames_.size();
530 
531  pointZones().setSize(nZones);
532 
533  for (/*nil*/; zonei < nZones; ++zonei)
534  {
535  pointZones().set
536  (
537  zonei,
538  new pointZone(pointZoneNames_[zonei], zonei, pointZones())
539  );
540  }
541  }
542  if (cellZoneNames_.size() > cellZones().size())
543  {
544  Info<< "Adding new cellZones." << endl;
545 
546  label zonei = cellZones().size(); // continue from here
547 
548  const label nZones = cellZoneNames_.size();
549 
550  cellZones().setSize(cellZoneNames_.size());
551 
552  for (/*nil*/; zonei < nZones; ++zonei)
553  {
554  cellZones().set
555  (
556  zonei,
557  new cellZone(cellZoneNames_[zonei], zonei, cellZones())
558  );
559  }
560  }
561  if (faceZoneNames_.size() > faceZones().size())
562  {
563  Info<< "Adding new faceZones." << endl;
564 
565  label zonei = faceZones().size(); // continue from here
566 
567  const label nZones = faceZoneNames_.size();
568 
569  faceZones().setSize(nZones);
570 
571  for (/*nil*/; zonei < nZones; ++zonei)
572  {
573  faceZones().set
574  (
575  zonei,
576  new faceZone(faceZoneNames_[zonei], zonei, faceZones())
577  );
578  }
579  }
580 
581 
582  // Shuffle the processor patches to be last
583  sortProcessorPatches();
584 
585  // Change mesh. No inflation
586  meshMod_.changeMesh(*this, false);
587 
588  // Clear topo change for the next operation
589  meshMod_.clear();
590 }
591 
592 
593 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:491
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
List< face > faceList
A List of faces.
Definition: faceListFwd.H:41
const word & name() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void addMesh(const polyMesh &m)
Add a mesh.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:558
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void merge()
Merge meshes.
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
T & last()
Access last element of the list, position [size()-1].
Definition: UList.H:812
List< word > wordList
A List of words.
Definition: fileName.H:58
const dimensionedScalar c
Speed of light in a vacuum.
messageStream Info
Information stream (stdout output on master, null elsewhere)
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return pointer to a new patch created on freestore from components.
Definition: polyPatchNew.C:28
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
const fileName & caseName() const
Return the Time::caseName()
Definition: IOobject.C:459
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:41
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
Namespace for OpenFOAM.