polyDualMesh.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  polyDualMesh
26 
27 Description
28  Calculates the dual of a polyMesh. Adheres to all the feature and patch
29  edges.
30 
31 Usage
32  \b polyDualMesh featureAngle
33 
34  Detects any boundary edge > angle and creates multiple boundary faces
35  for it. Normal behaviour is to have each point become a cell
36  (1.5 behaviour)
37 
38  Options:
39  - \par -concaveMultiCells
40  Creates multiple cells for each point on a concave edge. Might limit
41  the amount of distortion on some meshes.
42 
43  - \par -splitAllFaces
44  Normally only constructs a single face between two cells. This single
45  face might be too distorted. splitAllFaces will create a single face for
46  every original cell the face passes through. The mesh will thus have
47  multiple faces in between two cells! (so is not strictly
48  upper-triangular anymore - checkMesh will complain)
49 
50  - \par -doNotPreserveFaceZones:
51  By default all faceZones are preserved by marking all faces, edges and
52  points on them as features. The -doNotPreserveFaceZones disables this
53  behaviour.
54 
55  Note:
56  It is just a driver for meshDualiser. Substitute your own
57  simpleMarkFeatures to have different behaviour.
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #include "argList.H"
62 #include "Time.H"
63 #include "fvMesh.H"
64 #include "polyTopoChange.H"
65 #include "polyTopoChangeMap.H"
66 #include "PackedBoolList.H"
67 #include "meshTools.H"
68 #include "OFstream.H"
69 #include "meshDualiser.H"
70 
71 using namespace Foam;
72 
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75 // Naive feature detection. All boundary edges with angle > featureAngle become
76 // feature edges. All points on feature edges become feature points. All
77 // boundary faces become feature faces.
78 void simpleMarkFeatures
79 (
80  const polyMesh& mesh,
81  const PackedBoolList& isBoundaryEdge,
82  const scalar featureAngle,
83  const bool concaveMultiCells,
84  const bool doNotPreserveFaceZones,
85 
86  labelList& featureFaces,
87  labelList& featureEdges,
88  labelList& singleCellFeaturePoints,
89  labelList& multiCellFeaturePoints
90 )
91 {
92  const scalar minCos = Foam::cos(featureAngle);
93 
94  const polyBoundaryMesh& patches = mesh.boundaryMesh();
95 
96  // Working sets
97  labelHashSet featureEdgeSet;
98  labelHashSet singleCellFeaturePointSet;
99  labelHashSet multiCellFeaturePointSet;
100 
101 
102  // 1. Mark all edges between patches
103  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
104 
106  {
107  const polyPatch& pp = patches[patchi];
108  const labelList& meshEdges = pp.meshEdges();
109 
110  // All patch corner edges. These need to be feature points & edges!
111  for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
112  {
113  label meshEdgeI = meshEdges[edgeI];
114  featureEdgeSet.insert(meshEdgeI);
115  singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][0]);
116  singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][1]);
117  }
118  }
119 
120 
121 
122  // 2. Mark all geometric feature edges
123  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
124  // Make distinction between convex features where the boundary point becomes
125  // a single cell and concave features where the boundary point becomes
126  // multiple 'half' cells.
127 
128  // Addressing for all outside faces
129  primitivePatch allBoundary
130  (
132  (
133  mesh.faces(),
134  mesh.nFaces()-mesh.nInternalFaces(),
135  mesh.nInternalFaces()
136  ),
137  mesh.points()
138  );
139 
140  // Check for non-manifold points (surface pinched at point)
141  allBoundary.checkPointManifold(false, &singleCellFeaturePointSet);
142 
143  // Check for non-manifold edges (surface pinched at edge)
144  const labelListList& edgeFaces = allBoundary.edgeFaces();
145  const labelList& meshPoints = allBoundary.meshPoints();
146 
147  forAll(edgeFaces, edgeI)
148  {
149  const labelList& eFaces = edgeFaces[edgeI];
150 
151  if (eFaces.size() > 2)
152  {
153  const edge& e = allBoundary.edges()[edgeI];
154 
155  // Info<< "Detected non-manifold boundary edge:" << edgeI
156  // << " coords:"
157  // << allBoundary.points()[meshPoints[e[0]]]
158  // << allBoundary.points()[meshPoints[e[1]]] << endl;
159 
160  singleCellFeaturePointSet.insert(meshPoints[e[0]]);
161  singleCellFeaturePointSet.insert(meshPoints[e[1]]);
162  }
163  }
164 
165  // Check for features.
166  forAll(edgeFaces, edgeI)
167  {
168  const labelList& eFaces = edgeFaces[edgeI];
169 
170  if (eFaces.size() == 2)
171  {
172  label f0 = eFaces[0];
173  label f1 = eFaces[1];
174 
175  // check angle
176  const vector& n0 = allBoundary.faceNormals()[f0];
177  const vector& n1 = allBoundary.faceNormals()[f1];
178 
179  if ((n0 & n1) < minCos)
180  {
181  const edge& e = allBoundary.edges()[edgeI];
182  label v0 = meshPoints[e[0]];
183  label v1 = meshPoints[e[1]];
184 
185  label meshEdgeI = meshTools::findEdge(mesh, v0, v1);
186  featureEdgeSet.insert(meshEdgeI);
187 
188  // Check if convex or concave by looking at angle
189  // between face centres and normal
190  vector c1c0
191  (
192  allBoundary[f1].centre(allBoundary.points())
193  - allBoundary[f0].centre(allBoundary.points())
194  );
195 
196  if (concaveMultiCells && (c1c0 & n0) > small)
197  {
198  // Found concave edge. Make into multiCell features
199  Info<< "Detected concave feature edge:" << edgeI
200  << " cos:" << (c1c0 & n0)
201  << " coords:"
202  << allBoundary.points()[v0]
203  << allBoundary.points()[v1]
204  << endl;
205 
206  singleCellFeaturePointSet.erase(v0);
207  multiCellFeaturePointSet.insert(v0);
208  singleCellFeaturePointSet.erase(v1);
209  multiCellFeaturePointSet.insert(v1);
210  }
211  else
212  {
213  // Convex. singleCell feature.
214  if (!multiCellFeaturePointSet.found(v0))
215  {
216  singleCellFeaturePointSet.insert(v0);
217  }
218  if (!multiCellFeaturePointSet.found(v1))
219  {
220  singleCellFeaturePointSet.insert(v1);
221  }
222  }
223  }
224  }
225  }
226 
227 
228  // 3. Mark all feature faces
229  // ~~~~~~~~~~~~~~~~~~~~~~~~~
230 
231  // Face centres that need inclusion in the dual mesh
232  labelHashSet featureFaceSet(mesh.nFaces()-mesh.nInternalFaces());
233  // A. boundary faces.
234  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
235  {
236  featureFaceSet.insert(facei);
237  }
238 
239  // B. face zones.
240  const faceZoneList& faceZones = mesh.faceZones();
241 
242  if (doNotPreserveFaceZones)
243  {
244  if (faceZones.size() > 0)
245  {
247  << "Detected " << faceZones.size()
248  << " faceZones. These will not be preserved."
249  << endl;
250  }
251  }
252  else
253  {
254  if (faceZones.size() > 0)
255  {
256  Info<< "Detected " << faceZones.size()
257  << " faceZones. Preserving these by marking their"
258  << " points, edges and faces as features." << endl;
259  }
260 
261  forAll(faceZones, zoneI)
262  {
263  const faceZone& fz = faceZones[zoneI];
264 
265  Info<< "Inserting all faces in faceZone " << fz.name()
266  << " as features." << endl;
267 
268  forAll(fz, i)
269  {
270  label facei = fz[i];
271  const face& f = mesh.faces()[facei];
272  const labelList& fEdges = mesh.faceEdges()[facei];
273 
274  featureFaceSet.insert(facei);
275  forAll(f, fp)
276  {
277  // Mark point as multi cell point (since both sides of
278  // face should have different cells)
279  singleCellFeaturePointSet.erase(f[fp]);
280  multiCellFeaturePointSet.insert(f[fp]);
281 
282  // Make sure there are points on the edges.
283  featureEdgeSet.insert(fEdges[fp]);
284  }
285  }
286  }
287  }
288 
289  // Transfer to arguments
290  featureFaces = featureFaceSet.toc();
291  featureEdges = featureEdgeSet.toc();
292  singleCellFeaturePoints = singleCellFeaturePointSet.toc();
293  multiCellFeaturePoints = multiCellFeaturePointSet.toc();
294 }
295 
296 
297 // Dump features to .obj files
298 void dumpFeatures
299 (
300  const polyMesh& mesh,
301  const labelList& featureFaces,
302  const labelList& featureEdges,
303  const labelList& singleCellFeaturePoints,
304  const labelList& multiCellFeaturePoints
305 )
306 {
307  {
308  OFstream str("featureFaces.obj");
309  Info<< "Dumping centres of featureFaces to obj file " << str.name()
310  << endl;
311  forAll(featureFaces, i)
312  {
313  meshTools::writeOBJ(str, mesh.faceCentres()[featureFaces[i]]);
314  }
315  }
316  {
317  OFstream str("featureEdges.obj");
318  Info<< "Dumping featureEdges to obj file " << str.name() << endl;
319  label vertI = 0;
320 
321  forAll(featureEdges, i)
322  {
323  const edge& e = mesh.edges()[featureEdges[i]];
324  meshTools::writeOBJ(str, mesh.points()[e[0]]);
325  vertI++;
326  meshTools::writeOBJ(str, mesh.points()[e[1]]);
327  vertI++;
328  str<< "l " << vertI-1 << ' ' << vertI << nl;
329  }
330  }
331  {
332  OFstream str("singleCellFeaturePoints.obj");
333  Info<< "Dumping featurePoints that become a single cell to obj file "
334  << str.name() << endl;
335  forAll(singleCellFeaturePoints, i)
336  {
337  meshTools::writeOBJ(str, mesh.points()[singleCellFeaturePoints[i]]);
338  }
339  }
340  {
341  OFstream str("multiCellFeaturePoints.obj");
342  Info<< "Dumping featurePoints that become multiple cells to obj file "
343  << str.name() << endl;
344  forAll(multiCellFeaturePoints, i)
345  {
346  meshTools::writeOBJ(str, mesh.points()[multiCellFeaturePoints[i]]);
347  }
348  }
349 }
350 
351 
352 int main(int argc, char *argv[])
353 {
354  #include "addOverwriteOption.H"
356 
357  argList::validArgs.append("featureAngle [0-180]");
359  (
360  "splitAllFaces",
361  "have multiple faces in between cells"
362  );
364  (
365  "concaveMultiCells",
366  "split cells on concave boundary edges into multiple cells"
367  );
369  (
370  "doNotPreserveFaceZones",
371  "disable the default behaviour of preserving faceZones by having"
372  " multiple faces in between cells"
373  );
374 
375  #include "setRootCase.H"
376  #include "createTime.H"
377  #include "createMeshNoChangers.H"
378 
379  const word oldInstance = mesh.pointsInstance();
380 
381  // Mark boundary edges and points.
382  // (Note: in 1.4.2 we can use the built-in mesh point ordering
383  // facility instead)
384  PackedBoolList isBoundaryEdge(mesh.nEdges());
385  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
386  {
387  const labelList& fEdges = mesh.faceEdges()[facei];
388 
389  forAll(fEdges, i)
390  {
391  isBoundaryEdge.set(fEdges[i], 1);
392  }
393  }
394 
395  const scalar featureAngle = degToRad(args.argRead<scalar>(1));
396  const scalar minCos = Foam::cos(featureAngle);
397 
398  Info<< "Feature:" << radToDeg(featureAngle) << endl
399  << "minCos :" << minCos << endl
400  << endl;
401 
402 
403  const bool splitAllFaces = args.optionFound("splitAllFaces");
404  if (splitAllFaces)
405  {
406  Info<< "Splitting all internal faces to create multiple faces"
407  << " between two cells." << nl
408  << endl;
409  }
410 
411  const bool overwrite = args.optionFound("overwrite");
412  const bool doNotPreserveFaceZones = args.optionFound
413  (
414  "doNotPreserveFaceZones"
415  );
416  const bool concaveMultiCells = args.optionFound("concaveMultiCells");
417  if (concaveMultiCells)
418  {
419  Info<< "Generating multiple cells for points on concave feature edges."
420  << nl << endl;
421  }
422 
423 
424  mesh.cellZones().clear();
425  mesh.pointZones().clear();
426 
427  if (doNotPreserveFaceZones)
428  {
429  mesh.faceZones().clear();
430  }
431 
432  // Face(centre)s that need inclusion in the dual mesh
433  labelList featureFaces;
434  // Edge(centre)s ,,
435  labelList featureEdges;
436  // Points (that become a single cell) that need inclusion in the dual mesh
437  labelList singleCellFeaturePoints;
438  // Points (that become a multiple cells) ,,
439  labelList multiCellFeaturePoints;
440 
441  // Sample implementation of feature detection.
442  simpleMarkFeatures
443  (
444  mesh,
445  isBoundaryEdge,
446  featureAngle,
447  concaveMultiCells,
448  doNotPreserveFaceZones,
449 
450  featureFaces,
451  featureEdges,
452  singleCellFeaturePoints,
453  multiCellFeaturePoints
454  );
455 
456  mesh.faceZones().clear();
457 
458  // If we want to split all polyMesh faces into one dualface per cell
459  // we are passing through we also need a point
460  // at the polyMesh facecentre and edgemid of the faces we want to
461  // split.
462  if (splitAllFaces)
463  {
464  featureEdges = identityMap(mesh.nEdges());
465  featureFaces = identityMap(mesh.nFaces());
466  }
467 
468  // Write obj files for debugging
469  dumpFeatures
470  (
471  mesh,
472  featureFaces,
473  featureEdges,
474  singleCellFeaturePoints,
475  multiCellFeaturePoints
476  );
477 
478 
479  // Topo change container
480  polyTopoChange meshMod(mesh.boundaryMesh().size());
481 
482  // Mesh dualiser engine
483  meshDualiser dualMaker(mesh);
484 
485  // Insert all commands into polyTopoChange to create dual of mesh. This does
486  // all the hard work.
487  dualMaker.setRefinement
488  (
489  splitAllFaces,
490  featureFaces,
491  featureEdges,
492  singleCellFeaturePoints,
493  multiCellFeaturePoints,
494  meshMod
495  );
496 
497  // Create mesh, return map from old to new mesh.
498  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
499 
500  // Update mesh objects
501  mesh.topoChange(map);
502 
503  if (!overwrite)
504  {
505  runTime++;
506  }
507  else
508  {
509  mesh.setInstance(oldInstance);
510  }
511 
512  Info<< "Writing dual mesh to " << runTime.name() << endl;
513 
514  mesh.write();
515 
516  Info<< "End\n" << endl;
517 
518  return 0;
519 }
520 
521 
522 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:445
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:227
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Output to file stream.
Definition: OFstream.H:86
A bit-packed bool list.
void set(const PackedList< 1 > &)
Set specified bits.
label nEdges() const
Return number of edges in patch.
label nInternalEdges() const
Number of internal edges.
A List obtained as a section of another List.
Definition: SubList.H:56
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void clear()
Clear the zones.
Definition: ZoneList.C:236
const word & name() const
Return name.
Definition: Zone.H:184
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 void noParallel()
Remove the parallel options.
Definition: argList.C:175
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:183
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Creates dual of polyMesh. Every point becomes a cell (or multiple cells for feature points),...
Definition: meshDualiser.H:67
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const pointZoneList & pointZones() const
Return point zones.
Definition: polyMesh.H:440
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:965
const faceZoneList & faceZones() const
Return face zones.
Definition: polyMesh.H:446
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const labelList & meshEdges() const
Return global edge index for local edges.
Definition: polyPatch.C:327
Direct mesh changes based on v1.3 polyTopoChange syntax.
label nEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
label nInternalFaces() const
label nFaces() const
const labelListList & faceEdges() const
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
label findEdge(const edgeList &edges, const labelList &candidates, const label v0, const label v1)
Return edge among candidates that uses the two vertices.
Definition: meshTools.C:337
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
labelList f(nPoints)
Foam::argList args(argc, argv)