Loading...
Searching...
No Matches
polyDualMeshApp.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) 2016,2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 polyDualMesh
29
30Group
31 grpMeshManipulationUtilities
32
33Description
34 Creates the dual of a polyMesh, adhering to all the feature and patch edges.
35
36Usage
37 \b polyDualMesh featureAngle
38
39 Detects any boundary edge > angle and creates multiple boundary faces
40 for it. Normal behaviour is to have each point become a cell
41 (1.5 behaviour)
42
43 Options:
44 - \par -concaveMultiCells
45 Creates multiple cells for each point on a concave edge. Might limit
46 the amount of distortion on some meshes.
47
48 - \par -splitAllFaces
49 Normally only constructs a single face between two cells. This single
50 face might be too distorted. splitAllFaces will create a single face for
51 every original cell the face passes through. The mesh will thus have
52 multiple faces in between two cells! (so is not strictly
53 upper-triangular anymore - checkMesh will complain)
54
55 - \par -doNotPreserveFaceZones:
56 By default all faceZones are preserved by marking all faces, edges and
57 points on them as features. The -doNotPreserveFaceZones disables this
58 behaviour.
59
60Note
61 It is just a driver for meshDualiser. Substitute your own simpleMarkFeatures
62 to have different behaviour.
63
64\*---------------------------------------------------------------------------*/
65
66#include "argList.H"
67#include "Time.H"
68#include "fvMesh.H"
69#include "unitConversion.H"
70#include "polyTopoChange.H"
71#include "mapPolyMesh.H"
72#include "bitSet.H"
73#include "meshTools.H"
74#include "OFstream.H"
75#include "meshDualiser.H"
76#include "ReadFields.H"
77#include "volFields.H"
78#include "surfaceFields.H"
79#include "topoSet.H"
80#include "processorMeshes.H"
81
82using namespace Foam;
83
84// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85
86// Naive feature detection. All boundary edges with angle > featureAngle become
87// feature edges. All points on feature edges become feature points. All
88// boundary faces become feature faces.
89void simpleMarkFeatures
90(
91 const polyMesh& mesh,
92 const bitSet& isBoundaryEdge,
93 const scalar featureAngle,
94 const bool concaveMultiCells,
95 const bool doNotPreserveFaceZones,
96
97 labelList& featureFaces,
98 labelList& featureEdges,
99 labelList& singleCellFeaturePoints,
100 labelList& multiCellFeaturePoints
101)
102{
103 scalar minCos = Foam::cos(degToRad(featureAngle));
104
105 const polyBoundaryMesh& patches = mesh.boundaryMesh();
106
107 // Working sets
108 labelHashSet featureEdgeSet;
109 labelHashSet singleCellFeaturePointSet;
110 labelHashSet multiCellFeaturePointSet;
111
112
113 // 1. Mark all edges between patches
114 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
115
116 forAll(patches, patchi)
117 {
118 const polyPatch& pp = patches[patchi];
119 const labelList& meshEdges = pp.meshEdges();
120
121 // All patch corner edges. These need to be feature points & edges!
122 for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
123 {
124 label meshEdgeI = meshEdges[edgeI];
125 featureEdgeSet.insert(meshEdgeI);
126 singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][0]);
127 singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][1]);
128 }
129 }
130
131
132
133 // 2. Mark all geometric feature edges
134 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
135 // Make distinction between convex features where the boundary point becomes
136 // a single cell and concave features where the boundary point becomes
137 // multiple 'half' cells.
138
139 // Addressing for all outside faces
140 primitivePatch allBoundary(patches.faces(), mesh.points());
141
142 // Check for non-manifold points (surface pinched at point)
143 allBoundary.checkPointManifold(false, &singleCellFeaturePointSet);
144
145 // Check for non-manifold edges (surface pinched at edge)
146 const labelListList& edgeFaces = allBoundary.edgeFaces();
147 const labelList& meshPoints = allBoundary.meshPoints();
148
149 forAll(edgeFaces, edgeI)
150 {
151 const labelList& eFaces = edgeFaces[edgeI];
152
153 if (eFaces.size() > 2)
154 {
155 const edge& e = allBoundary.edges()[edgeI];
156
157 //Info<< "Detected non-manifold boundary edge:" << edgeI
158 // << " coords:"
159 // << allBoundary.points()[meshPoints[e[0]]]
160 // << allBoundary.points()[meshPoints[e[1]]] << endl;
161
162 singleCellFeaturePointSet.insert(meshPoints[e[0]]);
163 singleCellFeaturePointSet.insert(meshPoints[e[1]]);
164 }
165 }
166
167 // Check for features.
168 forAll(edgeFaces, edgeI)
169 {
170 const labelList& eFaces = edgeFaces[edgeI];
171
172 if (eFaces.size() == 2)
173 {
174 label f0 = eFaces[0];
175 label f1 = eFaces[1];
176
177 // check angle
178 const vector& n0 = allBoundary.faceNormals()[f0];
179 const vector& n1 = allBoundary.faceNormals()[f1];
180
181 if ((n0 & n1) < minCos)
182 {
183 const edge& e = allBoundary.edges()[edgeI];
184 label v0 = meshPoints[e[0]];
185 label v1 = meshPoints[e[1]];
186
187 label meshEdgeI = meshTools::findEdge(mesh, v0, v1);
188 featureEdgeSet.insert(meshEdgeI);
189
190 // Check if convex or concave by looking at angle
191 // between face centres and normal
192 vector c1c0
193 (
194 allBoundary[f1].centre(allBoundary.points())
195 - allBoundary[f0].centre(allBoundary.points())
196 );
197
198 if (concaveMultiCells && (c1c0 & n0) > SMALL)
199 {
200 // Found concave edge. Make into multiCell features
201 Info<< "Detected concave feature edge:" << edgeI
202 << " cos:" << (c1c0 & n0)
203 << " coords:"
204 << allBoundary.points()[v0]
205 << allBoundary.points()[v1]
206 << endl;
207
208 singleCellFeaturePointSet.erase(v0);
209 multiCellFeaturePointSet.insert(v0);
210 singleCellFeaturePointSet.erase(v1);
211 multiCellFeaturePointSet.insert(v1);
212 }
213 else
214 {
215 // Convex. singleCell feature.
216 if (!multiCellFeaturePointSet.found(v0))
217 {
218 singleCellFeaturePointSet.insert(v0);
219 }
220 if (!multiCellFeaturePointSet.found(v1))
221 {
222 singleCellFeaturePointSet.insert(v1);
223 }
224 }
225 }
226 }
227 }
228
229
230 // 3. Mark all feature faces
231 // ~~~~~~~~~~~~~~~~~~~~~~~~~
232
233 // Face centres that need inclusion in the dual mesh
234 labelHashSet featureFaceSet;
235 featureFaceSet.reserve(mesh.nBoundaryFaces());
236
237 // A. boundary faces.
238 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
239 {
240 featureFaceSet.insert(facei);
241 }
242
243 // B. face zones.
244 const faceZoneMesh& faceZones = mesh.faceZones();
245
246 if (doNotPreserveFaceZones)
247 {
248 if (faceZones.size() > 0)
249 {
251 << "Detected " << faceZones.size()
252 << " faceZones. These will not be preserved."
253 << endl;
254 }
255 }
256 else
257 {
258 if (faceZones.size() > 0)
259 {
260 Info<< "Detected " << faceZones.size()
261 << " faceZones. Preserving these by marking their"
262 << " points, edges and faces as features." << endl;
263 }
264
265 forAll(faceZones, zoneI)
266 {
267 const faceZone& fz = faceZones[zoneI];
268
269 Info<< "Inserting all faces in faceZone " << fz.name()
270 << " as features." << endl;
271
272 forAll(fz, i)
273 {
274 label facei = fz[i];
275 const face& f = mesh.faces()[facei];
276 const labelList& fEdges = mesh.faceEdges()[facei];
277
278 featureFaceSet.insert(facei);
279 forAll(f, fp)
280 {
281 // Mark point as multi cell point (since both sides of
282 // face should have different cells)
283 singleCellFeaturePointSet.erase(f[fp]);
284 multiCellFeaturePointSet.insert(f[fp]);
285
286 // Make sure there are points on the edges.
287 featureEdgeSet.insert(fEdges[fp]);
288 }
289 }
290 }
291 }
292
293 // Transfer to arguments
294 featureFaces = featureFaceSet.toc();
295 featureEdges = featureEdgeSet.toc();
296 singleCellFeaturePoints = singleCellFeaturePointSet.toc();
297 multiCellFeaturePoints = multiCellFeaturePointSet.toc();
298}
299
300
301// Dump features to .obj files
302void dumpFeatures
303(
304 const polyMesh& mesh,
305 const labelList& featureFaces,
306 const labelList& featureEdges,
307 const labelList& singleCellFeaturePoints,
308 const labelList& multiCellFeaturePoints
309)
310{
311 {
312 OFstream str("featureFaces.obj");
313 Info<< "Dumping centres of featureFaces to obj file " << str.name()
314 << endl;
315 forAll(featureFaces, i)
316 {
317 meshTools::writeOBJ(str, mesh.faceCentres()[featureFaces[i]]);
318 }
319 }
320 {
321 OFstream str("featureEdges.obj");
322 Info<< "Dumping featureEdges to obj file " << str.name() << endl;
323 label vertI = 0;
324
325 forAll(featureEdges, i)
326 {
327 const edge& e = mesh.edges()[featureEdges[i]];
328 meshTools::writeOBJ(str, mesh.points()[e[0]]);
329 vertI++;
330 meshTools::writeOBJ(str, mesh.points()[e[1]]);
331 vertI++;
332 str<< "l " << vertI-1 << ' ' << vertI << nl;
333 }
334 }
335 {
336 OFstream str("singleCellFeaturePoints.obj");
337 Info<< "Dumping featurePoints that become a single cell to obj file "
338 << str.name() << endl;
339 forAll(singleCellFeaturePoints, i)
340 {
341 meshTools::writeOBJ(str, mesh.points()[singleCellFeaturePoints[i]]);
342 }
343 }
344 {
345 OFstream str("multiCellFeaturePoints.obj");
346 Info<< "Dumping featurePoints that become multiple cells to obj file "
347 << str.name() << endl;
348 forAll(multiCellFeaturePoints, i)
349 {
350 meshTools::writeOBJ(str, mesh.points()[multiCellFeaturePoints[i]]);
351 }
352 }
353}
354
355
356int main(int argc, char *argv[])
357{
359 (
360 "Creates the dual of a polyMesh,"
361 " adhering to all the feature and patch edges."
362 );
363
364 #include "addOverwriteOption.H"
366
368 (
369 "featureAngle",
370 "in degrees [0-180]"
371 );
372
374 (
375 "splitAllFaces",
376 "Have multiple faces in between cells"
377 );
379 (
380 "concaveMultiCells",
381 "Split cells on concave boundary edges into multiple cells"
382 );
384 (
385 "doNotPreserveFaceZones",
386 "Disable the default behaviour of preserving faceZones by having"
387 " multiple faces in between cells"
388 );
389
390 #include "setRootCase.H"
391 #include "createTime.H"
392 #include "createNamedMesh.H"
393
394 const word oldInstance = mesh.pointsInstance();
395
396 // Mark boundary edges and points.
397 // (Note: in 1.4.2 we can use the built-in mesh point ordering
398 // facility instead)
399 bitSet isBoundaryEdge(mesh.nEdges());
400 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); facei++)
401 {
402 const labelList& fEdges = mesh.faceEdges()[facei];
403
404 forAll(fEdges, i)
405 {
406 isBoundaryEdge.set(fEdges[i]);
407 }
408 }
409
410 const scalar featureAngle = args.get<scalar>(1);
411 const scalar minCos = Foam::cos(degToRad(featureAngle));
412
413 Info<< "Feature:" << featureAngle << endl
414 << "minCos :" << minCos << endl
415 << endl;
416
417
418 const bool splitAllFaces = args.found("splitAllFaces");
419 if (splitAllFaces)
420 {
421 Info<< "Splitting all internal faces to create multiple faces"
422 << " between two cells." << nl
423 << endl;
424 }
425
426 const bool overwrite = args.found("overwrite");
427 const bool doNotPreserveFaceZones = args.found("doNotPreserveFaceZones");
428 const bool concaveMultiCells = args.found("concaveMultiCells");
429 if (concaveMultiCells)
430 {
431 Info<< "Generating multiple cells for points on concave feature edges."
432 << nl << endl;
433 }
434
435
436 // Face(centre)s that need inclusion in the dual mesh
437 labelList featureFaces;
438 // Edge(centre)s ,,
439 labelList featureEdges;
440 // Points (that become a single cell) that need inclusion in the dual mesh
441 labelList singleCellFeaturePoints;
442 // Points (that become a multiple cells) ,,
443 labelList multiCellFeaturePoints;
444
445 // Sample implementation of feature detection.
446 simpleMarkFeatures
447 (
448 mesh,
449 isBoundaryEdge,
450 featureAngle,
451 concaveMultiCells,
452 doNotPreserveFaceZones,
453
454 featureFaces,
455 featureEdges,
456 singleCellFeaturePoints,
457 multiCellFeaturePoints
458 );
459
460 // If we want to split all polyMesh faces into one dualface per cell
461 // we are passing through we also need a point
462 // at the polyMesh facecentre and edgemid of the faces we want to
463 // split.
464 if (splitAllFaces)
465 {
466 featureEdges = identity(mesh.nEdges());
467 featureFaces = identity(mesh.nFaces());
468 }
469
470 // Write obj files for debugging
471 dumpFeatures
472 (
473 mesh,
474 featureFaces,
475 featureEdges,
476 singleCellFeaturePoints,
477 multiCellFeaturePoints
478 );
479
480
481
482 // Read objects in time directory
483 IOobjectList objects(mesh, runTime.timeName());
484
485 // Read vol fields.
487 ReadFields(mesh, objects, vsFlds);
488
490 ReadFields(mesh, objects, vvFlds);
491
493 ReadFields(mesh, objects, vstFlds);
494
496 ReadFields(mesh, objects, vsymtFlds);
497
499 ReadFields(mesh, objects, vtFlds);
500
501 // Read surface fields.
503 ReadFields(mesh, objects, ssFlds);
504
506 ReadFields(mesh, objects, svFlds);
507
509 ReadFields(mesh, objects, sstFlds);
510
512 ReadFields(mesh, objects, ssymtFlds);
513
515 ReadFields(mesh, objects, stFlds);
516
517
518 // Topo change container
519 polyTopoChange meshMod(mesh.boundaryMesh().size());
520
521 // Mesh dualiser engine
522 meshDualiser dualMaker(mesh);
523
524 // Insert all commands into polyTopoChange to create dual of mesh. This does
525 // all the hard work.
526 dualMaker.setRefinement
527 (
528 splitAllFaces,
529 featureFaces,
530 featureEdges,
531 singleCellFeaturePoints,
532 multiCellFeaturePoints,
533 meshMod
534 );
535
536 // Create mesh, return map from old to new mesh.
537 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
538
539 // Update fields
540 mesh.updateMesh(map());
541
542 // Optionally inflate mesh
543 if (map().hasMotionPoints())
544 {
545 mesh.movePoints(map().preMotionPoints());
546 }
547
548 if (!overwrite)
549 {
550 ++runTime;
551 }
552 else
553 {
554 mesh.setInstance(oldInstance);
555 }
556
557 Info<< "Writing dual mesh to " << runTime.timeName() << endl;
558
559 mesh.write();
562
563 Info<< "End\n" << endl;
564
565 return 0;
566}
567
568
569// ************************************************************************* //
Field reading functions for post-processing utilities.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > toc() const
The table of contents (the keys) in unsorted order.
Definition HashTable.C:141
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition HashTable.C:489
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition argList.C:366
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition argList.C:389
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Creates dual of polyMesh. Every point becomes a cell (or multiple cells for feature points),...
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Direct mesh changes based on v1.3 polyTopoChange syntax.
static void removeFiles(const polyMesh &mesh)
Helper: remove all procAddressing files from mesh instance.
static void removeFiles(const polyMesh &)
Helper: remove all sets files from mesh instance.
Definition topoSet.C:693
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
const word & name() const noexcept
The zone name.
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
#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:352
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
Namespace for OpenFOAM.
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Read Geometric fields of templated type.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with faceZone content on a polyMesh.
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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...
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
Foam::argList args(argc, argv)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.
Unit conversion functions.