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 "addMeshOption.H"
65  #include "addRegionOption.H"
67  (
68  "excludeProcPatches",
69  "exclude processor patches"
70  );
72  (
73  "patches",
74  "(patch0 .. patchN)",
75  "only triangulate selected patches (wildcards supported)"
76  );
78  (
79  "faceZones",
80  "(fz0 .. fzN)",
81  "triangulate selected faceZones (wildcards supported)"
82  );
83 
84  #include "setRootCase.H"
85  #include "createTime.H"
86 
87  const fileName outFileName(args[1]);
88 
89  Info<< "Extracting surface from boundaryMesh ..."
90  << endl << endl;
91 
92  const bool includeProcPatches =
93  !(
94  args.optionFound("excludeProcPatches")
95  || Pstream::parRun()
96  );
97 
98  if (includeProcPatches)
99  {
100  Info<< "Including all processor patches." << nl << endl;
101  }
102  else if (Pstream::parRun())
103  {
104  Info<< "Excluding all processor patches." << nl << endl;
105  }
106 
107  Info<< "Reading mesh from time " << runTime.userTimeName() << endl;
108 
109  #include "createSpecifiedPolyMesh.H"
110 
111 
112  // Create local surface from:
113  // - explicitly named patches only (-patches (at your option)
114  // - all patches (default in sequential mode)
115  // - all non-processor patches (default in parallel mode)
116  // - all non-processor patches (sequential mode, -excludeProcPatches
117  // (at your option)
118 
119  // Construct table of patches to include.
121 
123 
124  if (args.optionFound("patches"))
125  {
127  (
128  wordReList(args.optionLookup("patches")())
129  );
130  }
131  else
132  {
134  {
135  const polyPatch& patch = bMesh[patchi];
136 
137  if (includeProcPatches || !isA<processorPolyPatch>(patch))
138  {
140  }
141  }
142  }
143 
144 
145 
146  const faceZoneList& mfz = mesh.faceZones();
147  labelHashSet includeFaceZones(mfz.size());
148 
149  if (args.optionFound("faceZones"))
150  {
151  wordReList zoneNames(args.optionLookup("faceZones")());
152  const wordList allZoneNames(mfz.toc());
153  forAll(zoneNames, i)
154  {
155  const wordRe& zoneName = zoneNames[i];
156 
157  labelList zoneIDs = findStrings(zoneName, allZoneNames);
158 
159  forAll(zoneIDs, j)
160  {
161  includeFaceZones.insert(zoneIDs[j]);
162  }
163 
164  if (zoneIDs.empty())
165  {
167  << "Cannot find any faceZone name matching "
168  << zoneName << endl;
169  }
170 
171  }
172  Info<< "Additionally triangulating faceZones "
173  << UIndirectList<word>(allZoneNames, includeFaceZones.sortedToc())
174  << endl;
175  }
176 
177 
178 
179  // From (name of) patch to compact 'zone' index
180  HashTable<label> compactZoneID(1000);
181  // Mesh face and compact zone indx
182  DynamicList<label> faceLabels;
183  DynamicList<label> compactZones;
184 
185  {
186  // Collect sizes. Hash on names to handle local-only patches (e.g.
187  // processor patches)
188  HashTable<label> patchSize(1000);
189  label nFaces = 0;
191  {
192  const polyPatch& pp = bMesh[iter.key()];
193  patchSize.insert(pp.name(), pp.size());
194  nFaces += pp.size();
195  }
196 
197  HashTable<label> zoneSize(1000);
198  forAllConstIter(labelHashSet, includeFaceZones, iter)
199  {
200  const faceZone& pp = mfz[iter.key()];
201  zoneSize.insert(pp.name(), pp.size());
202  nFaces += pp.size();
203  }
204 
205 
208 
209 
210  // Allocate compact numbering for all patches/faceZones
211  forAllConstIter(HashTable<label>, patchSize, iter)
212  {
213  label sz = compactZoneID.size();
214  compactZoneID.insert(iter.key(), sz);
215  }
216 
217  forAllConstIter(HashTable<label>, zoneSize, iter)
218  {
219  label sz = compactZoneID.size();
220  // Info<< "For faceZone " << iter.key() << " allocating zoneID "
221  // << sz << endl;
222  compactZoneID.insert(iter.key(), sz);
223  }
224 
225 
226  Pstream::mapCombineScatter(compactZoneID);
227 
228 
229  // Rework HashTable into labelList just for speed of conversion
230  labelList patchToCompactZone(bMesh.size(), -1);
231  labelList faceZoneToCompactZone(bMesh.size(), -1);
232  forAllConstIter(HashTable<label>, compactZoneID, iter)
233  {
234  label patchi = bMesh.findIndex(iter.key());
235  if (patchi != -1)
236  {
237  patchToCompactZone[patchi] = iter();
238  }
239  else
240  {
241  label zoneI = mfz.findIndex(iter.key());
242  faceZoneToCompactZone[zoneI] = iter();
243  }
244  }
245 
246 
247  faceLabels.setCapacity(nFaces);
248  compactZones.setCapacity(nFaces);
249 
250  // Collect faces on patches
252  {
253  const polyPatch& pp = bMesh[iter.key()];
254  forAll(pp, i)
255  {
256  faceLabels.append(pp.start()+i);
257  compactZones.append(patchToCompactZone[pp.index()]);
258  }
259  }
260  // Collect faces on faceZones
261  forAllConstIter(labelHashSet, includeFaceZones, iter)
262  {
263  const label fzi = iter.key();
264  const faceZone& fz = mfz[fzi];
265  forAll(fz, i)
266  {
267  faceLabels.append(fz[i]);
268  compactZones.append(faceZoneToCompactZone[fzi]);
269  }
270  }
271  }
272 
273 
274  // Addressing engine for all faces
275  uindirectPrimitivePatch allBoundary
276  (
277  UIndirectList<face>(mesh.faces(), faceLabels),
278  mesh.points()
279  );
280 
281 
282  // Find correspondence to master points
283  labelList pointToGlobal;
284  labelList uniqueMeshPoints;
286  (
287  allBoundary.meshPoints(),
288  allBoundary.meshPointMap(),
289  pointToGlobal,
290  uniqueMeshPoints
291  );
292 
293  // Gather all unique points on master
294  List<pointField> gatheredPoints(Pstream::nProcs());
295  gatheredPoints[Pstream::myProcNo()] = pointField
296  (
297  mesh.points(),
298  uniqueMeshPoints
299  );
300  Pstream::gatherList(gatheredPoints);
301 
302  // Gather all faces
303  List<faceList> gatheredFaces(Pstream::nProcs());
304  gatheredFaces[Pstream::myProcNo()] = allBoundary.localFaces();
305  forAll(gatheredFaces[Pstream::myProcNo()], i)
306  {
307  inplaceRenumber(pointToGlobal, gatheredFaces[Pstream::myProcNo()][i]);
308  }
309  Pstream::gatherList(gatheredFaces);
310 
311  // Gather all ZoneIDs
312  List<labelList> gatheredZones(Pstream::nProcs());
313  gatheredZones[Pstream::myProcNo()] = move(compactZones);
314  Pstream::gatherList(gatheredZones);
315 
316  // On master combine all points, faces, zones
317  if (Pstream::master())
318  {
319  pointField allPoints = ListListOps::combine<pointField>
320  (
321  gatheredPoints,
323  );
324  gatheredPoints.clear();
325 
326  faceList allFaces = ListListOps::combine<faceList>
327  (
328  gatheredFaces,
330  );
331  gatheredFaces.clear();
332 
333  labelList allZones = ListListOps::combine<labelList>
334  (
335  gatheredZones,
337  );
338  gatheredZones.clear();
339 
340 
341  // Zones
342  surfZoneIdentifierList surfZones(compactZoneID.size());
343  forAllConstIter(HashTable<label>, compactZoneID, iter)
344  {
345  surfZones[iter()] = surfZoneIdentifier(iter.key(), iter());
346  Info<< "surfZone " << iter() << " : " << surfZones[iter()].name()
347  << endl;
348  }
349 
350  UnsortedMeshedSurface<face> unsortedFace
351  (
352  move(allPoints),
353  move(allFaces),
354  move(allZones),
355  move(surfZones)
356  );
357 
358 
359  MeshedSurface<face> sortedFace(unsortedFace);
360 
361  Info<< "Writing merged surface to "
362  << runTime.globalPath()/outFileName << endl;
363 
364  sortedFace.write(runTime.globalPath()/outFileName);
365  }
366 
367  Info<< "End\n" << endl;
368 
369  return 0;
370 }
371 
372 
373 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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:171
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 >> &zonesOrientedIndices)
Insert given oriented indices.
Definition: faceZoneList.C:45
Named list of face indices representing a sub-set of the mesh faces.
Definition: faceZone.H:66
A class for handling file names.
Definition: fileName.H:82
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
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 patch set corresponding to the given names.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:443
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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:258
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:267
labelHashSet includePatches
Foam::argList args(argc, argv)