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-2022 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  "excludePatches",
149  "wordRes",
150  "Specify single patch or multiple patches to exclude from writing."
151  " Eg, 'outlet' or '( inlet \".*Wall\" )'",
152  true // mark as an advanced option
153  );
154 
155  #include "setRootCase.H"
156  #include "createTime.H"
157 
158  const auto userOutFileName = args.get<fileName>(1);
159 
160  if (!userOutFileName.has_ext())
161  {
163  << "Missing extension on output name " << userOutFileName
164  << exit(FatalError);
165  }
166 
167  Info<< "Extracting surface from boundaryMesh ..." << nl << nl;
168 
169  const bool includeProcPatches =
170  !(
171  args.found("excludeProcPatches")
172  || Pstream::parRun()
173  );
174 
175  if (includeProcPatches)
176  {
177  Info<< "Including all processor patches." << nl << endl;
178  }
179  else if (Pstream::parRun())
180  {
181  Info<< "Excluding all processor patches." << nl << endl;
182  }
183 
184  wordRes includePatches, excludePatches;
185  if (args.readListIfPresent<wordRe>("patches", includePatches))
186  {
187  Info<< "Including patches " << flatOutput(includePatches)
188  << nl << endl;
189  }
190  if (args.readListIfPresent<wordRe>("excludePatches", excludePatches))
191  {
192  Info<< "Excluding patches " << flatOutput(excludePatches)
193  << nl << endl;
194  }
195 
196  // Non-mandatory
197  const wordRes selectedFaceZones(args.getList<wordRe>("faceZones", false));
198  if (selectedFaceZones.size())
199  {
200  Info<< "Including faceZones " << flatOutput(selectedFaceZones)
201  << nl << endl;
202  }
203 
204  Info<< "Reading mesh from time " << runTime.value() << endl;
205 
206  #include "createNamedPolyMesh.H"
207 
208  // User specified times
210 
211  forAll(timeDirs, timeIndex)
212  {
213  runTime.setTime(timeDirs[timeIndex], timeIndex);
214  Info<< "Time [" << timeIndex << "] = " << runTime.timeName();
215 
216  fileName outFileName;
217  if (timeDirs.size() == 1)
218  {
219  outFileName = userOutFileName;
220  }
221  else
222  {
224  if (timeIndex && meshState == polyMesh::UNCHANGED)
225  {
226  Info<<" ... no mesh change." << nl;
227  continue;
228  }
229 
230  // The filename based on the original, but with additional
231  // time information. The extension was previously checked that
232  // it exists
233  const auto dot = userOutFileName.rfind('.');
234 
235  outFileName =
236  userOutFileName.substr(0, dot) + "_"
237  + Foam::name(runTime.value()) + "."
238  + userOutFileName.ext();
239  }
240 
241  Info<< nl;
242 
243  // Create local surface from:
244  // - explicitly named patches only (-patches (at your option)
245  // - all patches (default in sequential mode)
246  // - all non-processor patches (default in parallel mode)
247  // - all non-processor patches (sequential mode, -excludeProcPatches
248  // (at your option)
249 
250  // Construct table of patches to include.
252 
254  (
255  (includePatches.size() || excludePatches.size())
256  ? getSelectedPatches(bMesh, includePatches, excludePatches)
257  : includeProcPatches
258  ? identity(bMesh.size())
259  : identity(bMesh.nNonProcessor())
260  );
261 
262  labelList faceZoneIds;
263 
264  const faceZoneMesh& fzm = mesh.faceZones();
265 
266  if (selectedFaceZones.size())
267  {
268  faceZoneIds = fzm.indices(selectedFaceZones);
269 
270  Info<< "Additionally extracting faceZones "
271  << fzm.names(selectedFaceZones) << nl;
272  }
273 
274 
275  // From (name of) patch to compact 'zone' index
276  HashTable<label> compactZoneID(1024);
277  // Mesh face and compact zone indx
278  DynamicList<label> faceLabels;
279  DynamicList<label> compactZones;
280 
281  {
282  // Collect sizes. Hash on names to handle local-only patches (e.g.
283  // processor patches)
284  HashTable<label> patchSize(1024);
285  label nFaces = 0;
286  for (const label patchi : patchIds)
287  {
288  const polyPatch& pp = bMesh[patchi];
289  patchSize.insert(pp.name(), pp.size());
290  nFaces += pp.size();
291  }
292 
293  HashTable<label> zoneSize(1024);
294  for (const label zonei : faceZoneIds)
295  {
296  const faceZone& pp = fzm[zonei];
297  zoneSize.insert(pp.name(), pp.size());
298  nFaces += pp.size();
299  }
300 
301 
304 
305  if (Pstream::master())
306  {
307  // Allocate compact numbering for all patches/faceZones
308  forAllConstIters(patchSize, iter)
309  {
310  compactZoneID.insert(iter.key(), compactZoneID.size());
311  }
312 
313  forAllConstIters(zoneSize, iter)
314  {
315  compactZoneID.insert(iter.key(), compactZoneID.size());
316  }
317  }
318  Pstream::broadcast(compactZoneID);
319 
320 
321  // Rework HashTable into labelList just for speed of conversion
322  labelList patchToCompactZone(bMesh.size(), -1);
323  labelList faceZoneToCompactZone(bMesh.size(), -1);
324  forAllConstIters(compactZoneID, iter)
325  {
326  label patchi = bMesh.findPatchID(iter.key());
327  if (patchi != -1)
328  {
329  patchToCompactZone[patchi] = iter.val();
330  }
331  else
332  {
333  label zoneI = fzm.findZoneID(iter.key());
334  faceZoneToCompactZone[zoneI] = iter.val();
335  }
336  }
337 
338 
339  faceLabels.setCapacity(nFaces);
340  compactZones.setCapacity(nFaces);
341 
342  // Collect faces on patches
343  for (const label patchi : patchIds)
344  {
345  const polyPatch& pp = bMesh[patchi];
346  forAll(pp, i)
347  {
348  faceLabels.append(pp.start()+i);
349  compactZones.append(patchToCompactZone[pp.index()]);
350  }
351  }
352  // Collect faces on faceZones
353  for (const label zonei : faceZoneIds)
354  {
355  const faceZone& pp = fzm[zonei];
356  forAll(pp, i)
357  {
358  faceLabels.append(pp[i]);
359  compactZones.append(faceZoneToCompactZone[pp.index()]);
360  }
361  }
362  }
363 
364 
365  // Addressing engine for all faces
366  uindirectPrimitivePatch allBoundary
367  (
368  UIndirectList<face>(mesh.faces(), faceLabels),
369  mesh.points()
370  );
371 
372 
373  // Find correspondence to master points
374  labelList pointToGlobal;
375  labelList uniqueMeshPoints;
377  (
378  allBoundary.meshPoints(),
379  allBoundary.meshPointMap(),
380  pointToGlobal,
381  uniqueMeshPoints
382  );
383 
384  // Gather all unique points on master
385  List<pointField> gatheredPoints(Pstream::nProcs());
386  gatheredPoints[Pstream::myProcNo()] = pointField
387  (
388  mesh.points(),
389  uniqueMeshPoints
390  );
391  Pstream::gatherList(gatheredPoints);
392 
393  // Gather all faces
394  List<faceList> gatheredFaces(Pstream::nProcs());
395  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
396  for (face& f : gatheredFaces[Pstream::myProcNo()])
397  {
398  inplaceRenumber(pointToGlobal, f);
399  }
400  Pstream::gatherList(gatheredFaces);
401 
402  // Gather all ZoneIDs
403  List<labelList> gatheredZones(Pstream::nProcs());
404  gatheredZones[Pstream::myProcNo()].transfer(compactZones);
405  Pstream::gatherList(gatheredZones);
406 
407  // On master combine all points, faces, zones
408  if (Pstream::master())
409  {
410  pointField allPoints = ListListOps::combine<pointField>
411  (
412  gatheredPoints,
414  );
415  gatheredPoints.clear();
416 
417  faceList allFaces = ListListOps::combine<faceList>
418  (
419  gatheredFaces,
421  );
422  gatheredFaces.clear();
423 
424  labelList allZones = ListListOps::combine<labelList>
425  (
426  gatheredZones,
428  );
429  gatheredZones.clear();
430 
431 
432  // Zones
433  surfZoneIdentifierList surfZones(compactZoneID.size());
434  forAllConstIters(compactZoneID, iter)
435  {
436  surfZones[*iter] = surfZoneIdentifier(iter.key(), *iter);
437  Info<< "surfZone " << *iter
438  << " : " << surfZones[*iter].name()
439  << endl;
440  }
441 
442  UnsortedMeshedSurface<face> unsortedFace
443  (
444  std::move(allPoints),
445  std::move(allFaces),
446  std::move(allZones),
447  surfZones
448  );
449 
450 
451  MeshedSurface<face> sortedFace(unsortedFace);
452 
453  fileName globalCasePath
454  (
455  outFileName.isAbsolute()
456  ? outFileName
457  : (
459  ? runTime.globalPath()/outFileName
460  : runTime.path()/outFileName
461  )
462  );
463 
464  Info<< "Writing merged surface to " << globalCasePath << endl;
465 
466  sortedFace.write(globalCasePath);
467  }
468  }
469 
470  Info<< "End\n" << endl;
471 
472  return 0;
473 }
474 
475 
476 // ************************************************************************* //
A surface geometry mesh, in which the surface zone information is conveyed by the &#39;zoneId&#39; associated...
Definition: MeshedSurface.H:79
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
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:453
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
fileName path() const
Return path.
Definition: Time.H:449
A class for handling file names.
Definition: fileName.H:71
static void setAdvanced(const word &optName, bool advanced=true)
Set an existing option as being &#39;advanced&#39; or normal.
Definition: argList.C:395
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:578
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Object access operator or list access operator.
Definition: UList.H:945
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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:639
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:365
void setCapacity(const label len)
Alter the size of the underlying storage.
Definition: DynamicListI.H:296
Required Variables.
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)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:688
bool processorCase() const noexcept
Return true if this is a processor case.
Definition: TimePathsI.H:29
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 processes in communicator.
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:1066
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
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) is 1 for serial run.
Definition: UPstream.H:656
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 INVALID.
Definition: exprTraits.C:52
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:668
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) zone indices for all matches.
Definition: ZoneMesh.C:365
label findZoneID(const word &zoneName) const
Find zone index by name, return -1 if not found.
Definition: ZoneMesh.C:513
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:31
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:47
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:376
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1293
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:977
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:558
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1091
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...
const word & name() const noexcept
The patch name.
label index() const noexcept
The index of this zone in the zone list.
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:646
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:57
static bool master(const label communicator=worldComm)
Am I the master rank.
Definition: UPstream.H:672
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:441
const polyBoundaryMesh & patches
const word & name() const noexcept
The zone name.
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:342
wordList names() const
A list of the zone names.
Definition: ZoneMesh.C:292
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
label index() const noexcept
The index of this patch in the boundaryMesh.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:89
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
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
label timeIndex
Definition: getTimeIndex.H:24
fileName globalPath() const
Return global path for the case.
Definition: TimePathsI.H:73
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225