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-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Application
28  surfaceMeshExtract
29 
30 Group
31  grpSurfaceUtilities
32 
33 Description
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"
47 #include "UnsortedMeshedSurface.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"
54 #include "stringListOps.H" // For stringListOps::findMatching()
55 #include "indirectPrimitivePatch.H"
56 #include "globalMeshData.H"
57 #include "globalIndex.H"
58 #include "timeSelector.H"
59 
60 using namespace Foam;
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 labelList getSelectedPatches
65 (
67  const wordRes& allow,
68  const wordRes& deny
69 )
70 {
71  // Name-based selection
72  labelList indices
73  (
75  (
76  patches,
77  allow,
78  deny,
80  )
81  );
82 
83 
84  // Remove undesirable patches
85 
86  label count = 0;
87  for (const label patchi : indices)
88  {
89  const polyPatch& pp = patches[patchi];
90 
91  if (isType<emptyPolyPatch>(pp))
92  {
93  continue;
94  }
95  else if (Pstream::parRun() && bool(isA<processorPolyPatch>(pp)))
96  {
97  break; // No processor patches for parallel output
98  }
99 
100  indices[count] = patchi;
101  ++count;
102  }
103 
104  indices.resize(count);
105 
106  return indices;
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111 
112 int main(int argc, char *argv[])
113 {
115  (
116  "Extract patch or faceZone surfaces from a polyMesh."
117  " The name is historical, it only triangulates faces"
118  " when the output format requires it."
119  );
121 
122  // Less frequently used - reduce some clutter
123  argList::setAdvanced("decomposeParDict");
124 
125  argList::addArgument("output", "The output surface file");
126 
127  #include "addRegionOption.H"
129  (
130  "excludeProcPatches",
131  "Exclude processor patches"
132  );
134  (
135  "faceZones",
136  "wordRes",
137  "Specify single or multiple faceZones to extract\n"
138  "Eg, 'cells' or '( slice \"mfp-.*\" )'"
139  );
141  (
142  "patches",
143  "wordRes",
144  "Specify single patch or multiple patches to extract.\n"
145  "Eg, 'top' or '( front \".*back\" )'"
146  );
148  (
149  "exclude-patches",
150  "wordRes",
151  "Specify single patch or multiple patches to exclude from -patches."
152  " Eg, 'outlet' or '( inlet \".*Wall\" )'",
153  true // mark as an advanced option
154  );
155  argList::addOptionCompat("exclude-patches", {"excludePatches", 2306});
156 
157  #include "setRootCase.H"
158  #include "createTime.H"
159 
160  const auto userOutFileName = args.get<fileName>(1);
161 
162  if (!userOutFileName.has_ext())
163  {
165  << "Missing extension on output name " << userOutFileName
166  << exit(FatalError);
167  }
168 
169  Info<< "Extracting surface from boundaryMesh ..." << nl << nl;
170 
171  const bool includeProcPatches =
172  (!UPstream::parRun() && !args.found("excludeProcPatches"));
173 
174  if (includeProcPatches)
175  {
176  Info<< "Including all processor patches." << nl << endl;
177  }
178  else if (UPstream::parRun())
179  {
180  Info<< "Excluding all processor patches." << nl << endl;
181  }
182 
183  wordRes includePatches, excludePatches;
184  if (args.readListIfPresent<wordRe>("patches", includePatches))
185  {
186  Info<< "Including patches " << flatOutput(includePatches)
187  << nl << endl;
188  }
189  if (args.readListIfPresent<wordRe>("exclude-patches", excludePatches))
190  {
191  Info<< "Excluding patches " << flatOutput(excludePatches)
192  << nl << endl;
193  }
194 
195  // Non-mandatory
196  const wordRes selectedFaceZones(args.getList<wordRe>("faceZones", false));
197  if (selectedFaceZones.size())
198  {
199  Info<< "Including faceZones " << flatOutput(selectedFaceZones)
200  << nl << endl;
201  }
202 
203  Info<< "Reading mesh from time " << runTime.value() << endl;
204 
205  #include "createNamedPolyMesh.H"
206 
207  // User specified times
209 
210  forAll(timeDirs, timeIndex)
211  {
212  runTime.setTime(timeDirs[timeIndex], timeIndex);
213  Info<< "Time [" << timeIndex << "] = " << runTime.timeName();
214 
215  fileName outFileName;
216  if (timeDirs.size() == 1)
217  {
218  outFileName = userOutFileName;
219  }
220  else
221  {
224  {
225  Info<<" ... no mesh change." << nl;
226  continue;
227  }
228 
229  // The filename based on the original, but with additional
230  // time information. The extension was previously checked that
231  // it exists
232  const auto dot = userOutFileName.rfind('.');
233 
234  outFileName =
235  userOutFileName.substr(0, dot) + "_"
236  + Foam::name(runTime.value()) + "."
237  + userOutFileName.ext();
238  }
239 
240  Info<< nl;
241 
242  // Create local surface from:
243  // - explicitly named patches only (-patches (at your option)
244  // - all patches (default in sequential mode)
245  // - all non-processor patches (default in parallel mode)
246  // - all non-processor patches (sequential mode, -excludeProcPatches
247  // (at your option)
248 
249  // Construct table of patches to include.
251 
253  (
254  (includePatches.size() || excludePatches.size())
255  ? getSelectedPatches(bMesh, includePatches, excludePatches)
256  : includeProcPatches
257  ? identity(bMesh.size())
258  : identity(bMesh.nNonProcessor())
259  );
260 
261  labelList faceZoneIds;
262 
263  const faceZoneMesh& fzm = mesh.faceZones();
264 
265  if (selectedFaceZones.size())
266  {
267  faceZoneIds = fzm.indices(selectedFaceZones);
268 
269  Info<< "Additionally extracting faceZones "
270  << fzm.names(selectedFaceZones) << nl;
271  }
272 
273 
274  // From (name of) patch to compact 'zone' index
275  HashTable<label> compactZoneID(1024);
276  // Mesh face and compact zone indx
278  DynamicList<label> compactZones;
279 
280  {
281  // Collect sizes. Hash on names to handle local-only patches (e.g.
282  // processor patches)
283  HashTable<label> patchSize(1024);
284  label nFaces = 0;
285  for (const label patchi : patchIds)
286  {
287  const polyPatch& pp = bMesh[patchi];
288  patchSize.insert(pp.name(), pp.size());
289  nFaces += pp.size();
290  }
291 
292  HashTable<label> zoneSize(1024);
293  for (const label zonei : faceZoneIds)
294  {
295  const faceZone& pp = fzm[zonei];
296  zoneSize.insert(pp.name(), pp.size());
297  nFaces += pp.size();
298  }
299 
300 
303 
304  if (Pstream::master())
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  (
368  mesh.points()
369  );
370 
371 
372  // Find correspondence to master points
373  labelList pointToGlobal;
374  labelList uniqueMeshPoints;
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  for (face& f : gatheredFaces[Pstream::myProcNo()])
396  {
397  inplaceRenumber(pointToGlobal, f);
398  }
399  Pstream::gatherList(gatheredFaces);
400 
401  // Gather all ZoneIDs
402  List<labelList> gatheredZones(Pstream::nProcs());
403  gatheredZones[Pstream::myProcNo()].transfer(compactZones);
404  Pstream::gatherList(gatheredZones);
405 
406  // On master combine all points, faces, zones
407  if (Pstream::master())
408  {
409  pointField allPoints = ListListOps::combine<pointField>
410  (
411  gatheredPoints,
413  );
414  gatheredPoints.clear();
415 
416  faceList allFaces = ListListOps::combine<faceList>
417  (
418  gatheredFaces,
420  );
421  gatheredFaces.clear();
422 
423  labelList allZones = ListListOps::combine<labelList>
424  (
425  gatheredZones,
427  );
428  gatheredZones.clear();
429 
430 
431  // Zones
432  surfZoneIdentifierList surfZones(compactZoneID.size());
433  forAllConstIters(compactZoneID, iter)
434  {
435  surfZones[*iter] = surfZoneIdentifier(iter.key(), *iter);
436  Info<< "surfZone " << *iter
437  << " : " << surfZones[*iter].name()
438  << endl;
439  }
440 
441  UnsortedMeshedSurface<face> unsortedFace
442  (
443  std::move(allPoints),
444  std::move(allFaces),
445  std::move(allZones),
446  surfZones
447  );
448 
449 
450  MeshedSurface<face> sortedFace(unsortedFace);
451 
452  fileName globalCasePath
453  (
454  outFileName.isAbsolute()
455  ? outFileName
456  : (
458  ? runTime.globalPath()/outFileName
459  : runTime.path()/outFileName
460  )
461  );
462 
463  Info<< "Writing merged surface to " << globalCasePath << endl;
464 
465  sortedFace.write(globalCasePath);
466  }
467  }
468 
469  Info<< "End\n" << endl;
470 
471  return 0;
472 }
473 
474 
475 // ************************************************************************* //
A surface geometry mesh, in which the surface zone information is conveyed by the &#39;zoneId&#39; associated...
Definition: MeshedSurface.H:76
labelList patchIds
const Type & value() const noexcept
Return const reference to value.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static void mapCombineGather(Container &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combine Map elements.
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
A class for handling file names.
Definition: fileName.H:72
static void setAdvanced(const word &optName, bool advanced=true)
Set an existing option as being &#39;advanced&#39; or normal.
Definition: argList.C:404
Identifies a surface patch/zone by name and index, with optional geometric type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:521
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Object access operator or list access operator (default is pass-through)
Definition: UList.H:1013
Database for mesh data, solution data, solver performance and other reduced data. ...
Definition: meshState.H:53
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void gatherList(const UList< commsStruct > &comms, UList< T > &values, const int tag, const label comm)
Gather data, but keep individual values separate. Uses the specified communication schedule...
List< T > getList(const label index) const
Get a List of values from the argument at index.
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
static void addOptionCompat(const word &optName, std::pair< const char *, int > compat)
Specify an alias for the option name.
Definition: argList.C:418
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:303
Required Classes.
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
Operations on lists of strings.
static bool isAbsolute(const std::string &str)
Return true if filename starts with a &#39;/&#39; or &#39;\&#39; or (windows-only) with a filesystem-root.
Definition: fileNameI.H:129
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
labelList faceLabels(nFaceLabels)
bool processorCase() const noexcept
True if this is a processor case.
Definition: TimePathsI.H:52
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:340
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word ext() const
Return file name extension (part after last .)
Definition: wordI.H:171
A list of faces which address into the list of points.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1077
Required Classes.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:715
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) zone indices for all matches.
Definition: ZoneMesh.C:524
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:718
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
fileName globalPath() const
Return global path for the case = rootPath/globalCaseName. Same as TimePaths::globalPath() ...
Definition: Time.H:512
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
labelList findMatching(const StringListType &input, const wordRes::filter &pred, AccessOp aop=identityOp())
Return ids for items with matching names.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:385
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings...
Definition: wordRe.H:78
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:671
bool readListIfPresent(const word &optName, List< T > &list) const
If named option is present, get a List of values treating a single entry like a list of size 1...
Definition: argListI.H:387
labelList f(nPoints)
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:260
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
const polyBoundaryMesh & patches
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:351
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:451
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
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:69
Foam::argList args(argc, argv)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
label timeIndex
Definition: getTimeIndex.H:24
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225