Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
reconstructParMesh
Description
Chris Greenshields
committed
Reconstructs a mesh using geometric information only.
Writes point/face/cell procAddressing so afterwards reconstructPar can be
used to reconstruct fields.
- uses geometric matching tolerance (set with -mergeTol (at your option)
If the parallel case does not have correct procBoundaries use the
-fullMatch option which will check all boundary faces (bit slower).
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "timeSelector.H"
#include "IOobjectList.H"
#include "labelIOList.H"
#include "processorPolyPatch.H"
#include "mapAddedPolyMesh.H"
#include "polyMeshAdder.H"
#include "faceCoupleInfo.H"
#include "zeroGradientFvPatchFields.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Tolerance (as fraction of the bounding box). Needs to be fairly lax since
// usually meshes get written with limited precision (6 digits)
static const scalar defaultMergeTol = 1e-7;
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
static void renumber
(
const labelList& map,
labelList& elems
)
{
forAll(elems, i)
{
if (elems[i] >= 0)
{
elems[i] = map[elems[i]];
}
}
}
// Determine which faces are coupled. Uses geometric merge distance.
// Looks either at all boundaryFaces (fullMatch) or only at the
// procBoundaries for procI. Assumes that masterMesh contains already merged
// all the processors < procI.
autoPtr<faceCoupleInfo> determineCoupledFaces
(
const bool fullMatch,
const label procI,
const polyMesh& masterMesh,
const polyMesh& meshToAdd,
const scalar mergeDist
)
{
if (fullMatch || masterMesh.nCells() == 0)
{
return autoPtr<faceCoupleInfo>
(
new faceCoupleInfo
(
masterMesh,
meshToAdd,
mergeDist, // Absolute merging distance
true // Matching faces identical
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
)
);
}
else
{
// Pick up all patches on masterMesh ending in "toDDD" where DDD is
// the processor number procI.
const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
const string toProcString("to" + name(procI));
DynamicList<label> masterFaces
(
masterMesh.nFaces()
- masterMesh.nInternalFaces()
);
forAll(masterPatches, patchI)
{
const polyPatch& pp = masterPatches[patchI];
if
(
isA<processorPolyPatch>(pp)
&& (
pp.name().rfind(toProcString)
== (pp.name().size()-toProcString.size())
)
)
{
label meshFaceI = pp.start();
forAll(pp, i)
{
masterFaces.append(meshFaceI++);
}
}
}
masterFaces.shrink();
// Pick up all patches on meshToAdd ending in "procBoundaryDDDtoYYY"
// where DDD is the processor number procI and YYY is < procI.
const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
DynamicList<label> addFaces
(
meshToAdd.nFaces()
- meshToAdd.nInternalFaces()
);
forAll(addPatches, patchI)
{
const polyPatch& pp = addPatches[patchI];
if (isA<processorPolyPatch>(pp))
{
bool isConnected = false;
for (label mergedProcI = 0; mergedProcI < procI; mergedProcI++)
{
const string fromProcString
(
"procBoundary"
+ name(procI)
+ "to"
+ name(mergedProcI)
);
if (pp.name() == fromProcString)
{
isConnected = true;
break;
}
}
if (isConnected)
{
label meshFaceI = pp.start();
forAll(pp, i)
{
addFaces.append(meshFaceI++);
}
}
}
}
addFaces.shrink();
return autoPtr<faceCoupleInfo>
(
new faceCoupleInfo
(
masterMesh,
masterFaces,
meshToAdd,
addFaces,
mergeDist, // Absolute merging distance
true, // Matching faces identical?
false, // If perfect match are faces already ordered
// (e.g. processor patches)
false // are faces each on separate patch?
)
);
}
}
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
autoPtr<mapPolyMesh> mergeSharedPoints
(
const scalar mergeDist,
polyMesh& mesh,
labelListList& pointProcAddressing
)
{
// Find out which sets of points get merged and create a map from
// mesh point to unique point.
Map<label> pointToMaster
(
fvMeshAdder::findSharedPoints
(
mesh,
mergeDist
)
);
Info<< "mergeSharedPoints : detected " << pointToMaster.size()
<< " points that are to be merged." << endl;
if (returnReduce(pointToMaster.size(), sumOp<label>()) == 0)
{
return autoPtr<mapPolyMesh>(NULL);
}
polyTopoChange meshMod(mesh);
fvMeshAdder::mergePoints(mesh, pointToMaster, meshMod);
// Change the mesh (no inflation). Note: parallel comms allowed.
autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
// Update fields. No inflation, parallel sync.
mesh.updateMesh(map);
// pointProcAddressing give indices into the master mesh so adapt them
// for changed point numbering.
// Adapt constructMaps for merged points.
forAll(pointProcAddressing, procI)
{
labelList& constructMap = pointProcAddressing[procI];
forAll(constructMap, i)
{
label oldPointI = constructMap[i];
// New label of point after changeMesh.
label newPointI = map().reversePointMap()[oldPointI];
if (newPointI < -1)
{
constructMap[i] = -newPointI-2;
}
else if (newPointI >= 0)
{
constructMap[i] = newPointI;
}
else
{
FatalErrorIn("fvMeshDistribute::mergeSharedPoints()")
<< "Problem. oldPointI:" << oldPointI
<< " newPointI:" << newPointI << abort(FatalError);
}
}
}
return map;
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
boundBox procBounds
(
const argList& args,
const PtrList<Time>& databases,
const word& regionDir
)
{
boundBox bb = boundBox::invertedBox;
forAll(databases, procI)
{
fileName pointsInstance
(
databases[procI].findInstance
(
regionDir/polyMesh::meshSubDir,
"points"
)
);
if (pointsInstance != databases[procI].timeName())
{
FatalErrorIn(args.executable())
<< "Your time was specified as " << databases[procI].timeName()
<< " but there is no polyMesh/points in that time." << endl
<< "(there is a points file in " << pointsInstance
<< ")" << endl
<< "Please rerun with the correct time specified"
<< " (through the -constant, -time or -latestTime "
<< "(at your option)."
<< endl << exit(FatalError);
}
Info<< "Reading points from "
<< databases[procI].caseName()
<< " for time = " << databases[procI].timeName()
<< nl << endl;
pointIOField points
(
IOobject
(
"points",
databases[procI].findInstance
(
regionDir/polyMesh::meshSubDir,
"points"
),
regionDir/polyMesh::meshSubDir,
databases[procI],
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
)
);
boundBox domainBb(points, false);
bb.min() = min(bb.min(), domainBb.min());
bb.max() = max(bb.max(), domainBb.max());
}
return bb;
}
void writeCellDistance
(
Time& runTime,
const fvMesh& masterMesh,
const labelListList& cellProcAddressing
)
{
// Write the decomposition as labelList for use with 'manual'
// decomposition method.
labelIOList cellDecomposition
(
IOobject
(
"cellDecomposition",
masterMesh.facesInstance(),
masterMesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
masterMesh.nCells()
);
forAll(cellProcAddressing, procI)
{
const labelList& pCells = cellProcAddressing[procI];
UIndirectList<label>(cellDecomposition, pCells) = procI;
}
cellDecomposition.write();
Info<< nl << "Wrote decomposition to "
<< cellDecomposition.objectPath()
<< " for use in manual decomposition." << endl;
// Write as volScalarField for postprocessing. Change time to 0
// if was 'constant'
{
const scalar oldTime = runTime.value();
const label oldIndex = runTime.timeIndex();
if (runTime.timeName() == runTime.constant() && oldIndex == 0)
{
runTime.setTime(0, oldIndex+1);
}
volScalarField cellDist
(
IOobject
(
"cellDist",
runTime.timeName(),
masterMesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
masterMesh,
dimensionedScalar("cellDist", dimless, 0),
zeroGradientFvPatchScalarField::typeName
);
forAll(cellDecomposition, cellI)
{
cellDist[cellI] = cellDecomposition[cellI];
}
cellDist.write();
Info<< nl << "Wrote decomposition as volScalarField to "
<< cellDist.name() << " for use in postprocessing."
<< endl;
// Restore time
runTime.setTime(oldTime, oldIndex);
}
}
int main(int argc, char *argv[])
{
argList::addNote
(
"reconstruct a mesh using geometric information only"
);
// Enable -constant ... if someone really wants it
// Enable -zeroTime to prevent accidentally trashing the initial fields
timeSelector::addOptions(true, true);
argList::addOption
(
"mergeTol",
"scalar",
"specify the merge distance relative to the bounding box size "
);
argList::addBoolOption
(
"fullMatch",
"do (slower) geometric matching on all boundary faces"
);
argList::addBoolOption
(
"cellDist",
"write cell distribution as a labelList - for use with 'manual' "
"decomposition method or as a volScalarField for post-processing."
);
#include "addTimeOptions.H"
#include "addRegionOption.H"
#include "setRootCase.H"
#include "createTime.H"
<< " individual processor" << nl
<< "meshes back into one master mesh. Use it if the original"
<< " master mesh has" << nl
<< "been deleted or if the processor meshes have been modified"
<< " (topology change)." << nl
<< "This tool will write the resulting mesh to a new time step"
<< " and construct" << nl
<< "xxxxProcAddressing files in the processor meshes so"
<< " reconstructPar can be" << nl
<< "used to regenerate the fields on the master mesh." << nl
<< nl
<< "Not well tested & use at your own risk!" << nl
<< endl;
word regionName = polyMesh::defaultRegion;
word regionDir = word::null;
if
(
args.optionReadIfPresent("region", regionName)
&& regionName != polyMesh::defaultRegion
)
regionDir = regionName;
Info<< "Operating on region " << regionName << nl << endl;
}
scalar mergeTol = defaultMergeTol;
args.optionReadIfPresent("mergeTol", mergeTol);
scalar writeTol = Foam::pow(10.0, -scalar(IOstream::defaultPrecision()));
<< "Write tolerance : " << writeTol << endl;
if (runTime.writeFormat() == IOstream::ASCII && mergeTol < writeTol)
{
FatalErrorIn(args.executable())
<< "Your current settings specify ASCII writing with "
<< IOstream::defaultPrecision() << " digits precision." << endl
<< "Your merging tolerance (" << mergeTol << ") is finer than this."
<< endl
<< "Please change your writeFormat to binary"
<< " or increase the writePrecision" << endl
<< "or adjust the merge tolerance (-mergeTol)."
<< exit(FatalError);
}
const bool fullMatch = args.optionFound("fullMatch");
if (fullMatch)
{
Info<< "Doing geometric matching on all boundary faces." << nl << endl;
Info<< "Doing geometric matching on correct procBoundaries only."
<< nl << "This assumes a correct decomposition." << endl;
}
bool writeCellDist = args.optionFound("cellDist");
int nProcs = 0;
while
(
(
args.rootPath()
/ args.caseName()
/ fileName(word("processor") + name(nProcs))
)
)
{
nProcs++;
}
Info<< "Found " << nProcs << " processor directories" << nl << endl;
// Read all time databases
PtrList<Time> databases(nProcs);
forAll(databases, procI)
<< args.caseName()/fileName(word("processor") + name(procI))
<< endl;
databases.set
(
procI,
new Time
(
Time::controlDictName,
args.rootPath(),
args.caseName()/fileName(word("processor") + name(procI))
)
);
// Use the times list from the master processor
// and select a subset based on the command-line options
instantList Times = timeSelector::select
(
databases[0].times(),
args
);
// Set startTime and endTime depending on -time and -latestTime options
#include "checkTimeOptions.H"
if (Times.empty())
{
FatalErrorIn(args.executable())
<< "No times selected"
<< exit(FatalError);
// Loop over all times
for (label timeI = startTime; timeI < endTime; timeI++)
// Set time for global database
runTime.setTime(Times[timeI], timeI);
Info<< "Time = " << runTime.timeName() << endl << endl;
// Set time for all databases
forAll(databases, procI)
databases[procI].setTime(Times[timeI], timeI);
const fileName meshPath =
databases[0].path()
/databases[0].timeName()
/regionDir
/polyMesh::meshSubDir;
if (!isFile(meshPath/"faces"))
{
Info<< "No mesh." << nl << endl;
continue;
}
// Read point on individual processors to determine merge tolerance
// (otherwise single cell domains might give problems)
const boundBox bb = procBounds(args, databases, regionDir);
const scalar mergeDist = mergeTol*bb.mag();
Info<< "Overall mesh bounding box : " << bb << nl
<< "Relative tolerance : " << mergeTol << nl
<< "Absolute matching distance : " << mergeDist << nl
<< endl;
// Addressing from processor to reconstructed case
labelListList cellProcAddressing(nProcs);
labelListList faceProcAddressing(nProcs);
labelListList pointProcAddressing(nProcs);
labelListList boundaryProcAddressing(nProcs);
// Internal faces on the final reconstructed mesh
label masterInternalFaces;
// Owner addressing on the final reconstructed mesh
labelList masterOwner;
// Construct empty mesh.
Info<< "Constructing empty mesh to add to." << nl << endl;
fvMesh masterMesh
(
IOobject
(
regionName,
runTime.timeName(),
runTime,
IOobject::NO_READ
),
xferCopy(pointField()),
xferCopy(faceList()),
xferCopy(cellList())
for (label procI = 0; procI < nProcs; procI++)
Info<< "Reading mesh to add from "
<< databases[procI].caseName()
<< " for time = " << databases[procI].timeName()
<< nl << endl;
fvMesh meshToAdd
(
IOobject
(
regionName,
databases[procI].timeName(),
databases[procI]
)
);
// Initialize its addressing
cellProcAddressing[procI] = identity(meshToAdd.nCells());
faceProcAddressing[procI] = identity(meshToAdd.nFaces());
pointProcAddressing[procI] = identity(meshToAdd.nPoints());
boundaryProcAddressing[procI] =
identity(meshToAdd.boundaryMesh().size());
// Find geometrically shared points/faces.
autoPtr<faceCoupleInfo> couples = determineCoupledFaces
(
fullMatch,
procI,
masterMesh,
meshToAdd,
mergeDist
);
// Add elements to mesh
Info<< "Adding to master mesh" << nl << endl;
autoPtr<mapAddedPolyMesh> map = fvMeshAdder::add
(
masterMesh,
meshToAdd,
couples
);
// Update all addressing so xxProcAddressing points to correct
// item in masterMesh.
// Processors that were already in masterMesh
for (label mergedI = 0; mergedI < procI; mergedI++)
{
renumber(map().oldCellMap(), cellProcAddressing[mergedI]);
renumber(map().oldFaceMap(), faceProcAddressing[mergedI]);
renumber(map().oldPointMap(), pointProcAddressing[mergedI]);
// Note: boundary is special since can contain -1.
renumber
(
map().oldPatchMap(),
boundaryProcAddressing[mergedI]
);
}
// Added processor
renumber(map().addedCellMap(), cellProcAddressing[procI]);
renumber(map().addedFaceMap(), faceProcAddressing[procI]);
renumber(map().addedPointMap(), pointProcAddressing[procI]);
renumber(map().addedPatchMap(), boundaryProcAddressing[procI]);
Info<< endl;
}
// See if any points on the mastermesh have become connected
// because of connections through processor meshes.
mergeSharedPoints(mergeDist, masterMesh, pointProcAddressing);
// Save some properties on the reconstructed mesh
masterInternalFaces = masterMesh.nInternalFaces();
masterOwner = masterMesh.faceOwner();
Info<< "\nWriting merged mesh to "
<< runTime.path()/runTime.timeName()
<< nl << endl;
if (!masterMesh.write())
{
FatalErrorIn(args.executable())
<< "Failed writing polyMesh."
<< exit(FatalError);
}
if (writeCellDist)
{
writeCellDistance(runTime, masterMesh, cellProcAddressing);
// Write the addressing
Info<< "Reconstructing the addressing from the processor meshes"
<< " to the newly reconstructed mesh" << nl << endl;
forAll(databases, procI)
{
Info<< "Reading processor " << procI << " mesh from "
<< databases[procI].caseName() << endl;
polyMesh procMesh
IOobject
(
regionName,
databases[procI].timeName(),
databases[procI]
)
);
// From processor point to reconstructed mesh point
Info<< "Writing pointProcAddressing to "
<< databases[procI].caseName()
/procMesh.facesInstance()
/polyMesh::meshSubDir
<< endl;
labelIOList
IOobject
(
"pointProcAddressing",
procMesh.facesInstance(),
polyMesh::meshSubDir,
procMesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false // Do not register
),
pointProcAddressing[procI]
).write();
// From processor face to reconstructed mesh face
Info<< "Writing faceProcAddressing to "
<< databases[procI].caseName()
/procMesh.facesInstance()
/polyMesh::meshSubDir
<< endl;
labelIOList faceProcAddr
IOobject
(
"faceProcAddressing",
procMesh.facesInstance(),
polyMesh::meshSubDir,
procMesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false // Do not register
),
faceProcAddressing[procI]
);
// Now add turning index to faceProcAddressing.
// See reconstructPar for meaning of turning index.
forAll(faceProcAddr, procFaceI)
label masterFaceI = faceProcAddr[procFaceI];
if
(
!procMesh.isInternalFace(procFaceI)
&& masterFaceI < masterInternalFaces
)
// proc face is now external but used to be internal face.
// Check if we have owner or neighbour.
label procOwn = procMesh.faceOwner()[procFaceI];
label masterOwn = masterOwner[masterFaceI];
if (cellProcAddressing[procI][procOwn] == masterOwn)
{
// No turning. Offset by 1.
faceProcAddr[procFaceI]++;
}
else
{
// Turned face.
faceProcAddr[procFaceI] =
-1 - faceProcAddr[procFaceI];
}
// No turning. Offset by 1.
faceProcAddr[procFaceI]++;
faceProcAddr.write();
// From processor cell to reconstructed mesh cell
Info<< "Writing cellProcAddressing to "
<< databases[procI].caseName()
/procMesh.facesInstance()
/polyMesh::meshSubDir
<< endl;
labelIOList
IOobject
(
"cellProcAddressing",
procMesh.facesInstance(),
polyMesh::meshSubDir,
procMesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false // Do not register
),
cellProcAddressing[procI]
).write();
// From processor patch to reconstructed mesh patch
Info<< "Writing boundaryProcAddressing to "
<< databases[procI].caseName()
/procMesh.facesInstance()
/polyMesh::meshSubDir
<< endl;
labelIOList
IOobject
(
"boundaryProcAddressing",
procMesh.facesInstance(),
polyMesh::meshSubDir,
procMesh,
IOobject::NO_READ,
IOobject::NO_WRITE,
false // Do not register
),
boundaryProcAddressing[procI]
).write();
Info<< endl;
}
Info<< "End.\n" << endl;
return 0;
}
// ************************************************************************* //