surfaceMeshExtract.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-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2021 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
27Application
28 surfaceMeshExtract
29
30Group
31 grpSurfaceUtilities
32
33Description
34 Extract patch or faceZone surfaces from a polyMesh.
35 Depending on output surface format triangulates faces.
36
37 Region numbers on faces no guaranteed to be the same as the patch indices.
38
39 Optionally only extracts named patches.
40
41 If run in parallel, processor patches get filtered out by default and
42 the mesh is merged (based on topology).
43
44\*---------------------------------------------------------------------------*/
45
46#include "MeshedSurface.H"
48#include "argList.H"
49#include "Time.H"
50#include "polyMesh.H"
51#include "emptyPolyPatch.H"
52#include "processorPolyPatch.H"
53#include "ListListOps.H"
55#include "globalMeshData.H"
56#include "globalIndex.H"
57#include "timeSelector.H"
58
59using namespace Foam;
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63labelList getSelectedPatches
64(
66 const wordRes& allow,
67 const wordRes& deny
68)
69{
70 // Name-based selection
71 labelList indices
72 (
73 stringListOps::findMatching
74 (
75 patches,
76 allow,
77 deny,
79 )
80 );
81
82
83 // Remove undesirable patches
84
85 label count = 0;
86 for (const label patchi : indices)
87 {
88 const polyPatch& pp = patches[patchi];
89
90 if (isType<emptyPolyPatch>(pp))
91 {
92 continue;
93 }
94 else if (Pstream::parRun() && bool(isA<processorPolyPatch>(pp)))
95 {
96 break; // No processor patches for parallel output
97 }
98
99 indices[count] = patchi;
100 ++count;
101 }
102
103 indices.resize(count);
104
105 return indices;
106}
107
108
109// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
110
111int main(int argc, char *argv[])
112{
113 argList::addNote
114 (
115 "Extract patch or faceZone surfaces from a polyMesh."
116 " The name is historical, it only triangulates faces"
117 " when the output format requires it."
118 );
119 timeSelector::addOptions();
120
121 // Less frequently used - reduce some clutter
122 argList::setAdvanced("decomposeParDict");
123 argList::setAdvanced("noFunctionObjects");
124
125 argList::addArgument("output", "The output surface file");
126
127 #include "addRegionOption.H"
128 argList::addBoolOption
129 (
130 "excludeProcPatches",
131 "Exclude processor patches"
132 );
133 argList::addOption
134 (
135 "faceZones",
136 "wordRes",
137 "Specify single or multiple faceZones to extract\n"
138 "Eg, 'cells' or '( slice \"mfp-.*\" )'"
139 );
140 argList::addOption
141 (
142 "patches",
143 "wordRes",
144 "Specify single patch or multiple patches to extract.\n"
145 "Eg, 'top' or '( front \".*back\" )'"
146 );
147 argList::addOption
148 (
149 "excludePatches",
150 "wordRes",
151 "Specify single patch or multiple patches to exclude from writing."
152 " Eg, 'outlet' or '( inlet \".*Wall\" )'",
153 true // mark as an advanced option
154 );
155
156 #include "setRootCase.H"
157 #include "createTime.H"
158
159 const auto userOutFileName = args.get<fileName>(1);
160
161 if (!userOutFileName.hasExt())
162 {
164 << "Missing extension on output name " << userOutFileName
165 << exit(FatalError);
166 }
167
168 Info<< "Extracting surface from boundaryMesh ..." << nl << nl;
169
170 const bool includeProcPatches =
171 !(
172 args.found("excludeProcPatches")
173 || Pstream::parRun()
174 );
175
176 if (includeProcPatches)
177 {
178 Info<< "Including all processor patches." << nl << endl;
179 }
180 else if (Pstream::parRun())
181 {
182 Info<< "Excluding all processor patches." << nl << endl;
183 }
184
185 wordRes includePatches, excludePatches;
186 if (args.readListIfPresent<wordRe>("patches", includePatches))
187 {
188 Info<< "Including patches " << flatOutput(includePatches)
189 << nl << endl;
190 }
191 if (args.readListIfPresent<wordRe>("excludePatches", excludePatches))
192 {
193 Info<< "Excluding patches " << flatOutput(excludePatches)
194 << nl << endl;
195 }
196
197 // Non-mandatory
198 const wordRes selectedFaceZones(args.getList<wordRe>("faceZones", false));
199 if (selectedFaceZones.size())
200 {
201 Info<< "Including faceZones " << flatOutput(selectedFaceZones)
202 << nl << endl;
203 }
204
205 Info<< "Reading mesh from time " << runTime.value() << endl;
206
207 #include "createNamedPolyMesh.H"
208
209 // User specified times
210 instantList timeDirs = timeSelector::select0(runTime, args);
211
212 forAll(timeDirs, timeIndex)
213 {
214 runTime.setTime(timeDirs[timeIndex], timeIndex);
215 Info<< "Time [" << timeIndex << "] = " << runTime.timeName();
216
217 fileName outFileName;
218 if (timeDirs.size() == 1)
219 {
220 outFileName = userOutFileName;
221 }
222 else
223 {
224 polyMesh::readUpdateState meshState = mesh.readUpdate();
225 if (timeIndex && meshState == polyMesh::UNCHANGED)
226 {
227 Info<<" ... no mesh change." << nl;
228 continue;
229 }
230
231 // The filename based on the original, but with additional
232 // time information. The extension was previously checked that
233 // it exists
234 const auto dot = userOutFileName.rfind('.');
235
236 outFileName =
237 userOutFileName.substr(0, dot) + "_"
238 + Foam::name(runTime.value()) + "."
239 + userOutFileName.ext();
240 }
241
242 Info<< nl;
243
244 // Create local surface from:
245 // - explicitly named patches only (-patches (at your option)
246 // - all patches (default in sequential mode)
247 // - all non-processor patches (default in parallel mode)
248 // - all non-processor patches (sequential mode, -excludeProcPatches
249 // (at your option)
250
251 // Construct table of patches to include.
252 const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
253
255 (
256 (includePatches.size() || excludePatches.size())
257 ? getSelectedPatches(bMesh, includePatches, excludePatches)
258 : includeProcPatches
259 ? identity(bMesh.size())
261 );
262
263 labelList faceZoneIds;
264
265 const faceZoneMesh& fzm = mesh.faceZones();
266
267 if (selectedFaceZones.size())
268 {
269 faceZoneIds = fzm.indices(selectedFaceZones);
270
271 Info<< "Additionally extracting faceZones "
272 << fzm.names(selectedFaceZones) << nl;
273 }
274
275
276 // From (name of) patch to compact 'zone' index
277 HashTable<label> compactZoneID(1024);
278 // Mesh face and compact zone indx
279 DynamicList<label> faceLabels;
280 DynamicList<label> compactZones;
281
282 {
283 // Collect sizes. Hash on names to handle local-only patches (e.g.
284 // processor patches)
285 HashTable<label> patchSize(1024);
286 label nFaces = 0;
287 for (const label patchi : patchIds)
288 {
289 const polyPatch& pp = bMesh[patchi];
290 patchSize.insert(pp.name(), pp.size());
291 nFaces += pp.size();
292 }
293
294 HashTable<label> zoneSize(1024);
295 for (const label zonei : faceZoneIds)
296 {
297 const faceZone& pp = fzm[zonei];
298 zoneSize.insert(pp.name(), pp.size());
299 nFaces += pp.size();
300 }
301
302
303 Pstream::mapCombineGather(patchSize, plusEqOp<label>());
304 Pstream::mapCombineGather(zoneSize, plusEqOp<label>());
305
306 // Allocate compact numbering for all patches/faceZones
307 forAllConstIters(patchSize, iter)
308 {
309 compactZoneID.insert(iter.key(), compactZoneID.size());
310 }
311
312 forAllConstIters(zoneSize, iter)
313 {
314 compactZoneID.insert(iter.key(), compactZoneID.size());
315 }
316
317 Pstream::broadcast(compactZoneID);
318
319
320 // Rework HashTable into labelList just for speed of conversion
321 labelList patchToCompactZone(bMesh.size(), -1);
322 labelList faceZoneToCompactZone(bMesh.size(), -1);
323 forAllConstIters(compactZoneID, iter)
324 {
325 label patchi = bMesh.findPatchID(iter.key());
326 if (patchi != -1)
327 {
328 patchToCompactZone[patchi] = iter.val();
329 }
330 else
331 {
332 label zoneI = fzm.findZoneID(iter.key());
333 faceZoneToCompactZone[zoneI] = iter.val();
334 }
335 }
336
337
338 faceLabels.setCapacity(nFaces);
339 compactZones.setCapacity(nFaces);
340
341 // Collect faces on patches
342 for (const label patchi : patchIds)
343 {
344 const polyPatch& pp = bMesh[patchi];
345 forAll(pp, i)
346 {
347 faceLabels.append(pp.start()+i);
348 compactZones.append(patchToCompactZone[pp.index()]);
349 }
350 }
351 // Collect faces on faceZones
352 for (const label zonei : faceZoneIds)
353 {
354 const faceZone& pp = fzm[zonei];
355 forAll(pp, i)
356 {
357 faceLabels.append(pp[i]);
358 compactZones.append(faceZoneToCompactZone[pp.index()]);
359 }
360 }
361 }
362
363
364 // Addressing engine for all faces
365 uindirectPrimitivePatch allBoundary
366 (
367 UIndirectList<face>(mesh.faces(), faceLabels),
368 mesh.points()
369 );
370
371
372 // Find correspondence to master points
373 labelList pointToGlobal;
374 labelList uniqueMeshPoints;
375 autoPtr<globalIndex> globalNumbers = mesh.globalData().mergePoints
376 (
377 allBoundary.meshPoints(),
378 allBoundary.meshPointMap(),
379 pointToGlobal,
380 uniqueMeshPoints
381 );
382
383 // Gather all unique points on master
384 List<pointField> gatheredPoints(Pstream::nProcs());
385 gatheredPoints[Pstream::myProcNo()] = pointField
386 (
387 mesh.points(),
388 uniqueMeshPoints
389 );
390 Pstream::gatherList(gatheredPoints);
391
392 // Gather all faces
393 List<faceList> gatheredFaces(Pstream::nProcs());
394 gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
395 forAll(gatheredFaces[Pstream::myProcNo()], i)
396 {
398 (
399 pointToGlobal,
400 gatheredFaces[Pstream::myProcNo()][i]
401 );
402 }
403 Pstream::gatherList(gatheredFaces);
404
405 // Gather all ZoneIDs
406 List<labelList> gatheredZones(Pstream::nProcs());
407 gatheredZones[Pstream::myProcNo()].transfer(compactZones);
408 Pstream::gatherList(gatheredZones);
409
410 // On master combine all points, faces, zones
411 if (Pstream::master())
412 {
413 pointField allPoints = ListListOps::combine<pointField>
414 (
415 gatheredPoints,
417 );
418 gatheredPoints.clear();
419
420 faceList allFaces = ListListOps::combine<faceList>
421 (
422 gatheredFaces,
424 );
425 gatheredFaces.clear();
426
427 labelList allZones = ListListOps::combine<labelList>
428 (
429 gatheredZones,
431 );
432 gatheredZones.clear();
433
434
435 // Zones
436 surfZoneIdentifierList surfZones(compactZoneID.size());
437 forAllConstIters(compactZoneID, iter)
438 {
439 surfZones[*iter] = surfZoneIdentifier(iter.key(), *iter);
440 Info<< "surfZone " << *iter
441 << " : " << surfZones[*iter].name()
442 << endl;
443 }
444
445 UnsortedMeshedSurface<face> unsortedFace
446 (
447 std::move(allPoints),
448 std::move(allFaces),
449 std::move(allZones),
450 surfZones
451 );
452
453
454 MeshedSurface<face> sortedFace(unsortedFace);
455
456 fileName globalCasePath
457 (
458 outFileName.isAbsolute()
459 ? outFileName
460 : (
461 runTime.processorCase()
462 ? runTime.globalPath()/outFileName
463 : runTime.path()/outFileName
464 )
465 );
466
467 Info<< "Writing merged surface to " << globalCasePath << endl;
468
469 sortedFace.write(globalCasePath);
470 }
471 }
472
473 Info<< "End\n" << endl;
474
475 return 0;
476}
477
478
479// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:303
A HashTable similar to std::unordered_map.
Definition: HashTable.H:123
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Definition: MeshedSurface.H:99
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A list of faces which address into the list of points.
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A surface geometry mesh, in which the surface zone information is conveyed by the 'zoneId' associated...
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:525
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:304
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) zone indices for all matches.
Definition: ZoneMesh.C:377
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:278
bool readListIfPresent(const word &optName, List< T > &list) const
Definition: argListI.H:394
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178
List< T > getList(const label index) const
Get a List of values from the argument at index.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:67
A class for handling file names.
Definition: fileName.H:76
static bool isAbsolute(const std::string &str)
Definition: fileNameI.H:136
label index() const noexcept
The index of this patch in the boundaryMesh.
const word & name() const noexcept
The patch name.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
Identifies a surface patch/zone by name and index, with optional geometric type.
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
word ext() const
Return file name extension (part after last .)
Definition: word.C:126
label index() const noexcept
The index of this zone in the zone list.
const word & name() const noexcept
The zone name.
labelList patchIds
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const label nNonProcessor
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:78
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
void inplaceRenumber(const labelUList &oldToNew, IntListType &lists)
Inplace renumber the values (not the indices) of a list of lists.
Definition: ensightOutput.H:90
Namespace for OpenFOAM.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
messageStream Info
Information stream (stdout output on master, null elsewhere)
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:238