surfaceMeshTriangulate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  surfaceMeshTriangulate
26 
27 Description
28  Extracts surface from a polyMesh. Depending on output surface format
29  triangulates faces.
30 
31  Region numbers on faces cannot be guaranteed to be the same as the patch
32  indices.
33 
34  Optionally only triangulates named patches.
35 
36  If run in parallel the processor patches get filtered out by default and
37  the mesh gets merged (based on topology).
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "MeshedSurface.H"
42 #include "UnsortedMeshedSurface.H"
43 #include "argList.H"
44 #include "Time.H"
45 #include "polyMesh.H"
46 #include "processorPolyPatch.H"
47 #include "ListListOps.H"
49 #include "globalMeshData.H"
50 #include "globalIndex.H"
51 
52 using namespace Foam;
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 
57 int main(int argc, char *argv[])
58 {
60  (
61  "extract surface from a polyMesh"
62  );
63  argList::validArgs.append("output surface file");
64  #include "addRegionOption.H"
66  (
67  "excludeProcPatches",
68  "exclude processor patches"
69  );
71  (
72  "patches",
73  "(patch0 .. patchN)",
74  "only triangulate selected patches (wildcards supported)"
75  );
77  (
78  "faceZones",
79  "(fz0 .. fzN)",
80  "triangulate selected faceZones (wildcards supported)"
81  );
82 
83  #include "setRootCase.H"
84  #include "createTime.H"
85 
86  const fileName outFileName(args[1]);
87 
88  Info<< "Extracting surface from boundaryMesh ..."
89  << endl << endl;
90 
91  const bool includeProcPatches =
92  !(
93  args.optionFound("excludeProcPatches")
94  || Pstream::parRun()
95  );
96 
97  if (includeProcPatches)
98  {
99  Info<< "Including all processor patches." << nl << endl;
100  }
101  else if (Pstream::parRun())
102  {
103  Info<< "Excluding all processor patches." << nl << endl;
104  }
105 
106  Info<< "Reading mesh from time " << runTime.userTimeName() << endl;
107 
108  #include "createNamedPolyMesh.H"
109 
110 
111  // Create local surface from:
112  // - explicitly named patches only (-patches (at your option)
113  // - all patches (default in sequential mode)
114  // - all non-processor patches (default in parallel mode)
115  // - all non-processor patches (sequential mode, -excludeProcPatches
116  // (at your option)
117 
118  // Construct table of patches to include.
119  const polyBoundaryMesh& bMesh = mesh.boundaryMesh();
120 
122 
123  if (args.optionFound("patches"))
124  {
126  (
127  wordReList(args.optionLookup("patches")())
128  );
129  }
130  else
131  {
133  {
134  const polyPatch& patch = bMesh[patchi];
135 
136  if (includeProcPatches || !isA<processorPolyPatch>(patch))
137  {
139  }
140  }
141  }
142 
143 
144 
145  const faceZoneList& mfz = mesh.faceZones();
146  labelHashSet includeFaceZones(mfz.size());
147 
148  if (args.optionFound("faceZones"))
149  {
150  wordReList zoneNames(args.optionLookup("faceZones")());
151  const wordList allZoneNames(mfz.toc());
152  forAll(zoneNames, i)
153  {
154  const wordRe& zoneName = zoneNames[i];
155 
156  labelList zoneIDs = findStrings(zoneName, allZoneNames);
157 
158  forAll(zoneIDs, j)
159  {
160  includeFaceZones.insert(zoneIDs[j]);
161  }
162 
163  if (zoneIDs.empty())
164  {
166  << "Cannot find any faceZone name matching "
167  << zoneName << endl;
168  }
169 
170  }
171  Info<< "Additionally triangulating faceZones "
172  << UIndirectList<word>(allZoneNames, includeFaceZones.sortedToc())
173  << endl;
174  }
175 
176 
177 
178  // From (name of) patch to compact 'zone' index
179  HashTable<label> compactZoneID(1000);
180  // Mesh face and compact zone indx
181  DynamicList<label> faceLabels;
182  DynamicList<label> compactZones;
183 
184  {
185  // Collect sizes. Hash on names to handle local-only patches (e.g.
186  // processor patches)
187  HashTable<label> patchSize(1000);
188  label nFaces = 0;
190  {
191  const polyPatch& pp = bMesh[iter.key()];
192  patchSize.insert(pp.name(), pp.size());
193  nFaces += pp.size();
194  }
195 
196  HashTable<label> zoneSize(1000);
197  forAllConstIter(labelHashSet, includeFaceZones, iter)
198  {
199  const faceZone& pp = mfz[iter.key()];
200  zoneSize.insert(pp.name(), pp.size());
201  nFaces += pp.size();
202  }
203 
204 
207 
208 
209  // Allocate compact numbering for all patches/faceZones
210  forAllConstIter(HashTable<label>, patchSize, iter)
211  {
212  label sz = compactZoneID.size();
213  compactZoneID.insert(iter.key(), sz);
214  }
215 
216  forAllConstIter(HashTable<label>, zoneSize, iter)
217  {
218  label sz = compactZoneID.size();
219  // Info<< "For faceZone " << iter.key() << " allocating zoneID "
220  // << sz << endl;
221  compactZoneID.insert(iter.key(), sz);
222  }
223 
224 
225  Pstream::mapCombineScatter(compactZoneID);
226 
227 
228  // Rework HashTable into labelList just for speed of conversion
229  labelList patchToCompactZone(bMesh.size(), -1);
230  labelList faceZoneToCompactZone(bMesh.size(), -1);
231  forAllConstIter(HashTable<label>, compactZoneID, iter)
232  {
233  label patchi = bMesh.findIndex(iter.key());
234  if (patchi != -1)
235  {
236  patchToCompactZone[patchi] = iter();
237  }
238  else
239  {
240  label zoneI = mfz.findIndex(iter.key());
241  faceZoneToCompactZone[zoneI] = iter();
242  }
243  }
244 
245 
246  faceLabels.setCapacity(nFaces);
247  compactZones.setCapacity(nFaces);
248 
249  // Collect faces on patches
251  {
252  const polyPatch& pp = bMesh[iter.key()];
253  forAll(pp, i)
254  {
255  faceLabels.append(pp.start()+i);
256  compactZones.append(patchToCompactZone[pp.index()]);
257  }
258  }
259  // Collect faces on faceZones
260  forAllConstIter(labelHashSet, includeFaceZones, iter)
261  {
262  const label fzi = iter.key();
263  const faceZone& fz = mfz[fzi];
264  forAll(fz, i)
265  {
266  faceLabels.append(fz[i]);
267  compactZones.append(faceZoneToCompactZone[fzi]);
268  }
269  }
270  }
271 
272 
273  // Addressing engine for all faces
274  uindirectPrimitivePatch allBoundary
275  (
276  UIndirectList<face>(mesh.faces(), faceLabels),
277  mesh.points()
278  );
279 
280 
281  // Find correspondence to master points
282  labelList pointToGlobal;
283  labelList uniqueMeshPoints;
284  autoPtr<globalIndex> globalNumbers = mesh.globalData().mergePoints
285  (
286  allBoundary.meshPoints(),
287  allBoundary.meshPointMap(),
288  pointToGlobal,
289  uniqueMeshPoints
290  );
291 
292  // Gather all unique points on master
293  List<pointField> gatheredPoints(Pstream::nProcs());
294  gatheredPoints[Pstream::myProcNo()] = pointField
295  (
296  mesh.points(),
297  uniqueMeshPoints
298  );
299  Pstream::gatherList(gatheredPoints);
300 
301  // Gather all faces
302  List<faceList> gatheredFaces(Pstream::nProcs());
303  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
304  forAll(gatheredFaces[Pstream::myProcNo()], i)
305  {
306  inplaceRenumber(pointToGlobal, gatheredFaces[Pstream::myProcNo()][i]);
307  }
308  Pstream::gatherList(gatheredFaces);
309 
310  // Gather all ZoneIDs
311  List<labelList> gatheredZones(Pstream::nProcs());
312  gatheredZones[Pstream::myProcNo()] = move(compactZones);
313  Pstream::gatherList(gatheredZones);
314 
315  // On master combine all points, faces, zones
316  if (Pstream::master())
317  {
318  pointField allPoints = ListListOps::combine<pointField>
319  (
320  gatheredPoints,
322  );
323  gatheredPoints.clear();
324 
325  faceList allFaces = ListListOps::combine<faceList>
326  (
327  gatheredFaces,
329  );
330  gatheredFaces.clear();
331 
332  labelList allZones = ListListOps::combine<labelList>
333  (
334  gatheredZones,
336  );
337  gatheredZones.clear();
338 
339 
340  // Zones
341  surfZoneIdentifierList surfZones(compactZoneID.size());
342  forAllConstIter(HashTable<label>, compactZoneID, iter)
343  {
344  surfZones[iter()] = surfZoneIdentifier(iter.key(), iter());
345  Info<< "surfZone " << iter() << " : " << surfZones[iter()].name()
346  << endl;
347  }
348 
349  UnsortedMeshedSurface<face> unsortedFace
350  (
351  move(allPoints),
352  move(allFaces),
353  move(allZones),
354  move(surfZones)
355  );
356 
357 
358  MeshedSurface<face> sortedFace(unsortedFace);
359 
360  Info<< "Writing merged surface to "
361  << runTime.globalPath()/outFileName << endl;
362 
363  sortedFace.write(runTime.globalPath()/outFileName);
364  }
365 
366  Info<< "End\n" << endl;
367 
368  return 0;
369 }
370 
371 
372 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
wordList toc() const
Return the table of contents.
void setCapacity(const label)
Alter the size of the underlying storage.
Definition: DynamicListI.H:130
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
A list of faces which address into the list of points.
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A surface geometry mesh, in which the surface zone information is conveyed by the 'zoneId' associated...
const word & name() const
Return name.
Definition: Zone.H:184
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void insert(const List< Map< bool >> &zonesIndices)
Insert given indices and corresponding face flips into zones.
Definition: faceZoneList.C:57
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:59
A class for handling file names.
Definition: fileName.H:82
label index() const
Return the index of this patch in the boundaryMesh.
const word & name() const
Return name.
Foam::polyBoundaryMesh.
label findIndex(const word &patchName) const
Find patch index given a name.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch indices corresponding to the given names.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
An identifier for a surface zone on a meshed surface.
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:77
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const polyBoundaryMesh & bMesh
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
static const char nl
Definition: Ostream.H:266
labelHashSet includePatches
Foam::argList args(argc, argv)