Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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{
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 OCharStream os;
77 p.write(os);
78 ISpanStream is(os.view());
79
80 patchDicts_.push_back(dictionary(is));
81 }
82
83 if (nameFound)
84 {
85 // Duplicate name is not allowed. Create a composite name from the
86 // patch name and case name
87 const word& caseName = p.boundaryMesh().mesh().time().caseName();
88
89 patchNames_.append(pName + "_" + caseName);
90
91 Info<< "label patchIndex(const polyPatch& p) : "
92 << "Patch " << p.index() << " named "
93 << pName << " in mesh " << caseName
94 << " already exists, but patch types "
95 << " do not match.\nCreating a composite name as "
96 << patchNames_.last() << endl;
97 }
98 else
99 {
100 patchNames_.append(pName);
101 }
102
103 return patchNames_.size() - 1;
104}
105
106
107Foam::label Foam::mergePolyMesh::zoneIndex
108(
110 const word& curName
111)
112{
113 forAll(names, zonei)
114 {
115 if (names[zonei] == curName)
116 {
117 return zonei;
118 }
119 }
120
121 // Not found. Add new name to the list
122 names.append(curName);
123
124 return names.size() - 1;
125}
126
127
128void Foam::mergePolyMesh::sortProcessorPatches()
129{
130 Info<< "Reordering processor patches last" << endl;
131
132 // Updates boundaryMesh() and meshMod_ to guarantee processor patches
133 // are last. This could be done inside the merge() but it is far easier
134 // to do separately.
135
136
137 // 1. Shuffle the patches in the boundaryMesh
138
139 const polyBoundaryMesh& oldPatches = boundaryMesh();
140
141 polyPatchList newPatches(oldPatches.size());
142
143 labelList oldToSorted(oldPatches.size());
144
145 label nPatches = 0;
146
147 forAll(oldPatches, patchi)
148 {
149 const polyPatch& pp = oldPatches[patchi];
150
152 {
153 newPatches.set
154 (
155 nPatches,
156 pp.clone
157 (
158 oldPatches,
159 nPatches,
160 0,
161 nInternalFaces()
162 )
163 );
164
165 oldToSorted[patchi] = nPatches;
166 ++nPatches;
167 }
168 }
169
170 forAll(oldPatches, patchi)
171 {
172 const polyPatch& pp = oldPatches[patchi];
173
175 {
176 newPatches.set
177 (
178 nPatches,
179 pp.clone
180 (
181 oldPatches,
182 oldToSorted[patchi],
183 0,
184 nInternalFaces()
185 )
186 );
187
188 oldToSorted[patchi] = nPatches;
189 ++nPatches;
190 }
191 }
192
193
194 removeBoundary();
195 addPatches(newPatches);
196
197
198 // Update the polyTopoChange
200 (
201 meshMod_.region()
202 );
203
204 forAll(patchID, facei)
205 {
206 label patchi = patchID[facei];
207 if (patchi != -1)
208 {
209 patchID[facei] = oldToSorted[patchID[facei]];
210 }
211 }
212}
213
214
215// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
216
217Foam::mergePolyMesh::mergePolyMesh(const IOobject& io)
218:
219 polyMesh(io),
220 meshMod_(*this),
221 patchNames_(2*boundaryMesh().size()),
222 patchDicts_(2*boundaryMesh().size()),
223 pointZoneNames_(),
224 faceZoneNames_(),
225 cellZoneNames_()
226{
227 // Insert the original patches into the list
228 wordList curPatchNames = boundaryMesh().names();
229
230 OCharStream os;
231 forAll(boundaryMesh(), patchi)
232 {
233 os.rewind();
234 boundaryMesh()[patchi].write(os);
235 ISpanStream is(os.view());
236
237 patchNames_.push_back(boundaryMesh()[patchi].name());
238 patchDicts_.push_back(dictionary(is));
239 }
240
241 // Insert point, face and cell zones into the list
242
243 // Point zones
244 wordList curPointZoneNames = pointZones().names();
245 if (curPointZoneNames.size())
246 {
247 pointZoneNames_.setCapacity(2*curPointZoneNames.size());
248 pointZoneNames_.append(curPointZoneNames);
249 }
250
251 // Face zones
252 wordList curFaceZoneNames = faceZones().names();
253 if (curFaceZoneNames.size())
254 {
255 faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
256 faceZoneNames_.append(curFaceZoneNames);
257 }
258
259 // Cell zones
260 wordList curCellZoneNames = cellZones().names();
261 if (curCellZoneNames.size())
262 {
263 cellZoneNames_.setCapacity(2*curCellZoneNames.size());
264 cellZoneNames_.append(curCellZoneNames);
265 }
266}
267
268
269// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
270
271void Foam::mergePolyMesh::addMesh(const polyMesh& m)
272{
273 // Add all the points, faces and cells of the new mesh
274
275 // Add points
276
277 label zoneID = -1;
278
279 const pointField& p = m.points();
280 labelList renumberPoints(p.size());
281
282 const pointZoneMesh& pz = m.pointZones();
283 labelList pointZoneIndices(pz.size());
284
285 forAll(pz, zoneI)
286 {
287 pointZoneIndices[zoneI] = zoneIndex(pointZoneNames_, pz[zoneI].name());
288 }
289
290 forAll(p, pointi)
291 {
292 // Grab zone ID. If a point is not in a zone, it will return -1
293 zoneID = pz.whichZone(pointi);
294
295 if (zoneID >= 0)
296 {
297 // Translate zone ID into the new index
298 zoneID = pointZoneIndices[zoneID];
299 }
300
301 renumberPoints[pointi] =
302 meshMod_.setAction
303 (
304 polyAddPoint
305 (
306 p[pointi], // Point to add
307 -1, // Master point (straight addition)
308 zoneID, // Zone for point
309 pointi < m.nPoints() // Is in cell?
310 )
311 );
312 }
313
314 // Add cells
315
316 const cellList& c = m.cells();
317 labelList renumberCells(c.size());
318
319 const cellZoneMesh& cz = m.cellZones();
320 labelList cellZoneIndices(cz.size());
321
322 forAll(cz, zoneI)
323 {
324 cellZoneIndices[zoneI] = zoneIndex(cellZoneNames_, cz[zoneI].name());
325 }
326
327 forAll(c, celli)
328 {
329 // Grab zone ID. If a cell is not in a zone, it will return -1
330 zoneID = cz.whichZone(celli);
331
332 if (zoneID >= 0)
333 {
334 // Translate zone ID into the new index
335 zoneID = cellZoneIndices[zoneID];
336 }
337
338 renumberCells[celli] =
339 meshMod_.setAction
340 (
341 polyAddCell
342 (
343 -1, // Master point
344 -1, // Master edge
345 -1, // Master face
346 -1, // Master cell
347 zoneID // Zone for cell
348 )
349 );
350 }
351
352 // Add faces
353 const polyBoundaryMesh& bm = m.boundaryMesh();
354
355 // Gather the patch indices
356 labelList patchIndices(bm.size());
357
358 forAll(patchIndices, patchi)
359 {
360 patchIndices[patchi] = patchIndex(bm[patchi]);
361 }
362
363 // Temporary: update number of allowable patches. This should be
364 // determined at the top - before adding anything.
365 meshMod_.setNumPatches(patchNames_.size());
366
367
368
369 const faceZoneMesh& fz = m.faceZones();
370 labelList faceZoneIndices(fz.size());
371
372 forAll(fz, zoneI)
373 {
374 faceZoneIndices[zoneI] = zoneIndex(faceZoneNames_, fz[zoneI].name());
375 }
376
377 const faceList& f = m.faces();
378 labelList renumberFaces(f.size());
379
380 const labelList& own = m.faceOwner();
381 const labelList& nei = m.faceNeighbour();
382
383 label newOwn, newNei, newPatch, newZone;
384 bool newZoneFlip;
385
386 forAll(f, facei)
387 {
388 const face& curFace = f[facei];
389
390 face newFace(curFace.size());
391
392 forAll(curFace, pointi)
393 {
394 newFace[pointi] = renumberPoints[curFace[pointi]];
395 }
396
397 if (debug)
398 {
399 // Check that the face is valid
400 if (min(newFace) < 0)
401 {
403 << "Error in point mapping for face " << facei
404 << ". Old face: " << curFace << " New face: " << newFace
405 << abort(FatalError);
406 }
407 }
408
409 if (facei < m.nInternalFaces() || facei >= m.nFaces())
410 {
411 newPatch = -1;
412 }
413 else
414 {
415 newPatch = patchIndices[bm.whichPatch(facei)];
416 }
417
418 newOwn = own[facei];
419 if (newOwn > -1) newOwn = renumberCells[newOwn];
420
421 if (newPatch > -1)
422 {
423 newNei = -1;
424 }
425 else
426 {
427 newNei = nei[facei];
428 newNei = renumberCells[newNei];
429 }
430
431
432 newZone = fz.whichZone(facei);
433 newZoneFlip = false;
434
435 if (newZone >= 0)
436 {
437 newZoneFlip = fz[newZone].flipMap()[fz[newZone].whichFace(facei)];
438
439 // Grab the new zone
440 newZone = faceZoneIndices[newZone];
441 }
442
443 renumberFaces[facei] =
444 meshMod_.setAction
445 (
446 polyAddFace
447 (
448 newFace,
449 newOwn,
450 newNei,
451 -1,
452 -1,
453 -1,
454 false,
455 newPatch,
456 newZone,
457 newZoneFlip
458 )
459 );
460 }
461}
462
463
465{
466 Info<< "patch names: " << patchNames_ << nl
467 << "patch dicts: " << patchDicts_ << nl
468 << "point zone names: " << pointZoneNames_ << nl
469 << "face zone names: " << faceZoneNames_ << nl
470 << "cell zone names: " << cellZoneNames_ << endl;
471
472 // Add the patches if necessary
473 if (patchNames_.size() != boundaryMesh().size())
474 {
475 Info<< "Copying old patches" << endl;
476
477 polyPatchList newPatches(patchNames_.size());
478
479 const polyBoundaryMesh& oldPatches = boundaryMesh();
480
481 // Note. Re-using counter in two for loops
482 label patchi = 0;
483
484 for (patchi = 0; patchi < oldPatches.size(); patchi++)
485 {
486 newPatches.set
487 (
488 patchi,
489 oldPatches[patchi].clone(oldPatches)
490 );
491 }
492
493 Info<< "Adding new patches. " << endl;
494
495 label endOfLastPatch =
496 patchi == 0
497 ? 0
498 : oldPatches[patchi - 1].start() + oldPatches[patchi - 1].size();
499
500 for (; patchi < patchNames_.size(); patchi++)
501 {
502 // Add a patch
503 dictionary dict(patchDicts_[patchi]);
504 dict.set("nFaces", 0);
505 dict.set("startFace", endOfLastPatch);
506
507 newPatches.set
508 (
509 patchi,
510 (
512 (
513 patchNames_[patchi],
514 dict,
515 patchi,
516 oldPatches
517 )
518 )
519 );
520 }
521
522 removeBoundary();
523 addPatches(newPatches);
524 }
525
526 // Add the zones if necessary
527 if (pointZoneNames_.size() > pointZones().size())
528 {
529 Info<< "Adding new pointZones." << endl;
530
531 label zonei = pointZones().size(); // continue from here
532
533 const label nZones = pointZoneNames_.size();
534
535 pointZones().setSize(nZones);
536
537 for (/*nil*/; zonei < nZones; ++zonei)
538 {
539 pointZones().set
540 (
541 zonei,
542 new pointZone(pointZoneNames_[zonei], zonei, pointZones())
543 );
544 }
545 }
546 if (cellZoneNames_.size() > cellZones().size())
547 {
548 Info<< "Adding new cellZones." << endl;
549
550 label zonei = cellZones().size(); // continue from here
551
552 const label nZones = cellZoneNames_.size();
553
554 cellZones().setSize(cellZoneNames_.size());
555
556 for (/*nil*/; zonei < nZones; ++zonei)
557 {
558 cellZones().set
559 (
560 zonei,
561 new cellZone(cellZoneNames_[zonei], zonei, cellZones())
562 );
563 }
564 }
565 if (faceZoneNames_.size() > faceZones().size())
566 {
567 Info<< "Adding new faceZones." << endl;
568
569 label zonei = faceZones().size(); // continue from here
570
571 const label nZones = faceZoneNames_.size();
572
573 faceZones().setSize(nZones);
574
575 for (/*nil*/; zonei < nZones; ++zonei)
576 {
577 faceZones().set
578 (
579 zonei,
580 new faceZone(faceZoneNames_[zonei], zonei, faceZones())
581 );
582 }
583 }
584
585
586 // Shuffle the processor patches to be last
587 sortProcessorPatches();
588
589 // Change mesh. No inflation
590 meshMod_.changeMesh(*this, false);
591
592 // Clear topo change for the next operation
593 meshMod_.clear();
594}
595
596
597// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const fileName & caseName() const noexcept
Return the Time::caseName().
Definition IOobject.C:468
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A subset of mesh cells.
Definition cellZone.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
Add a given mesh to the original mesh to create a single new mesh.
void addMesh(const polyMesh &m)
Add a mesh.
void merge()
Merge meshes.
A subset of mesh points.
Definition pointZone.H:64
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
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.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const auto & io
auto & name
auto & names
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition polyPatch.H:61
List< word > wordList
List of word.
Definition fileName.H:60
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with pointZone content on a polyMesh.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
List< cell > cellList
List of cell.
Definition cellListFwd.H:41
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with cellZone content on a polyMesh.
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
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
label nPatches
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299