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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
40namespace Foam
41{
42 defineTypeNameAndDebug(mergePolyMesh, 1);
43}
44
45
46// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47
48Foam::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
105Foam::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
126void 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 DynamicList<polyPatch*> newPatches(oldPatches.size());
140
141 labelList oldToSorted(oldPatches.size());
142
143 forAll(oldPatches, patchi)
144 {
145 const polyPatch& pp = oldPatches[patchi];
146
147 if (!isA<processorPolyPatch>(pp))
148 {
149 oldToSorted[patchi] = newPatches.size();
150 newPatches.append
151 (
152 pp.clone
153 (
154 oldPatches,
155 oldToSorted[patchi],
156 0,
157 nInternalFaces()
158 ).ptr()
159 );
160 }
161 }
162 forAll(oldPatches, patchi)
163 {
164 const polyPatch& pp = oldPatches[patchi];
165
166 if (isA<processorPolyPatch>(pp))
167 {
168 oldToSorted[patchi] = newPatches.size();
169 newPatches.append
170 (
171 pp.clone
172 (
173 oldPatches,
174 oldToSorted[patchi],
175 0,
176 nInternalFaces()
177 ).ptr()
178 );
179 }
180 }
181
182
183 removeBoundary();
184 addPatches(newPatches);
185
186
187 // Update the polyTopoChange
188 DynamicList<label>& patchID = const_cast<DynamicList<label>&>
189 (
190 meshMod_.region()
191 );
192
193 forAll(patchID, facei)
194 {
195 label patchi = patchID[facei];
196 if (patchi != -1)
197 {
198 patchID[facei] = oldToSorted[patchID[facei]];
199 }
200 }
201}
202
203
204// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
205
207:
208 polyMesh(io),
209 meshMod_(*this),
210 patchNames_(2*boundaryMesh().size()),
211 patchDicts_(2*boundaryMesh().size()),
212 pointZoneNames_(),
213 faceZoneNames_(),
214 cellZoneNames_()
215{
216 // Insert the original patches into the list
217 wordList curPatchNames = boundaryMesh().names();
218
219 forAll(boundaryMesh(), patchi)
220 {
221 patchNames_.append(boundaryMesh()[patchi].name());
222
223 OStringStream os;
224 boundaryMesh()[patchi].write(os);
225 patchDicts_.append(dictionary(IStringStream(os.str())()));
226 }
227
228 // Insert point, face and cell zones into the list
229
230 // Point zones
231 wordList curPointZoneNames = pointZones().names();
232 if (curPointZoneNames.size())
233 {
234 pointZoneNames_.setCapacity(2*curPointZoneNames.size());
235 pointZoneNames_.append(curPointZoneNames);
236 }
237
238 // Face zones
239 wordList curFaceZoneNames = faceZones().names();
240 if (curFaceZoneNames.size())
241 {
242 faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
243 faceZoneNames_.append(curFaceZoneNames);
244 }
245
246 // Cell zones
247 wordList curCellZoneNames = cellZones().names();
248 if (curCellZoneNames.size())
249 {
250 cellZoneNames_.setCapacity(2*curCellZoneNames.size());
251 cellZoneNames_.append(curCellZoneNames);
252 }
253}
254
255
256// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
257
258void Foam::mergePolyMesh::addMesh(const polyMesh& m)
259{
260 // Add all the points, faces and cells of the new mesh
261
262 // Add points
263
264 label zoneID = -1;
265
266 const pointField& p = m.points();
267 labelList renumberPoints(p.size());
268
269 const pointZoneMesh& pz = m.pointZones();
270 labelList pointZoneIndices(pz.size());
271
272 forAll(pz, zoneI)
273 {
274 pointZoneIndices[zoneI] = zoneIndex(pointZoneNames_, pz[zoneI].name());
275 }
276
277 forAll(p, pointi)
278 {
279 // Grab zone ID. If a point is not in a zone, it will return -1
280 zoneID = pz.whichZone(pointi);
281
282 if (zoneID >= 0)
283 {
284 // Translate zone ID into the new index
285 zoneID = pointZoneIndices[zoneID];
286 }
287
288 renumberPoints[pointi] =
289 meshMod_.setAction
290 (
291 polyAddPoint
292 (
293 p[pointi], // Point to add
294 -1, // Master point (straight addition)
295 zoneID, // Zone for point
296 pointi < m.nPoints() // Is in cell?
297 )
298 );
299 }
300
301 // Add cells
302
303 const cellList& c = m.cells();
304 labelList renumberCells(c.size());
305
306 const cellZoneMesh& cz = m.cellZones();
307 labelList cellZoneIndices(cz.size());
308
309 forAll(cz, zoneI)
310 {
311 cellZoneIndices[zoneI] = zoneIndex(cellZoneNames_, cz[zoneI].name());
312 }
313
314 forAll(c, celli)
315 {
316 // Grab zone ID. If a cell is not in a zone, it will return -1
317 zoneID = cz.whichZone(celli);
318
319 if (zoneID >= 0)
320 {
321 // Translate zone ID into the new index
322 zoneID = cellZoneIndices[zoneID];
323 }
324
325 renumberCells[celli] =
326 meshMod_.setAction
327 (
328 polyAddCell
329 (
330 -1, // Master point
331 -1, // Master edge
332 -1, // Master face
333 -1, // Master cell
334 zoneID // Zone for cell
335 )
336 );
337 }
338
339 // Add faces
340 const polyBoundaryMesh& bm = m.boundaryMesh();
341
342 // Gather the patch indices
343 labelList patchIndices(bm.size());
344
345 forAll(patchIndices, patchi)
346 {
347 patchIndices[patchi] = patchIndex(bm[patchi]);
348 }
349
350 // Temporary: update number of allowable patches. This should be
351 // determined at the top - before adding anything.
352 meshMod_.setNumPatches(patchNames_.size());
353
354
355
356 const faceZoneMesh& fz = m.faceZones();
357 labelList faceZoneIndices(fz.size());
358
359 forAll(fz, zoneI)
360 {
361 faceZoneIndices[zoneI] = zoneIndex(faceZoneNames_, fz[zoneI].name());
362 }
363
364 const faceList& f = m.faces();
365 labelList renumberFaces(f.size());
366
367 const labelList& own = m.faceOwner();
368 const labelList& nei = m.faceNeighbour();
369
370 label newOwn, newNei, newPatch, newZone;
371 bool newZoneFlip;
372
373 forAll(f, facei)
374 {
375 const face& curFace = f[facei];
376
377 face newFace(curFace.size());
378
379 forAll(curFace, pointi)
380 {
381 newFace[pointi] = renumberPoints[curFace[pointi]];
382 }
383
384 if (debug)
385 {
386 // Check that the face is valid
387 if (min(newFace) < 0)
388 {
390 << "Error in point mapping for face " << facei
391 << ". Old face: " << curFace << " New face: " << newFace
392 << abort(FatalError);
393 }
394 }
395
396 if (facei < m.nInternalFaces() || facei >= m.nFaces())
397 {
398 newPatch = -1;
399 }
400 else
401 {
402 newPatch = patchIndices[bm.whichPatch(facei)];
403 }
404
405 newOwn = own[facei];
406 if (newOwn > -1) newOwn = renumberCells[newOwn];
407
408 if (newPatch > -1)
409 {
410 newNei = -1;
411 }
412 else
413 {
414 newNei = nei[facei];
415 newNei = renumberCells[newNei];
416 }
417
418
419 newZone = fz.whichZone(facei);
420 newZoneFlip = false;
421
422 if (newZone >= 0)
423 {
424 newZoneFlip = fz[newZone].flipMap()[fz[newZone].whichFace(facei)];
425
426 // Grab the new zone
427 newZone = faceZoneIndices[newZone];
428 }
429
430 renumberFaces[facei] =
431 meshMod_.setAction
432 (
433 polyAddFace
434 (
435 newFace,
436 newOwn,
437 newNei,
438 -1,
439 -1,
440 -1,
441 false,
442 newPatch,
443 newZone,
444 newZoneFlip
445 )
446 );
447 }
448}
449
450
452{
453 Info<< "patch names: " << patchNames_ << nl
454 << "patch dicts: " << patchDicts_ << nl
455 << "point zone names: " << pointZoneNames_ << nl
456 << "face zone names: " << faceZoneNames_ << nl
457 << "cell zone names: " << cellZoneNames_ << endl;
458
459 // Add the patches if necessary
460 if (patchNames_.size() != boundaryMesh().size())
461 {
462 Info<< "Copying old patches" << endl;
463
464 List<polyPatch*> newPatches(patchNames_.size());
465
466 const polyBoundaryMesh& oldPatches = boundaryMesh();
467
468 // Note. Re-using counter in two for loops
469 label patchi = 0;
470
471 for (patchi = 0; patchi < oldPatches.size(); patchi++)
472 {
473 newPatches[patchi] = oldPatches[patchi].clone(oldPatches).ptr();
474 }
475
476 Info<< "Adding new patches. " << endl;
477
478 label endOfLastPatch =
479 patchi == 0
480 ? 0
481 : oldPatches[patchi - 1].start() + oldPatches[patchi - 1].size();
482
483 for (; patchi < patchNames_.size(); patchi++)
484 {
485 // Add a patch
486 dictionary dict(patchDicts_[patchi]);
487 dict.set("nFaces", 0);
488 dict.set("startFace", endOfLastPatch);
489
490 newPatches[patchi] =
491 (
493 (
494 patchNames_[patchi],
495 dict,
496 patchi,
497 oldPatches
498 ).ptr()
499 );
500 }
501
502 removeBoundary();
503 addPatches(newPatches);
504 }
505
506 // Add the zones if necessary
507 if (pointZoneNames_.size() > pointZones().size())
508 {
509 Info<< "Adding new pointZones." << endl;
510
511 label zonei = pointZones().size(); // continue from here
512
513 const label nZones = pointZoneNames_.size();
514
515 pointZones().setSize(nZones);
516
517 for (/*nil*/; zonei < nZones; ++zonei)
518 {
519 pointZones().set
520 (
521 zonei,
522 new pointZone(pointZoneNames_[zonei], zonei, pointZones())
523 );
524 }
525 }
526 if (cellZoneNames_.size() > cellZones().size())
527 {
528 Info<< "Adding new cellZones." << endl;
529
530 label zonei = cellZones().size(); // continue from here
531
532 const label nZones = cellZoneNames_.size();
533
534 cellZones().setSize(cellZoneNames_.size());
535
536 for (/*nil*/; zonei < nZones; ++zonei)
537 {
538 cellZones().set
539 (
540 zonei,
541 new cellZone(cellZoneNames_[zonei], zonei, cellZones())
542 );
543 }
544 }
545 if (faceZoneNames_.size() > faceZones().size())
546 {
547 Info<< "Adding new faceZones." << endl;
548
549 label zonei = faceZones().size(); // continue from here
550
551 const label nZones = faceZoneNames_.size();
552
553 faceZones().setSize(nZones);
554
555 for (/*nil*/; zonei < nZones; ++zonei)
556 {
557 faceZones().set
558 (
559 zonei,
560 new faceZone(faceZoneNames_[zonei], zonei, faceZones())
561 );
562 }
563 }
564
565
566 // Shuffle the processor patches to be last
567 sortProcessorPatches();
568
569 // Change mesh. No inflation
570 meshMod_.changeMesh(*this, false);
571
572 // Clear topo change for the next operation
573 meshMod_.clear();
574}
575
576
577// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
const fileName & caseName() const
Return the Time::caseName()
Definition: IOobject.C:518
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
const fileName & caseName() const
Return case name.
Definition: TimePathsI.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:216
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:780
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
Add a given mesh to the original mesh to create a single new mesh.
Definition: mergePolyMesh.H:57
void addMesh(const polyMesh &m)
Add a mesh.
void merge()
Merge meshes.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
void setNumPatches(const label nPatches)
label patchIndex() const
Return access to the patch index.
virtual bool write(const bool valid=true) const
Write using setting from DB.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
label nZones
const labelIOList & zoneID
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
List< label > labelList
A List of labels.
Definition: List.H:66
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
List< cell > cellList
A List of cells.
Definition: cellListFwd.H:47
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333