1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2004-2011 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
25 Calculate the dual of a polyMesh. Adheres to all the feature and patch
30 - polyDualMesh featureAngle
32 Detects any boundary edge > angle and creates multiple boundary faces
33 for it. Normal behaviour is to have each point become a cell
36 \param -concaveMultiCells
37 Creates multiple cells for each point on a concave edge. Might limit
38 the amount of distortion on some meshes.
41 Normally only constructs a single face between two cells. This single face
42 might be too distorted. splitAllFaces will create a single face for every
43 original cell the face passes through. The mesh will thus have
44 multiple faces inbetween two cells! (so is not strictly upper-triangular
45 anymore - checkMesh will complain)
47 \param -doNotPreserveFaceZones:
48 By default all faceZones are preserved by marking all faces, edges and
49 points on them as features. The -doNotPreserveFaceZones disables this
52 Note: is just a driver for meshDualiser. Substitute your own
53 simpleMarkFeatures to have different behaviour.
55 \*---------------------------------------------------------------------------*/
59 #include "timeSelector.H"
61 #include "unitConversion.H"
62 #include "polyTopoChange.H"
63 #include "mapPolyMesh.H"
64 #include "PackedBoolList.H"
65 #include "meshTools.H"
67 #include "meshDualiser.H"
68 #include "ReadFields.H"
69 #include "volFields.H"
70 #include "surfaceFields.H"
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 // Naive feature detection. All boundary edges with angle > featureAngle become
77 // feature edges. All points on feature edges become feature points. All
78 // boundary faces become feature faces.
79 void simpleMarkFeatures
82 const PackedBoolList& isBoundaryEdge,
83 const scalar featureAngle,
84 const bool concaveMultiCells,
85 const bool doNotPreserveFaceZones,
87 labelList& featureFaces,
88 labelList& featureEdges,
89 labelList& singleCellFeaturePoints,
90 labelList& multiCellFeaturePoints
93 scalar minCos = Foam::cos(degToRad(featureAngle));
95 const polyBoundaryMesh& patches = mesh.boundaryMesh();
98 labelHashSet featureEdgeSet;
99 labelHashSet singleCellFeaturePointSet;
100 labelHashSet multiCellFeaturePointSet;
103 // 1. Mark all edges between patches
104 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
106 forAll(patches, patchI)
108 const polyPatch& pp = patches[patchI];
109 const labelList& meshEdges = pp.meshEdges();
111 // All patch corner edges. These need to be feature points & edges!
112 for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
114 label meshEdgeI = meshEdges[edgeI];
115 featureEdgeSet.insert(meshEdgeI);
116 singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][0]);
117 singleCellFeaturePointSet.insert(mesh.edges()[meshEdgeI][1]);
123 // 2. Mark all geometric feature edges
124 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
125 // Make distinction between convex features where the boundary point becomes
126 // a single cell and concave features where the boundary point becomes
127 // multiple 'half' cells.
129 // Addressing for all outside faces
130 primitivePatch allBoundary
135 mesh.nFaces()-mesh.nInternalFaces(),
136 mesh.nInternalFaces()
141 // Check for non-manifold points (surface pinched at point)
142 allBoundary.checkPointManifold(false, &singleCellFeaturePointSet);
144 // Check for non-manifold edges (surface pinched at edge)
145 const labelListList& edgeFaces = allBoundary.edgeFaces();
146 const labelList& meshPoints = allBoundary.meshPoints();
148 forAll(edgeFaces, edgeI)
150 const labelList& eFaces = edgeFaces[edgeI];
152 if (eFaces.size() > 2)
154 const edge& e = allBoundary.edges()[edgeI];
156 //Info<< "Detected non-manifold boundary edge:" << edgeI
158 // << allBoundary.points()[meshPoints[e[0]]]
159 // << allBoundary.points()[meshPoints[e[1]]] << endl;
161 singleCellFeaturePointSet.insert(meshPoints[e[0]]);
162 singleCellFeaturePointSet.insert(meshPoints[e[1]]);
166 // Check for features.
167 forAll(edgeFaces, edgeI)
169 const labelList& eFaces = edgeFaces[edgeI];
171 if (eFaces.size() == 2)
173 label f0 = eFaces[0];
174 label f1 = eFaces[1];
177 const vector& n0 = allBoundary.faceNormals()[f0];
178 const vector& n1 = allBoundary.faceNormals()[f1];
180 if ((n0 & n1) < minCos)
182 const edge& e = allBoundary.edges()[edgeI];
183 label v0 = meshPoints[e[0]];
184 label v1 = meshPoints[e[1]];
186 label meshEdgeI = meshTools::findEdge(mesh, v0, v1);
187 featureEdgeSet.insert(meshEdgeI);
189 // Check if convex or concave by looking at angle
190 // between face centres and normal
193 allBoundary[f1].centre(allBoundary.points())
194 - allBoundary[f0].centre(allBoundary.points())
197 if (concaveMultiCells && (c1c0 & n0) > SMALL)
199 // Found concave edge. Make into multiCell features
200 Info<< "Detected concave feature edge:" << edgeI
201 << " cos:" << (c1c0 & n0)
203 << allBoundary.points()[v0]
204 << allBoundary.points()[v1]
207 singleCellFeaturePointSet.erase(v0);
208 multiCellFeaturePointSet.insert(v0);
209 singleCellFeaturePointSet.erase(v1);
210 multiCellFeaturePointSet.insert(v1);
214 // Convex. singleCell feature.
215 if (!multiCellFeaturePointSet.found(v0))
217 singleCellFeaturePointSet.insert(v0);
219 if (!multiCellFeaturePointSet.found(v1))
221 singleCellFeaturePointSet.insert(v1);
229 // 3. Mark all feature faces
230 // ~~~~~~~~~~~~~~~~~~~~~~~~~
232 // Face centres that need inclusion in the dual mesh
233 labelHashSet featureFaceSet(mesh.nFaces()-mesh.nInternalFaces());
234 // A. boundary faces.
235 for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
237 featureFaceSet.insert(faceI);
241 const faceZoneMesh& faceZones = mesh.faceZones();
243 if (doNotPreserveFaceZones)
245 if (faceZones.size() > 0)
247 WarningIn("simpleMarkFeatures(..)")
248 << "Detected " << faceZones.size()
249 << " faceZones. These will not be preserved."
255 if (faceZones.size() > 0)
257 Info<< "Detected " << faceZones.size()
258 << " faceZones. Preserving these by marking their"
259 << " points, edges and faces as features." << endl;
262 forAll(faceZones, zoneI)
264 const faceZone& fz = faceZones[zoneI];
266 Info<< "Inserting all faces in faceZone " << fz.name()
267 << " as features." << endl;
272 const face& f = mesh.faces()[faceI];
273 const labelList& fEdges = mesh.faceEdges()[faceI];
275 featureFaceSet.insert(faceI);
278 // Mark point as multi cell point (since both sides of
279 // face should have different cells)
280 singleCellFeaturePointSet.erase(f[fp]);
281 multiCellFeaturePointSet.insert(f[fp]);
283 // Make sure there are points on the edges.
284 featureEdgeSet.insert(fEdges[fp]);
290 // Transfer to arguments
291 featureFaces = featureFaceSet.toc();
292 featureEdges = featureEdgeSet.toc();
293 singleCellFeaturePoints = singleCellFeaturePointSet.toc();
294 multiCellFeaturePoints = multiCellFeaturePointSet.toc();
298 // Dump features to .obj files
301 const polyMesh& mesh,
302 const labelList& featureFaces,
303 const labelList& featureEdges,
304 const labelList& singleCellFeaturePoints,
305 const labelList& multiCellFeaturePoints
309 OFstream str("featureFaces.obj");
310 Info<< "Dumping centres of featureFaces to obj file " << str.name()
312 forAll(featureFaces, i)
314 meshTools::writeOBJ(str, mesh.faceCentres()[featureFaces[i]]);
318 OFstream str("featureEdges.obj");
319 Info<< "Dumping featureEdges to obj file " << str.name() << endl;
322 forAll(featureEdges, i)
324 const edge& e = mesh.edges()[featureEdges[i]];
325 meshTools::writeOBJ(str, mesh.points()[e[0]]);
327 meshTools::writeOBJ(str, mesh.points()[e[1]]);
329 str<< "l " << vertI-1 << ' ' << vertI << nl;
333 OFstream str("singleCellFeaturePoints.obj");
334 Info<< "Dumping featurePoints that become a single cell to obj file "
335 << str.name() << endl;
336 forAll(singleCellFeaturePoints, i)
338 meshTools::writeOBJ(str, mesh.points()[singleCellFeaturePoints[i]]);
342 OFstream str("multiCellFeaturePoints.obj");
343 Info<< "Dumping featurePoints that become multiple cells to obj file "
344 << str.name() << endl;
345 forAll(multiCellFeaturePoints, i)
347 meshTools::writeOBJ(str, mesh.points()[multiCellFeaturePoints[i]]);
353 int main(int argc, char *argv[])
355 # include "addOverwriteOption.H"
356 argList::noParallel();
357 timeSelector::addOptions(true, false);
359 argList::validArgs.append("featureAngle [0-180]");
360 argList::addBoolOption
363 "have multiple faces inbetween cells"
365 argList::addBoolOption
368 "split cells on concave boundary edges into multiple cells"
370 argList::addBoolOption
372 "doNotPreserveFaceZones",
373 "disable the default behaviour of preserving faceZones by having"
374 " multiple faces inbetween cells"
377 # include "setRootCase.H"
378 # include "createTime.H"
380 instantList timeDirs = timeSelector::select0(runTime, args);
382 # include "createMesh.H"
384 const word oldInstance = mesh.pointsInstance();
386 // Mark boundary edges and points.
387 // (Note: in 1.4.2 we can use the built-in mesh point ordering
389 PackedBoolList isBoundaryEdge(mesh.nEdges());
390 for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
392 const labelList& fEdges = mesh.faceEdges()[faceI];
396 isBoundaryEdge.set(fEdges[i], 1);
400 const scalar featureAngle = args.argRead<scalar>(1);
401 const scalar minCos = Foam::cos(degToRad(featureAngle));
403 Info<< "Feature:" << featureAngle << endl
404 << "minCos :" << minCos << endl
408 const bool splitAllFaces = args.optionFound("splitAllFaces");
411 Info<< "Splitting all internal faces to create multiple faces"
412 << " between two cells." << nl
416 const bool overwrite = args.optionFound("overwrite");
417 const bool doNotPreserveFaceZones = args.optionFound
419 "doNotPreserveFaceZones"
421 const bool concaveMultiCells = args.optionFound("concaveMultiCells");
422 if (concaveMultiCells)
424 Info<< "Generating multiple cells for points on concave feature edges."
429 // Face(centre)s that need inclusion in the dual mesh
430 labelList featureFaces;
432 labelList featureEdges;
433 // Points (that become a single cell) that need inclusion in the dual mesh
434 labelList singleCellFeaturePoints;
435 // Points (that become a multiple cells) ,,
436 labelList multiCellFeaturePoints;
438 // Sample implementation of feature detection.
445 doNotPreserveFaceZones,
449 singleCellFeaturePoints,
450 multiCellFeaturePoints
453 // If we want to split all polyMesh faces into one dualface per cell
454 // we are passing through we also need a point
455 // at the polyMesh facecentre and edgemid of the faces we want to
459 featureEdges = identity(mesh.nEdges());
460 featureFaces = identity(mesh.nFaces());
463 // Write obj files for debugging
469 singleCellFeaturePoints,
470 multiCellFeaturePoints
475 // Read objects in time directory
476 IOobjectList objects(mesh, runTime.timeName());
479 PtrList<volScalarField> vsFlds;
480 ReadFields(mesh, objects, vsFlds);
482 PtrList<volVectorField> vvFlds;
483 ReadFields(mesh, objects, vvFlds);
485 PtrList<volSphericalTensorField> vstFlds;
486 ReadFields(mesh, objects, vstFlds);
488 PtrList<volSymmTensorField> vsymtFlds;
489 ReadFields(mesh, objects, vsymtFlds);
491 PtrList<volTensorField> vtFlds;
492 ReadFields(mesh, objects, vtFlds);
494 // Read surface fields.
495 PtrList<surfaceScalarField> ssFlds;
496 ReadFields(mesh, objects, ssFlds);
498 PtrList<surfaceVectorField> svFlds;
499 ReadFields(mesh, objects, svFlds);
501 PtrList<surfaceSphericalTensorField> sstFlds;
502 ReadFields(mesh, objects, sstFlds);
504 PtrList<surfaceSymmTensorField> ssymtFlds;
505 ReadFields(mesh, objects, ssymtFlds);
507 PtrList<surfaceTensorField> stFlds;
508 ReadFields(mesh, objects, stFlds);
511 // Topo change container
512 polyTopoChange meshMod(mesh.boundaryMesh().size());
514 // Mesh dualiser engine
515 meshDualiser dualMaker(mesh);
517 // Insert all commands into polyTopoChange to create dual of mesh. This does
518 // all the hard work.
519 dualMaker.setRefinement
524 singleCellFeaturePoints,
525 multiCellFeaturePoints,
529 // Create mesh, return map from old to new mesh.
530 autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
533 mesh.updateMesh(map);
535 // Optionally inflate mesh
536 if (map().hasMotionPoints())
538 mesh.movePoints(map().preMotionPoints());
547 mesh.setInstance(oldInstance);
550 Info<< "Writing dual mesh to " << runTime.timeName() << endl;
554 Info<< "End\n" << endl;
560 // ************************************************************************* //