combinePatchFaces.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 | Copyright (C) 2011-2015 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  combinePatchFaces
26 
27 Description
28  Checks for multiple patch faces on same cell and combines them.
29  Multiple patch faces can result from e.g. removal of refined
30  neighbouring cells, leaving 4 exposed faces with same owner.
31 
32  Rules for merging:
33  - only boundary faces (since multiple internal faces between two cells
34  not allowed anyway)
35  - faces have to have same owner
36  - faces have to be connected via edge which are not features (so angle
37  between them < feature angle)
38  - outside of faces has to be single loop
39  - outside of face should not be (or just slightly) concave (so angle
40  between consecutive edges < concaveangle
41 
42  E.g. to allow all faces on same patch to be merged:
43 
44  combinePatchFaces 180 -concaveAngle 90
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "PstreamReduceOps.H"
49 #include "argList.H"
50 #include "Time.H"
51 #include "polyTopoChange.H"
52 #include "polyModifyFace.H"
53 #include "polyAddFace.H"
54 #include "combineFaces.H"
55 #include "removePoints.H"
56 #include "polyMesh.H"
57 #include "mapPolyMesh.H"
58 #include "unitConversion.H"
59 #include "motionSmoother.H"
60 
61 using namespace Foam;
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 // Merge faces on the same patch (usually from exposing refinement)
66 // Can undo merges if these cause problems.
68 (
69  const scalar minCos,
70  const scalar concaveSin,
71  const autoPtr<IOdictionary>& qualDictPtr,
72  const Time& runTime,
73  polyMesh& mesh
74 )
75 {
76  // Patch face merging engine
77  combineFaces faceCombiner(mesh);
78 
79  // Get all sets of faces that can be merged
80  labelListList allFaceSets(faceCombiner.getMergeSets(minCos, concaveSin));
81 
82  label nFaceSets = returnReduce(allFaceSets.size(), sumOp<label>());
83 
84  Info<< "Merging " << nFaceSets << " sets of faces." << endl;
85 
86  if (nFaceSets > 0)
87  {
88  // Store the faces of the face sets
89  List<faceList> allFaceSetsFaces(allFaceSets.size());
90  forAll(allFaceSets, setI)
91  {
92  allFaceSetsFaces[setI] = UIndirectList<face>
93  (
94  mesh.faces(),
95  allFaceSets[setI]
96  );
97  }
98 
100  {
101  // Topology changes container
102  polyTopoChange meshMod(mesh);
103 
104  // Merge all faces of a set into the first face of the set.
105  faceCombiner.setRefinement(allFaceSets, meshMod);
106 
107  // Change the mesh (no inflation)
108  map = meshMod.changeMesh(mesh, false, true);
109 
110  // Update fields
111  mesh.updateMesh(map);
112 
113  // Move mesh (since morphing does not do this)
114  if (map().hasMotionPoints())
115  {
116  mesh.movePoints(map().preMotionPoints());
117  }
118  else
119  {
120  // Delete mesh volumes. No other way to do this?
121  mesh.clearOut();
122  }
123  }
124 
125 
126  // Check for errors and undo
127  // ~~~~~~~~~~~~~~~~~~~~~~~~~
128 
129  // Faces in error.
130  labelHashSet errorFaces;
131 
132  if (qualDictPtr.valid())
133  {
134  motionSmoother::checkMesh(false, mesh, qualDictPtr(), errorFaces);
135  }
136  else
137  {
138  mesh.checkFacePyramids(false, -SMALL, &errorFaces);
139  }
140 
141  // Sets where the master is in error
142  labelHashSet errorSets;
143 
144  forAll(allFaceSets, setI)
145  {
146  label newMasterI = map().reverseFaceMap()[allFaceSets[setI][0]];
147 
148  if (errorFaces.found(newMasterI))
149  {
150  errorSets.insert(setI);
151  }
152  }
153  label nErrorSets = returnReduce(errorSets.size(), sumOp<label>());
154 
155  Info<< "Detected " << nErrorSets
156  << " error faces on boundaries that have been merged."
157  << " These will be restored to their original faces."
158  << endl;
159 
160  if (nErrorSets > 0)
161  {
162  // Renumber stored faces to new vertex numbering.
163  forAllConstIter(labelHashSet, errorSets, iter)
164  {
165  label setI = iter.key();
166 
167  faceList& setFaceVerts = allFaceSetsFaces[setI];
168 
169  forAll(setFaceVerts, i)
170  {
171  inplaceRenumber(map().reversePointMap(), setFaceVerts[i]);
172 
173  // Debug: check that all points are still there.
174  forAll(setFaceVerts[i], j)
175  {
176  label newVertI = setFaceVerts[i][j];
177 
178  if (newVertI < 0)
179  {
181  << "In set:" << setI << " old face labels:"
182  << allFaceSets[setI] << " new face vertices:"
183  << setFaceVerts[i] << " are unmapped vertices!"
184  << abort(FatalError);
185  }
186  }
187  }
188  }
189 
190 
191  // Topology changes container
192  polyTopoChange meshMod(mesh);
193 
194 
195  // Restore faces
196  forAllConstIter(labelHashSet, errorSets, iter)
197  {
198  label setI = iter.key();
199 
200  const labelList& setFaces = allFaceSets[setI];
201  const faceList& setFaceVerts = allFaceSetsFaces[setI];
202 
203  label newMasterI = map().reverseFaceMap()[setFaces[0]];
204 
205  // Restore. Get face properties.
206 
207  label own = mesh.faceOwner()[newMasterI];
208  label zoneID = mesh.faceZones().whichZone(newMasterI);
209  bool zoneFlip = false;
210  if (zoneID >= 0)
211  {
212  const faceZone& fZone = mesh.faceZones()[zoneID];
213  zoneFlip = fZone.flipMap()[fZone.whichFace(newMasterI)];
214  }
215  label patchID = mesh.boundaryMesh().whichPatch(newMasterI);
216 
217  Pout<< "Restoring new master face " << newMasterI
218  << " to vertices " << setFaceVerts[0] << endl;
219 
220  // Modify the master face.
221  meshMod.setAction
222  (
224  (
225  setFaceVerts[0], // original face
226  newMasterI, // label of face
227  own, // owner
228  -1, // neighbour
229  false, // face flip
230  patchID, // patch for face
231  false, // remove from zone
232  zoneID, // zone for face
233  zoneFlip // face flip in zone
234  )
235  );
236 
237 
238  // Add the previously removed faces
239  for (label i = 1; i < setFaces.size(); i++)
240  {
241  Pout<< "Restoring removed face " << setFaces[i]
242  << " with vertices " << setFaceVerts[i] << endl;
243 
244  meshMod.setAction
245  (
247  (
248  setFaceVerts[i], // vertices
249  own, // owner,
250  -1, // neighbour,
251  -1, // masterPointID,
252  -1, // masterEdgeID,
253  newMasterI, // masterFaceID,
254  false, // flipFaceFlux,
255  patchID, // patchID,
256  zoneID, // zoneID,
257  zoneFlip // zoneFlip
258  )
259  );
260  }
261  }
262 
263  // Change the mesh (no inflation)
264  map = meshMod.changeMesh(mesh, false, true);
265 
266  // Update fields
267  mesh.updateMesh(map);
268 
269  // Move mesh (since morphing does not do this)
270  if (map().hasMotionPoints())
271  {
272  mesh.movePoints(map().preMotionPoints());
273  }
274  else
275  {
276  // Delete mesh volumes. No other way to do this?
277  mesh.clearOut();
278  }
279  }
280  }
281  else
282  {
283  Info<< "No faces merged ..." << endl;
284  }
285 
286  return nFaceSets;
287 }
288 
289 
290 // Remove points not used by any face or points used by only two faces where
291 // the edges are in line
292 label mergeEdges(const scalar minCos, polyMesh& mesh)
293 {
294  Info<< "Merging all points on surface that" << nl
295  << "- are used by only two boundary faces and" << nl
296  << "- make an angle with a cosine of more than " << minCos
297  << "." << nl << endl;
298 
299  // Point removal analysis engine
300  removePoints pointRemover(mesh);
301 
302  // Count usage of points
303  boolList pointCanBeDeleted;
304  label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
305 
306  if (nRemove > 0)
307  {
308  Info<< "Removing " << nRemove
309  << " straight edge points ..." << endl;
310 
311  // Topology changes container
312  polyTopoChange meshMod(mesh);
313 
314  pointRemover.setRefinement(pointCanBeDeleted, meshMod);
315 
316  // Change the mesh (no inflation)
317  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
318 
319  // Update fields
320  mesh.updateMesh(map);
321 
322  // Move mesh (since morphing does not do this)
323  if (map().hasMotionPoints())
324  {
325  mesh.movePoints(map().preMotionPoints());
326  }
327  else
328  {
329  // Delete mesh volumes. No other way to do this?
330  mesh.clearOut();
331  }
332  }
333  else
334  {
335  Info<< "No straight edges simplified and no points removed ..." << endl;
336  }
337 
338  return nRemove;
339 }
340 
341 
342 
343 int main(int argc, char *argv[])
344 {
345  #include "addOverwriteOption.H"
346 
347  argList::validArgs.append("featureAngle [0..180]");
349  (
350  "concaveAngle",
351  "degrees",
352  "specify concave angle [0..180] (default: 30 degrees)"
353  );
355  (
356  "meshQuality",
357  "read user-defined mesh quality criterions from system/meshQualityDict"
358  );
359 
360  #include "setRootCase.H"
361  #include "createTime.H"
362  runTime.functionObjects().off();
363  #include "createPolyMesh.H"
364  const word oldInstance = mesh.pointsInstance();
365 
366  const scalar featureAngle = args.argRead<scalar>(1);
367  const scalar minCos = Foam::cos(degToRad(featureAngle));
368 
369  // Sin of angle between two consecutive edges on a face.
370  // If sin(angle) larger than this the face will be considered concave.
371  scalar concaveAngle = args.optionLookupOrDefault("concaveAngle", 30.0);
372  scalar concaveSin = Foam::sin(degToRad(concaveAngle));
373 
374  const bool overwrite = args.optionFound("overwrite");
375  const bool meshQuality = args.optionFound("meshQuality");
376 
377  Info<< "Merging all faces of a cell" << nl
378  << " - which are on the same patch" << nl
379  << " - which make an angle < " << featureAngle << " degrees"
380  << nl
381  << " (cos:" << minCos << ')' << nl
382  << " - even when resulting face becomes concave by more than "
383  << concaveAngle << " degrees" << nl
384  << " (sin:" << concaveSin << ')' << nl
385  << endl;
386 
387  autoPtr<IOdictionary> qualDict;
388  if (meshQuality)
389  {
390  Info<< "Enabling user-defined geometry checks." << nl << endl;
391 
392  qualDict.reset
393  (
394  new IOdictionary
395  (
396  IOobject
397  (
398  "meshQualityDict",
399  mesh.time().system(),
400  mesh,
403  )
404  )
405  );
406  }
407 
408 
409  if (!overwrite)
410  {
411  runTime++;
412  }
413 
414 
415 
416  // Merge faces on same patch
417  label nChanged = mergePatchFaces
418  (
419  minCos,
420  concaveSin,
421  qualDict,
422  runTime,
423  mesh
424  );
425 
426  // Merge points on straight edges and remove unused points
427  if (qualDict.valid())
428  {
429  Info<< "Merging all 'loose' points on surface edges, "
430  << "regardless of the angle they make." << endl;
431 
432  // Surface bnound to be used to extrude. Merge all loose points.
433  nChanged += mergeEdges(-1, mesh);
434  }
435  else
436  {
437  nChanged += mergeEdges(minCos, mesh);
438  }
439 
440  if (nChanged > 0)
441  {
442  if (overwrite)
443  {
444  mesh.setInstance(oldInstance);
445  }
446 
447  Info<< "Writing morphed mesh to time " << runTime.timeName() << endl;
448 
449  mesh.write();
450  }
451  else
452  {
453  Info<< "Mesh unchanged." << endl;
454  }
455 
456  Info<< "\nEnd\n" << endl;
457 
458  return 0;
459 }
460 
461 
462 // ************************************************************************* //
Foam::argList::validArgs
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:143
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::combineFaces
Combines boundary faces into single face. The faces get the patch of the first face ('the master')
Definition: combineFaces.H:55
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::argList::addOption
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:108
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::removePoints
Removes selected points from mesh and updates faces using these points.
Definition: removePoints.H:58
polyAddFace.H
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
addOverwriteOption.H
mapPolyMesh.H
motionSmoother.H
polyTopoChange.H
combineFaces.H
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::Time::functionObjects
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:435
unitConversion.H
Unit conversion functions.
Foam::polyModifyFace
Class describing modification of a face.
Definition: polyModifyFace.H:48
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:726
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::help::mergePatchFaces
faceList mergePatchFaces(const List< DynList< label > > &pfcs, const pointField &polyPoints)
Definition: helperFunctionsGeometryQueriesI.H:232
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
Definition: ListOpsTemplates.C:57
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
removePoints.H
Foam::label
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
Foam::TimePaths::system
const word & system() const
Return system name.
Definition: TimePaths.H:120
Foam::functionObjectList::off
virtual void off()
Switch the function objects off.
Definition: functionObjectList.C:178
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:802
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:703
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::argList::optionLookupOrDefault
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:237
Foam::ZoneMesh::whichZone
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:226
Foam::FatalError
error FatalError
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:314
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
PstreamReduceOps.H
polyModifyFace.H
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::motionSmootherAlgo::checkMesh
static bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, labelHashSet &wrongFaces)
Check mesh with mesh settings in dict. Collects incorrect faces.
Definition: motionSmootherAlgoCheck.C:415
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::autoPtr::valid
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set).
Definition: autoPtrI.H:83
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::polyAddFace
A face addition data class. A face can be inflated either from a point or from another face and can e...
Definition: polyAddFace.H:51
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::argList::argRead
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:177
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
args
Foam::argList args(argc, argv)
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Foam::primitiveMesh::checkFacePyramids
bool checkFacePyramids(const bool report=false, const scalar minPyrVol=-SMALL, labelHashSet *setPtr=NULL) const
Check face pyramid volume.
Definition: primitiveMeshCheck.C:585
createPolyMesh.H
Foam::fvMesh::clearOut
void clearOut()
Clear all geometry and addressing.
Definition: fvMesh.C:231
Foam::degToRad
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
Definition: unitConversion.H:45
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:256