mergePolyMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "mergePolyMesh.H"
27 #include "Time.H"
28 #include "polyTopoChangeMap.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(mergePolyMesh, 0);
35 }
36 
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 Foam::label Foam::mergePolyMesh::patchIndex(const polyPatch& p)
41 {
42  // Find the patch name on the list. If the patch is already there
43  // and patch types match, return index
44  const word& pType = p.type();
45  const word& pName = p.name();
46 
47  bool nameFound = false;
48 
49  forAll(patchNames_, patchi)
50  {
51  if (patchNames_[patchi] == pName)
52  {
53  if (word(patchDicts_[patchi]["type"]) == pType)
54  {
55  // Found name and types match
56  return patchi;
57  }
58  else
59  {
60  // Found the name, but type is different
61  nameFound = true;
62  }
63  }
64  }
65 
66  // Patch not found. Append to the list
67  {
68  OStringStream os;
69  p.write(os);
70  patchDicts_.append(dictionary(IStringStream(os.str())()));
71  }
72 
73  if (nameFound)
74  {
75  // Duplicate name is not allowed. Create a composite name from the
76  // patch name and case name
77  const word& caseName = p.boundaryMesh().mesh().time().caseName();
78 
79  patchNames_.append(pName + "_" + caseName);
80 
81  Info<< "label patchIndex(const polyPatch& p) : "
82  << "Patch " << p.index() << " named "
83  << pName << " in mesh " << caseName
84  << " already exists, but patch types "
85  << " do not match.\nCreating a composite name as "
86  << patchNames_.last() << endl;
87  }
88  else
89  {
90  patchNames_.append(pName);
91  }
92 
93  return patchNames_.size() - 1;
94 }
95 
96 
97 Foam::label Foam::mergePolyMesh::zoneIndex
98 (
99  DynamicList<word>& names,
100  const word& curName
101 )
102 {
103  forAll(names, zonei)
104  {
105  if (names[zonei] == curName)
106  {
107  return zonei;
108  }
109  }
110 
111  // Not found. Add new name to the list
112  names.append(curName);
113 
114  return names.size() - 1;
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
121 :
122  mesh_(mesh),
123  meshMod_(mesh_),
124  patchNames_(2*mesh_.boundaryMesh().size()),
125  patchDicts_(2*mesh_.boundaryMesh().size()),
126  pointZoneNames_(),
127  faceZoneNames_(),
128  cellZoneNames_()
129 {
130  // Insert the original patches into the list
131  wordList curPatchNames = mesh_.boundaryMesh().names();
132 
133  forAll(mesh_.boundaryMesh(), patchi)
134  {
135  patchNames_.append(mesh_.boundaryMesh()[patchi].name());
136 
137  OStringStream os;
138  mesh_.boundaryMesh()[patchi].write(os);
139  patchDicts_.append(dictionary(IStringStream(os.str())()));
140  }
141 
142  // Insert point, face and cell zones into the list
143 
144  // Point zones
145  wordList curPointZoneNames = mesh_.pointZones().toc();
146  if (curPointZoneNames.size())
147  {
148  pointZoneNames_.setCapacity(2*curPointZoneNames.size());
149  }
150 
151  forAll(curPointZoneNames, zonei)
152  {
153  pointZoneNames_.append(curPointZoneNames[zonei]);
154  }
155 
156  pointZonesAddedPoints_.setSize(pointZoneNames_.size());
157 
158  // Face zones
159  wordList curFaceZoneNames = mesh_.faceZones().toc();
160 
161  if (curFaceZoneNames.size())
162  {
163  faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
164  }
165 
166  forAll(curFaceZoneNames, zonei)
167  {
168  faceZoneNames_.append(curFaceZoneNames[zonei]);
169  }
170 
171  faceZonesAddedOrientedFaces_.setSize(faceZoneNames_.size());
172 
173  faceZonesAddedFaces_.setSize(faceZoneNames_.size());
174 
175  // Cell zones
176  wordList curCellZoneNames = mesh_.cellZones().toc();
177 
178  if (curCellZoneNames.size())
179  {
180  cellZoneNames_.setCapacity(2*curCellZoneNames.size());
181  }
182 
183  forAll(curCellZoneNames, zonei)
184  {
185  cellZoneNames_.append(curCellZoneNames[zonei]);
186  }
187 
188  cellZonesAddedCells_.setSize(cellZoneNames_.size());
189 }
190 
191 
192 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
193 
194 void Foam::mergePolyMesh::addMesh(const polyMesh& m)
195 {
196  // Add all the points, faces and cells of the new mesh
197 
198  // Add points
199 
200  const pointField& p = m.points();
201  labelList renumberPoints(p.size());
202 
203  const pointZoneList& pz = m.pointZones();
204  labelList pointZoneIndices(pz.size());
205 
206  forAll(pz, zonei)
207  {
208  pointZoneIndices[zonei] = zoneIndex(pointZoneNames_, pz[zonei].name());
209  pointZonesAddedPoints_.setSize(pointZoneNames_.size());
210  }
211 
212  forAll(p, pointi)
213  {
214  renumberPoints[pointi] = meshMod_.addPoint
215  (
216  p[pointi], // Point to add
217  -1, // Master point (straight addition)
218  pointi < m.nPoints() // Is in cell?
219  );
220 
221  const labelList zones(pz.whichZones(pointi));
222  forAll(zones, zonei)
223  {
224  pointZonesAddedPoints_[pointZoneIndices[zones[zonei]]]
225  .insert(renumberPoints[pointi]);
226  }
227  }
228 
229  // Add cells
230 
231  const cellList& c = m.cells();
232  labelList renumberCells(c.size());
233 
234  const cellZoneList& cz = m.cellZones();
235  labelList cellZoneIndices(cz.size());
236 
237  forAll(cz, zonei)
238  {
239  cellZoneIndices[zonei] = zoneIndex(cellZoneNames_, cz[zonei].name());
240  cellZonesAddedCells_.setSize(cellZoneNames_.size());
241  }
242 
243  forAll(c, celli)
244  {
245  renumberCells[celli] = meshMod_.addCell(-1);
246 
247  const labelList zones(cz.whichZones(celli));
248  forAll(zones, zonei)
249  {
250  cellZonesAddedCells_[cellZoneIndices[zones[zonei]]]
251  .insert(renumberCells[celli]);
252  }
253  }
254 
255  // Add faces
256  const polyBoundaryMesh& bm = m.boundaryMesh();
257 
258  // Gather the patch indices
259  labelList patchIndices(bm.size());
260 
261  forAll(patchIndices, patchi)
262  {
263  patchIndices[patchi] = patchIndex(bm[patchi]);
264  }
265 
266  // Temporary: update number of allowable patches. This should be
267  // determined at the top - before adding anything.
268  meshMod_.setNumPatches(patchNames_.size());
269 
270 
271 
272  const faceZoneList& fz = m.faceZones();
273  labelList faceZoneIndices(fz.size());
274 
275  forAll(fz, zonei)
276  {
277  faceZoneIndices[zonei] = zoneIndex(faceZoneNames_, fz[zonei].name());
278  }
279 
280  const faceList& f = m.faces();
281  labelList renumberFaces(f.size());
282 
283  const labelList& own = m.faceOwner();
284  const labelList& nei = m.faceNeighbour();
285 
286  label newOwn, newNei, newPatch;
287 
288  forAll(f, facei)
289  {
290  const face& curFace = f[facei];
291 
292  face newFace(curFace.size());
293 
294  forAll(curFace, pointi)
295  {
296  newFace[pointi] = renumberPoints[curFace[pointi]];
297  }
298 
299  if (debug)
300  {
301  // Check that the face is valid
302  if (min(newFace) < 0)
303  {
305  << "Error in point mapping for face " << facei
306  << ". Old face: " << curFace << " New face: " << newFace
307  << abort(FatalError);
308  }
309  }
310 
311  if (facei < m.nInternalFaces() || facei >= m.nFaces())
312  {
313  newPatch = -1;
314  }
315  else
316  {
317  newPatch = patchIndices[bm.whichPatch(facei)];
318  }
319 
320  newOwn = own[facei];
321  if (newOwn > -1) newOwn = renumberCells[newOwn];
322 
323  if (newPatch > -1)
324  {
325  newNei = -1;
326  }
327  else
328  {
329  newNei = nei[facei];
330  newNei = renumberCells[newNei];
331  }
332 
333  renumberFaces[facei] = meshMod_.addFace
334  (
335  newFace,
336  newOwn,
337  newNei,
338  -1,
339  false,
340  newPatch
341  );
342 
343  const labelList zones(fz.whichZones(facei));
344  forAll(zones, zonei)
345  {
346  const faceZone& fzi = fz[zones[zonei]];
347  const bool flip = fzi.flipMap()[fzi.localIndex(facei)];
348 
349  if (fzi.oriented())
350  {
351  faceZonesAddedOrientedFaces_[faceZoneIndices[zones[zonei]]]
352  .insert(renumberFaces[facei], flip);
353  }
354  else
355  {
356  faceZonesAddedFaces_[faceZoneIndices[zones[zonei]]]
357  .insert(renumberFaces[facei]);
358  }
359  }
360  }
361 }
362 
363 
365 {
366  if (debug)
367  {
368  Info<< "patch names: " << patchNames_ << nl
369  << "patch dicts: " << patchDicts_ << nl
370  << "point zone names: " << pointZoneNames_ << nl
371  << "face zone names: " << faceZoneNames_ << nl
372  << "cell zone names: " << cellZoneNames_ << endl;
373  }
374 
375  // Add the patches if necessary
376  if (patchNames_.size() != mesh_.boundaryMesh().size())
377  {
378  if (debug)
379  {
380  Info<< "Copying old patches" << endl;
381  }
382 
383  List<polyPatch*> newPatches(patchNames_.size());
384 
385  const polyBoundaryMesh& oldPatches = mesh_.boundaryMesh();
386 
387  // Note. Reusing counter in two for loops
388  label patchi = 0;
389 
390  for (patchi = 0; patchi < oldPatches.size(); patchi++)
391  {
392  newPatches[patchi] = oldPatches[patchi].clone(oldPatches).ptr();
393  }
394 
395  if (debug)
396  {
397  Info<< "Adding new patches. " << endl;
398  }
399 
400  label endOfLastPatch =
401  patchi == 0
402  ? 0
403  : oldPatches[patchi - 1].start() + oldPatches[patchi - 1].size();
404 
405  for (; patchi < patchNames_.size(); patchi++)
406  {
407  // Add a patch
408  dictionary dict(patchDicts_[patchi]);
409  dict.set("nFaces", 0);
410  dict.set("startFace", endOfLastPatch);
411 
412  newPatches[patchi] =
413  (
415  (
416  patchNames_[patchi],
417  dict,
418  patchi,
419  oldPatches
420  ).ptr()
421  );
422  }
423 
424  mesh_.removeBoundary();
425  mesh_.addPatches(newPatches);
426  }
427 
428  // Add the zones if necessary
429  if (pointZoneNames_.size() > mesh_.pointZones().size())
430  {
431  if (debug)
432  {
433  Info<< "Adding new pointZones. " << endl;
434  }
435 
436  label nZones = mesh_.pointZones().size();
437 
438  mesh_.pointZones().setSize(pointZoneNames_.size());
439 
440  for (label zonei = nZones; zonei < pointZoneNames_.size(); zonei++)
441  {
442  mesh_.pointZones().set
443  (
444  zonei,
445  new pointZone
446  (
447  pointZoneNames_[zonei],
448  labelList(),
449  mesh_.pointZones()
450  )
451  );
452  }
453  }
454 
455  if (cellZoneNames_.size() > mesh_.cellZones().size())
456  {
457  if (debug)
458  {
459  Info<< "Adding new cellZones. " << endl;
460  }
461 
462  label nZones = mesh_.cellZones().size();
463 
464  mesh_.cellZones().setSize(cellZoneNames_.size());
465 
466  for (label zonei = nZones; zonei < cellZoneNames_.size(); zonei++)
467  {
468  mesh_.cellZones().set
469  (
470  zonei,
471  new cellZone
472  (
473  cellZoneNames_[zonei],
474  labelList(),
475  mesh_.cellZones()
476  )
477  );
478  }
479  }
480 
481  if (faceZoneNames_.size() > mesh_.faceZones().size())
482  {
483  if (debug)
484  {
485  Info<< "Adding new faceZones. " << endl;
486  }
487 
488  label nZones = mesh_.faceZones().size();
489 
490  mesh_.faceZones().setSize(faceZoneNames_.size());
491 
492  for (label zonei = nZones; zonei < faceZoneNames_.size(); zonei++)
493  {
494  mesh_.faceZones().set
495  (
496  zonei,
497  new faceZone
498  (
499  faceZoneNames_[zonei],
500  labelList(),
501  boolList(),
502  mesh_.faceZones()
503  )
504  );
505  }
506  }
507 
508  // Change mesh
509  autoPtr<polyTopoChangeMap> map(meshMod_.changeMesh(mesh_));
510 
511  // Add the new points to the pointZones in the merged mesh
512  mesh_.pointZones().insert(pointZonesAddedPoints_);
513 
514  // Add the new oriented faces to the faceZones in the merged mesh
515  mesh_.faceZones().insert(faceZonesAddedOrientedFaces_);
516 
517  // Add the new faces to the faceZones in the merged mesh
518  mesh_.faceZones().insert(faceZonesAddedFaces_);
519 
520  // Add the new cells to the cellZones in the merged mesh
521  mesh_.cellZones().insert(cellZonesAddedCells_);
522 
523  mesh_.topoChange(map);
524 
525  // Clear topo change for the next operation
526  meshMod_.clear();
527 }
528 
529 
530 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
T & last()
Return the last element of the list.
Definition: UListI.H:128
void addMesh(const polyMesh &m)
Add a mesh.
mergePolyMesh(polyMesh &mesh)
Construct from polyMesh.
void merge()
Merge meshes.
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
Definition: polyPatchNew.C:32
void setNumPatches(const label nPatches)
Explicitly set the number of patches if construct-without-mesh.
label addCell(const label masterCellID)
Add cell and return new cell index.
label addPoint(const point &, const label masterPointID, const bool inCell)
Add point and return new point index.
label addFace(const face &f, const label own, const label nei, const label masterFaceID, const bool flipFaceFlux, const label patchID)
Add face to cells and return new face index.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
List< cell > cellList
list of cells
Definition: cellList.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
List< face > faceList
Definition: faceListFwd.H:41
static const char nl
Definition: Ostream.H:267
labelList f(nPoints)
dictionary dict
volScalarField & p