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 | Copyright (C) 2011-2016 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 "polyTopoChanger.H"
29 #include "mapPolyMesh.H"
30 #include "polyAddPoint.H"
31 #include "polyAddCell.H"
32 #include "polyAddFace.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 defineTypeNameAndDebug(mergePolyMesh, 1);
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::label Foam::mergePolyMesh::patchIndex(const polyPatch& p)
45 {
46  // Find the patch name on the list. If the patch is already there
47  // and patch types match, return index
48  const word& pType = p.type();
49  const word& pName = p.name();
50 
51  bool nameFound = false;
52 
53  forAll(patchNames_, patchi)
54  {
55  if (patchNames_[patchi] == pName)
56  {
57  if (word(patchDicts_[patchi]["type"]) == pType)
58  {
59  // Found name and types match
60  return patchi;
61  }
62  else
63  {
64  // Found the name, but type is different
65  nameFound = true;
66  }
67  }
68  }
69 
70  // Patch not found. Append to the list
71  {
72  OStringStream os;
73  p.write(os);
74  patchDicts_.append(dictionary(IStringStream(os.str())()));
75  }
76 
77  if (nameFound)
78  {
79  // Duplicate name is not allowed. Create a composite name from the
80  // patch name and case name
81  const word& caseName = p.boundaryMesh().mesh().time().caseName();
82 
83  patchNames_.append(pName + "_" + caseName);
84 
85  Info<< "label patchIndex(const polyPatch& p) : "
86  << "Patch " << p.index() << " named "
87  << pName << " in mesh " << caseName
88  << " already exists, but patch types "
89  << " do not match.\nCreating a composite name as "
90  << patchNames_.last() << endl;
91  }
92  else
93  {
94  patchNames_.append(pName);
95  }
96 
97  return patchNames_.size() - 1;
98 }
99 
100 
101 Foam::label Foam::mergePolyMesh::zoneIndex
102 (
103  DynamicList<word>& names,
104  const word& curName
105 )
106 {
107  forAll(names, zoneI)
108  {
109  if (names[zoneI] == curName)
110  {
111  return zoneI;
112  }
113  }
114 
115  // Not found. Add new name to the list
116  names.append(curName);
117 
118  return names.size() - 1;
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
124 Foam::mergePolyMesh::mergePolyMesh(const IOobject& io)
125 :
126  polyMesh(io),
127  meshMod_(*this),
128  patchNames_(2*boundaryMesh().size()),
129  patchDicts_(2*boundaryMesh().size()),
130  pointZoneNames_(),
131  faceZoneNames_(),
132  cellZoneNames_()
133 {
134  // Insert the original patches into the list
135  wordList curPatchNames = boundaryMesh().names();
136 
137  forAll(boundaryMesh(), patchi)
138  {
139  patchNames_.append(boundaryMesh()[patchi].name());
140 
141  OStringStream os;
142  boundaryMesh()[patchi].write(os);
143  patchDicts_.append(dictionary(IStringStream(os.str())()));
144  }
145 
146  // Insert point, face and cell zones into the list
147 
148  // Point zones
149  wordList curPointZoneNames = pointZones().names();
150  if (curPointZoneNames.size())
151  {
152  pointZoneNames_.setCapacity(2*curPointZoneNames.size());
153  }
154 
155  forAll(curPointZoneNames, zoneI)
156  {
157  pointZoneNames_.append(curPointZoneNames[zoneI]);
158  }
159 
160  // Face zones
161  wordList curFaceZoneNames = faceZones().names();
162 
163  if (curFaceZoneNames.size())
164  {
165  faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
166  }
167  forAll(curFaceZoneNames, zoneI)
168  {
169  faceZoneNames_.append(curFaceZoneNames[zoneI]);
170  }
171 
172  // Cell zones
173  wordList curCellZoneNames = cellZones().names();
174 
175  if (curCellZoneNames.size())
176  {
177  cellZoneNames_.setCapacity(2*curCellZoneNames.size());
178  }
179  forAll(curCellZoneNames, zoneI)
180  {
181  cellZoneNames_.append(curCellZoneNames[zoneI]);
182  }
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
191 void Foam::mergePolyMesh::addMesh(const polyMesh& m)
192 {
193  // Add all the points, faces and cells of the new mesh
194 
195  // Add points
196 
197  label zoneID = -1;
198 
199  const pointField& p = m.points();
200  labelList renumberPoints(p.size());
201 
202  const pointZoneMesh& pz = m.pointZones();
203  labelList pointZoneIndices(pz.size());
204 
205  forAll(pz, zoneI)
206  {
207  pointZoneIndices[zoneI] = zoneIndex(pointZoneNames_, pz[zoneI].name());
208  }
209 
210  forAll(p, pointi)
211  {
212  // Grab zone ID. If a point is not in a zone, it will return -1
213  zoneID = pz.whichZone(pointi);
214 
215  if (zoneID >= 0)
216  {
217  // Translate zone ID into the new index
218  zoneID = pointZoneIndices[zoneID];
219  }
220 
221  renumberPoints[pointi] =
222  meshMod_.setAction
223  (
224  polyAddPoint
225  (
226  p[pointi], // Point to add
227  -1, // Master point (straight addition)
228  zoneID, // Zone for point
229  pointi < m.nPoints() // Is in cell?
230  )
231  );
232  }
233 
234  // Add cells
235 
236  const cellList& c = m.cells();
237  labelList renumberCells(c.size());
238 
239  const cellZoneMesh& cz = m.cellZones();
240  labelList cellZoneIndices(cz.size());
241 
242  forAll(cz, zoneI)
243  {
244  cellZoneIndices[zoneI] = zoneIndex(cellZoneNames_, cz[zoneI].name());
245  }
246 
247  forAll(c, celli)
248  {
249  // Grab zone ID. If a cell is not in a zone, it will return -1
250  zoneID = cz.whichZone(celli);
251 
252  if (zoneID >= 0)
253  {
254  // Translate zone ID into the new index
255  zoneID = cellZoneIndices[zoneID];
256  }
257 
258  renumberCells[celli] =
259  meshMod_.setAction
260  (
261  polyAddCell
262  (
263  -1, // Master point
264  -1, // Master edge
265  -1, // Master face
266  -1, // Master cell
267  zoneID // Zone for cell
268  )
269  );
270  }
271 
272  // Add faces
273  const polyBoundaryMesh& bm = m.boundaryMesh();
274 
275  // Gather the patch indices
276  labelList patchIndices(bm.size());
277 
278  forAll(patchIndices, patchi)
279  {
280  patchIndices[patchi] = patchIndex(bm[patchi]);
281  }
282 
283  // Temporary: update number of allowable patches. This should be
284  // determined at the top - before adding anything.
285  meshMod_.setNumPatches(patchNames_.size());
286 
287 
288 
289  const faceZoneMesh& fz = m.faceZones();
290  labelList faceZoneIndices(fz.size());
291 
292  forAll(fz, zoneI)
293  {
294  faceZoneIndices[zoneI] = zoneIndex(faceZoneNames_, fz[zoneI].name());
295  }
296 
297  const faceList& f = m.faces();
298  labelList renumberFaces(f.size());
299 
300  const labelList& own = m.faceOwner();
301  const labelList& nei = m.faceNeighbour();
302 
303  label newOwn, newNei, newPatch, newZone;
304  bool newZoneFlip;
305 
306  forAll(f, facei)
307  {
308  const face& curFace = f[facei];
309 
310  face newFace(curFace.size());
311 
312  forAll(curFace, pointi)
313  {
314  newFace[pointi] = renumberPoints[curFace[pointi]];
315  }
316 
317  if (debug)
318  {
319  // Check that the face is valid
320  if (min(newFace) < 0)
321  {
323  << "Error in point mapping for face " << facei
324  << ". Old face: " << curFace << " New face: " << newFace
325  << abort(FatalError);
326  }
327  }
328 
329  if (facei < m.nInternalFaces() || facei >= m.nFaces())
330  {
331  newPatch = -1;
332  }
333  else
334  {
335  newPatch = patchIndices[bm.whichPatch(facei)];
336  }
337 
338  newOwn = own[facei];
339  if (newOwn > -1) newOwn = renumberCells[newOwn];
340 
341  if (newPatch > -1)
342  {
343  newNei = -1;
344  }
345  else
346  {
347  newNei = nei[facei];
348  newNei = renumberCells[newNei];
349  }
350 
351 
352  newZone = fz.whichZone(facei);
353  newZoneFlip = false;
354 
355  if (newZone >= 0)
356  {
357  newZoneFlip = fz[newZone].flipMap()[fz[newZone].whichFace(facei)];
358 
359  // Grab the new zone
360  newZone = faceZoneIndices[newZone];
361  }
362 
363  renumberFaces[facei] =
364  meshMod_.setAction
365  (
366  polyAddFace
367  (
368  newFace,
369  newOwn,
370  newNei,
371  -1,
372  -1,
373  -1,
374  false,
375  newPatch,
376  newZone,
377  newZoneFlip
378  )
379  );
380  }
381 
382 }
383 
384 
386 {
387  Info<< "patch names: " << patchNames_ << nl
388  << "patch dicts: " << patchDicts_ << nl
389  << "point zone names: " << pointZoneNames_ << nl
390  << "face zone names: " << faceZoneNames_ << nl
391  << "cell zone names: " << cellZoneNames_ << endl;
392 
393  // Add the patches if necessary
394  if (patchNames_.size() != boundaryMesh().size())
395  {
396  Info<< "Copying old patches" << endl;
397 
398  List<polyPatch*> newPatches(patchNames_.size());
399 
400  const polyBoundaryMesh& oldPatches = boundaryMesh();
401 
402  // Note. Re-using counter in two for loops
403  label patchi = 0;
404 
405  for (patchi = 0; patchi < oldPatches.size(); patchi++)
406  {
407  newPatches[patchi] = oldPatches[patchi].clone(oldPatches).ptr();
408  }
409 
410  Info<< "Adding new patches. " << endl;
411 
412  label endOfLastPatch =
413  oldPatches[patchi - 1].start() + oldPatches[patchi - 1].size();
414 
415  for (; patchi < patchNames_.size(); patchi++)
416  {
417  // Add a patch
418  dictionary dict(patchDicts_[patchi]);
419  dict.set("nFaces", 0);
420  dict.set("startFace", endOfLastPatch);
421 
422  newPatches[patchi] =
423  (
425  (
426  patchNames_[patchi],
427  dict,
428  patchi,
429  oldPatches
430  ).ptr()
431  );
432  }
433 
434  removeBoundary();
435  addPatches(newPatches);
436  }
437 
438  // Add the zones if necessary
439  if (pointZoneNames_.size() > pointZones().size())
440  {
441  Info<< "Adding new pointZones. " << endl;
442  label nZones = pointZones().size();
443 
444  pointZones().setSize(pointZoneNames_.size());
445 
446  for (label zoneI = nZones; zoneI < pointZoneNames_.size(); zoneI++)
447  {
448  pointZones().set
449  (
450  zoneI,
451  new pointZone
452  (
453  pointZoneNames_[zoneI],
454  labelList(),
455  zoneI,
456  pointZones()
457  )
458  );
459  }
460  }
461  if (cellZoneNames_.size() > cellZones().size())
462  {
463  Info<< "Adding new cellZones. " << endl;
464 
465  label nZones = cellZones().size();
466 
467  cellZones().setSize(cellZoneNames_.size());
468 
469  for (label zoneI = nZones; zoneI < cellZoneNames_.size(); zoneI++)
470  {
471  cellZones().set
472  (
473  zoneI,
474  new cellZone
475  (
476  cellZoneNames_[zoneI],
477  labelList(),
478  zoneI,
479  cellZones()
480  )
481  );
482  }
483  }
484  if (faceZoneNames_.size() > faceZones().size())
485  {
486  Info<< "Adding new faceZones. " << endl;
487 
488  label nZones = faceZones().size();
489 
490  faceZones().setSize(faceZoneNames_.size());
491 
492  for (label zoneI = nZones; zoneI < faceZoneNames_.size(); zoneI++)
493  {
494  faceZones().set
495  (
496  zoneI,
497  new faceZone
498  (
499  faceZoneNames_[zoneI],
500  labelList(),
501  boolList(),
502  zoneI,
503  faceZones()
504  )
505  );
506  }
507  }
508 
509  // Change mesh. No inflation
510  meshMod_.changeMesh(*this, false);
511 
512  // Clear topo change for the next operation
513  meshMod_.clear();
514 }
515 
516 
517 // ************************************************************************* //
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const dictionary & dict() const
const word & name() const
void addMesh(const polyMesh &m)
Add a mesh.
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< label > labelList
A List of labels.
Definition: labelList.H:56
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void merge()
Merge meshes.
static const char nl
Definition: Ostream.H:262
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:864
messageStream Info
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
List< cell > cellList
list of cells
Definition: cellList.H:42
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
Namespace for OpenFOAM.