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 "indirectPrimitivePatch.H"
55 #include "globalMeshData.H"
56 #include "globalIndex.H"
57 #include "timeSelector.H"
58 
59 using namespace Foam;
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 labelList getSelectedPatches
64 (
66  const wordRes& allow,
67  const wordRes& deny
68 )
69 {
70  // Name-based selection
71  labelList indices
72  (
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 
111 int main(int argc, char *argv[])
112 {
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  );
120 
121  // Less frequently used - reduce some clutter
122  argList::setAdvanced("decomposeParDict");
123 
124  argList::addArgument("output", "The output surface file");
125 
126  #include "addRegionOption.H"
128  (
129  "excludeProcPatches",
130  "Exclude processor patches"
131  );
133  (
134  "faceZones",
135  "wordRes",
136  "Specify single or multiple faceZones to extract\n"
137  "Eg, 'cells' or '( slice \"mfp-.*\" )'"
138  );
140  (
141  "patches",
142  "wordRes",
143  "Specify single patch or multiple patches to extract.\n"
144  "Eg, 'top' or '( front \".*back\" )'"
145  );
147  (
148  "exclude-patches",
149  "wordRes",
150  "Specify single patch or multiple patches to exclude from -patches."
151  " Eg, 'outlet' or '( inlet \".*Wall\" )'",
152  true // mark as an advanced option
153  );
154  argList::addOptionCompat("exclude-patches", {"excludePatches", 2306});
155 
156  #include "setRootCase.H"
157  #include "createTime.H"
158 
159  const auto userOutFileName = args.get<fileName>(1);
160 
161  if (!userOutFileName.has_ext())
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  (!UPstream::parRun() && !args.found("excludeProcPatches"));
172 
173  if (includeProcPatches)
174  {
175  Info<< "Including all processor patches." << nl << endl;
176  }
177  else if (UPstream::parRun())
178  {
179  Info<< "Excluding all processor patches." << nl << endl;
180  }
181 
182  wordRes includePatches, excludePatches;
183  if (args.readListIfPresent<wordRe>("patches", includePatches))
184  {
185  Info<< "Including patches " << flatOutput(includePatches)
186  << nl << endl;
187  }
188  if (args.readListIfPresent<wordRe>("exclude-patches", excludePatches))
189  {
190  Info<< "Excluding patches " << flatOutput(excludePatches)
191  << nl << endl;
192  }
193 
194  // Non-mandatory
195  const wordRes selectedFaceZones(args.getList<wordRe>("faceZones", false));
196  if (selectedFaceZones.size())
197  {
198  Info<< "Including faceZones " << flatOutput(selectedFaceZones)
199  << nl << endl;
200  }
201 
202  Info<< "Reading mesh from time " << runTime.value() << endl;
203 
204  #include "createNamedPolyMesh.H"
205 
206  // User specified times
208 
209  forAll(timeDirs, timeIndex)
210  {
211  runTime.setTime(timeDirs[timeIndex], timeIndex);
212  Info<< "Time [" << timeIndex << "] = " << runTime.timeName();
213 
214  fileName outFileName;
215  if (timeDirs.size() == 1)
216  {
217  outFileName = userOutFileName;
218  }
219  else
220  {
223  {
224  Info<<" ... no mesh change." << nl;
225  continue;
226  }
227 
228  // The filename based on the original, but with additional
229  // time information. The extension was previously checked that
230  // it exists
231  const auto dot = userOutFileName.rfind('.');
232 
233  outFileName =
234  userOutFileName.substr(0, dot) + "_"
235  + Foam::name(runTime.value()) + "."
236  + userOutFileName.ext();
237  }
238 
239  Info<< nl;
240 
241  // Create local surface from:
242  // - explicitly named patches only (-patches (at your option)
243  // - all patches (default in sequential mode)
244  // - all non-processor patches (default in parallel mode)
245  // - all non-processor patches (sequential mode, -excludeProcPatches
246  // (at your option)
247 
248  // Construct table of patches to include.
250 
252  (
253  (includePatches.size() || excludePatches.size())
254  ? getSelectedPatches(bMesh, includePatches, excludePatches)
255  : includeProcPatches
256  ? identity(bMesh.size())
257  : identity(bMesh.nNonProcessor())
258  );
259 
260  labelList faceZoneIds;
261 
262  const faceZoneMesh& fzm = mesh.faceZones();
263 
264  if (selectedFaceZones.size())
265  {
266  faceZoneIds = fzm.indices(selectedFaceZones);
267 
268  Info<< "Additionally extracting faceZones "
269  << fzm.names(selectedFaceZones) << nl;
270  }
271 
272 
273  // From (name of) patch to compact 'zone' index
274  HashTable<label> compactZoneID(1024);
275  // Mesh face and compact zone indx
277  DynamicList<label> compactZones;
278 
279  {
280  // Collect sizes. Hash on names to handle local-only patches (e.g.
281  // processor patches)
282  HashTable<label> patchSize(1024);
283  label nFaces = 0;
284  for (const label patchi : patchIds)
285  {
286  const polyPatch& pp = bMesh[patchi];
287  patchSize.insert(pp.name(), pp.size());
288  nFaces += pp.size();
289  }
290 
291  HashTable<label> zoneSize(1024);
292  for (const label zonei : faceZoneIds)
293  {
294  const faceZone& pp = fzm[zonei];
295  zoneSize.insert(pp.name(), pp.size());
296  nFaces += pp.size();
297  }
298 
299 
302 
303  if (Pstream::master())
304  {
305  // Allocate compact numbering for all patches/faceZones
306  forAllConstIters(patchSize, iter)
307  {
308  compactZoneID.insert(iter.key(), compactZoneID.size());
309  }
310 
311  forAllConstIters(zoneSize, iter)
312  {
313  compactZoneID.insert(iter.key(), compactZoneID.size());
314  }
315  }
316  Pstream::broadcast(compactZoneID);
317 
318 
319  // Rework HashTable into labelList just for speed of conversion
320  labelList patchToCompactZone(bMesh.size(), -1);
321  labelList faceZoneToCompactZone(bMesh.size(), -1);
322  forAllConstIters(compactZoneID, iter)
323  {
324  label patchi = bMesh.findPatchID(iter.key());
325  if (patchi != -1)
326  {
327  patchToCompactZone[patchi] = iter.val();
328  }
329  else
330  {
331  label zoneI = fzm.findZoneID(iter.key());
332  faceZoneToCompactZone[zoneI] = iter.val();
333  }
334  }
335 
336 
337  faceLabels.setCapacity(nFaces);
338  compactZones.setCapacity(nFaces);
339 
340  // Collect faces on patches
341  for (const label patchi : patchIds)
342  {
343  const polyPatch& pp = bMesh[patchi];
344  forAll(pp, i)
345  {
346  faceLabels.append(pp.start()+i);
347  compactZones.append(patchToCompactZone[pp.index()]);
348  }
349  }
350  // Collect faces on faceZones
351  for (const label zonei : faceZoneIds)
352  {
353  const faceZone& pp = fzm[zonei];
354  forAll(pp, i)
355  {
356  faceLabels.append(pp[i]);
357  compactZones.append(faceZoneToCompactZone[pp.index()]);
358  }
359  }
360  }
361 
362 
363  // Addressing engine for all faces
364  uindirectPrimitivePatch allBoundary
365  (
367  mesh.points()
368  );
369 
370 
371  // Find correspondence to master points
372  labelList pointToGlobal;
373  labelList uniqueMeshPoints;
375  (
376  allBoundary.meshPoints(),
377  allBoundary.meshPointMap(),
378  pointToGlobal,
379  uniqueMeshPoints
380  );
381 
382  // Gather all unique points on master
383  List<pointField> gatheredPoints(Pstream::nProcs());
384  gatheredPoints[Pstream::myProcNo()] = pointField
385  (
386  mesh.points(),
387  uniqueMeshPoints
388  );
389  Pstream::gatherList(gatheredPoints);
390 
391  // Gather all faces
392  List<faceList> gatheredFaces(Pstream::nProcs());
393  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
394  for (face& f : gatheredFaces[Pstream::myProcNo()])
395  {
396  inplaceRenumber(pointToGlobal, f);
397  }
398  Pstream::gatherList(gatheredFaces);
399 
400  // Gather all ZoneIDs
401  List<labelList> gatheredZones(Pstream::nProcs());
402  gatheredZones[Pstream::myProcNo()].transfer(compactZones);
403  Pstream::gatherList(gatheredZones);
404 
405  // On master combine all points, faces, zones
406  if (Pstream::master())
407  {
408  pointField allPoints = ListListOps::combine<pointField>
409  (
410  gatheredPoints,
412  );
413  gatheredPoints.clear();
414 
415  faceList allFaces = ListListOps::combine<faceList>
416  (
417  gatheredFaces,
419  );
420  gatheredFaces.clear();
421 
422  labelList allZones = ListListOps::combine<labelList>
423  (
424  gatheredZones,
426  );
427  gatheredZones.clear();
428 
429 
430  // Zones
431  surfZoneIdentifierList surfZones(compactZoneID.size());
432  forAllConstIters(compactZoneID, iter)
433  {
434  surfZones[*iter] = surfZoneIdentifier(iter.key(), *iter);
435  Info<< "surfZone " << *iter
436  << " : " << surfZones[*iter].name()
437  << endl;
438  }
439 
440  UnsortedMeshedSurface<face> unsortedFace
441  (
442  std::move(allPoints),
443  std::move(allFaces),
444  std::move(allZones),
445  surfZones
446  );
447 
448 
449  MeshedSurface<face> sortedFace(unsortedFace);
450 
451  fileName globalCasePath
452  (
453  outFileName.isAbsolute()
454  ? outFileName
455  : (
457  ? runTime.globalPath()/outFileName
458  : runTime.path()/outFileName
459  )
460  );
461 
462  Info<< "Writing merged surface to " << globalCasePath << endl;
463 
464  sortedFace.write(globalCasePath);
465  }
466  }
467 
468  Info<< "End\n" << endl;
469 
470  return 0;
471 }
472 
473 
474 // ************************************************************************* //
A surface geometry mesh, in which the surface zone information is conveyed by the &#39;zoneId&#39; associated...
Definition: MeshedSurface.H:76
static void mapCombineGather(const List< commsStruct > &comms, Container &values, const CombineOp &cop, const int tag, const label comm)
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
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:598
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
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:1004
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
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:1049
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.
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:1074
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:1065
Required Classes.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
static void gatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
Gather data, but keep individual values separate. Uses the specified communication schedule...
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:704
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) zone indices for all matches.
Definition: ZoneMesh.C:435
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:629
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:608
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:935
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:670
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:234
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of 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:56
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
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:362
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:90
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