combinePatchFaces.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  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 "combineFaces.H"
53 #include "removePoints.H"
54 #include "meshCheck.H"
55 #include "polyTopoChangeMap.H"
56 
57 using namespace Foam;
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 // Merge faces on the same patch (usually from exposing refinement)
62 // Can undo merges if these cause problems.
63 label mergePatchFaces
64 (
65  const scalar minCos,
66  const scalar concaveSin,
67  const autoPtr<IOdictionary>& qualDictPtr,
68  const Time& runTime,
69  polyMesh& mesh
70 )
71 {
72  // Patch face merging engine
73  combineFaces faceCombiner(mesh);
74 
75  // Get all sets of faces that can be merged
76  labelListList allFaceSets(faceCombiner.getMergeSets(minCos, concaveSin));
77 
78  label nFaceSets = returnReduce(allFaceSets.size(), sumOp<label>());
79 
80  Info<< "Merging " << nFaceSets << " sets of faces." << endl;
81 
82  if (nFaceSets > 0)
83  {
84  // Store the faces of the face sets
85  List<faceList> allFaceSetsFaces(allFaceSets.size());
86  forAll(allFaceSets, setI)
87  {
88  allFaceSetsFaces[setI] = UIndirectList<face>
89  (
90  mesh.faces(),
91  allFaceSets[setI]
92  );
93  }
94 
96  {
97  // Topology changes container
98  polyTopoChange meshMod(mesh);
99 
100  // Merge all faces of a set into the first face of the set.
101  faceCombiner.setRefinement(allFaceSets, meshMod);
102 
103  // Change the mesh
104  map = meshMod.changeMesh(mesh, true);
105 
106  // Update fields
107  mesh.topoChange(map);
108  }
109 
110 
111  // Check for errors and undo
112  // ~~~~~~~~~~~~~~~~~~~~~~~~~
113 
114  // Faces in error.
115  labelHashSet errorFaces;
116 
117  if (qualDictPtr.valid())
118  {
119  meshCheck::checkMesh(false, mesh, qualDictPtr(), errorFaces);
120  }
121  else
122  {
123  meshCheck::checkFacePyramids(mesh, false, -small, &errorFaces);
124  }
125 
126  // Sets where the master is in error
127  labelHashSet errorSets;
128 
129  forAll(allFaceSets, setI)
130  {
131  label newMasterI = map().reverseFaceMap()[allFaceSets[setI][0]];
132 
133  if (errorFaces.found(newMasterI))
134  {
135  errorSets.insert(setI);
136  }
137  }
138  label nErrorSets = returnReduce(errorSets.size(), sumOp<label>());
139 
140  Info<< "Detected " << nErrorSets
141  << " error faces on boundaries that have been merged."
142  << " These will be restored to their original faces."
143  << endl;
144 
145  if (nErrorSets > 0)
146  {
147  // Renumber stored faces to new vertex numbering.
148  forAllConstIter(labelHashSet, errorSets, iter)
149  {
150  label setI = iter.key();
151 
152  faceList& setFaceVerts = allFaceSetsFaces[setI];
153 
154  forAll(setFaceVerts, i)
155  {
156  inplaceRenumber(map().reversePointMap(), setFaceVerts[i]);
157 
158  // Debug: check that all points are still there.
159  forAll(setFaceVerts[i], j)
160  {
161  label newVertI = setFaceVerts[i][j];
162 
163  if (newVertI < 0)
164  {
166  << "In set:" << setI << " old face labels:"
167  << allFaceSets[setI] << " new face vertices:"
168  << setFaceVerts[i] << " are unmapped vertices!"
169  << abort(FatalError);
170  }
171  }
172  }
173  }
174 
175 
176  // Topology changes container
177  polyTopoChange meshMod(mesh);
178 
179 
180  // Restore faces
181  forAllConstIter(labelHashSet, errorSets, iter)
182  {
183  label setI = iter.key();
184 
185  const labelList& setFaces = allFaceSets[setI];
186  const faceList& setFaceVerts = allFaceSetsFaces[setI];
187 
188  label newMasterI = map().reverseFaceMap()[setFaces[0]];
189 
190  // Restore. Get face properties.
191 
192  label own = mesh.faceOwner()[newMasterI];
193  label patchID = mesh.boundaryMesh().whichPatch(newMasterI);
194 
195  Pout<< "Restoring new master face " << newMasterI
196  << " to vertices " << setFaceVerts[0] << endl;
197 
198  // Modify the master face.
199  meshMod.modifyFace
200  (
201  setFaceVerts[0], // original face
202  newMasterI, // label of face
203  own, // owner
204  -1, // neighbour
205  false, // face flip
206  patchID // patch for face
207  );
208 
209  // Add the previously removed faces
210  for (label i = 1; i < setFaces.size(); i++)
211  {
212  Pout<< "Restoring removed face " << setFaces[i]
213  << " with vertices " << setFaceVerts[i] << endl;
214 
215  meshMod.addFace
216  (
217  setFaceVerts[i], // vertices
218  own, // owner,
219  -1, // neighbour,
220  newMasterI, // masterFaceID,
221  false, // flipFaceFlux,
222  patchID // patchID,
223  );
224  }
225  }
226 
227  // Change the mesh
228  map = meshMod.changeMesh(mesh, true);
229 
230  // Update fields
231  mesh.topoChange(map);
232  }
233  }
234  else
235  {
236  Info<< "No faces merged ..." << endl;
237  }
238 
239  return nFaceSets;
240 }
241 
242 
243 // Remove points not used by any face or points used by only two faces where
244 // the edges are in line
245 label mergeEdges(const scalar minCos, polyMesh& mesh)
246 {
247  Info<< "Merging all points on surface that" << nl
248  << "- are used by only two boundary faces and" << nl
249  << "- make an angle with a cosine of more than " << minCos
250  << "." << nl << endl;
251 
252  // Point removal analysis engine
253  removePoints pointRemover(mesh);
254 
255  // Count usage of points
256  boolList pointCanBeDeleted;
257  label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
258 
259  if (nRemove > 0)
260  {
261  Info<< "Removing " << nRemove
262  << " straight edge points ..." << endl;
263 
264  // Topology changes container
265  polyTopoChange meshMod(mesh);
266 
267  pointRemover.setRefinement(pointCanBeDeleted, meshMod);
268 
269  // Change the mesh
270  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
271 
272  // Update fields
273  mesh.topoChange(map);
274  }
275  else
276  {
277  Info<< "No straight edges simplified and no points removed ..." << endl;
278  }
279 
280  return nRemove;
281 }
282 
283 
284 
285 int main(int argc, char *argv[])
286 {
287  #include "addOverwriteOption.H"
288 
289  argList::validArgs.append("featureAngle [0..180]");
291  (
292  "concaveAngle",
293  "degrees",
294  "specify concave angle [0..180] (default: 30 degrees)"
295  );
297  (
298  "meshQuality",
299  "read user-defined mesh quality criterions from system/meshQualityDict"
300  );
301 
302  #include "setRootCase.H"
304  #include "createPolyMesh.H"
305  const word oldInstance = mesh.pointsInstance();
306 
307  const scalar featureAngle = degToRad(args.argRead<scalar>(1));
308  const scalar minCos = Foam::cos(featureAngle);
309 
310  // Sin of angle between two consecutive edges on a face.
311  // If sin(angle) larger than this the face will be considered concave.
312  const scalar concaveAngle =
313  degToRad(args.optionLookupOrDefault("concaveAngle", 30.0));
314  const scalar concaveSin = Foam::sin(concaveAngle);
315 
316  const bool overwrite = args.optionFound("overwrite");
317  const bool meshQuality = args.optionFound("meshQuality");
318 
319  Info<< "Merging all faces of a cell" << nl
320  << " - which are on the same patch" << nl
321  << " - which make an angle < " << radToDeg(featureAngle)
322  << " degrees" << nl
323  << " (cos:" << minCos << ')' << nl
324  << " - even when resulting face becomes concave by more than "
325  << radToDeg(concaveAngle) << " degrees" << nl
326  << " (sin:" << concaveSin << ')' << nl
327  << endl;
328 
329  autoPtr<IOdictionary> qualDict;
330  if (meshQuality)
331  {
332  Info<< "Enabling user-defined geometry checks." << nl << endl;
333 
334  qualDict.reset
335  (
336  new IOdictionary
337  (
338  IOobject
339  (
340  "meshQualityDict",
341  mesh.time().system(),
342  mesh,
345  )
346  )
347  );
348  }
349 
350 
351  if (!overwrite)
352  {
353  runTime++;
354  }
355 
356 
357 
358  // Merge faces on same patch
359  label nChanged = mergePatchFaces
360  (
361  minCos,
362  concaveSin,
363  qualDict,
364  runTime,
365  mesh
366  );
367 
368  // Merge points on straight edges and remove unused points
369  if (qualDict.valid())
370  {
371  Info<< "Merging all 'loose' points on surface edges, "
372  << "regardless of the angle they make." << endl;
373 
374  // Surface bnound to be used to extrude. Merge all loose points.
375  nChanged += mergeEdges(-1, mesh);
376  }
377  else
378  {
379  nChanged += mergeEdges(minCos, mesh);
380  }
381 
382  if (nChanged > 0)
383  {
384  if (overwrite)
385  {
386  mesh.setInstance(oldInstance);
387  }
388 
389  Info<< "Writing morphed mesh to time " << runTime.name() << endl;
390 
391  mesh.write();
392  }
393  else
394  {
395  Info<< "Mesh unchanged." << endl;
396  }
397 
398  Info<< "\nEnd\n" << endl;
399 
400  return 0;
401 }
402 
403 
404 // ************************************************************************* //
Inter-processor communication reduction functions.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const word & system()
Return system name.
Definition: TimePaths.H:112
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A List with indirect addressing.
Definition: UIndirectList.H:60
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 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
T argRead(const label index) const
Read a value from the argument at index.
Definition: argListI.H:183
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Combines boundary faces into single face. The faces get the patch of the first face ('the master')
Definition: combineFaces.H:56
const word & name() const
Return const reference to name.
const Time & time() const
Return time.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
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 polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual bool write(const bool write=true) const
Write using setting from DB.
Removes selected points from mesh and updates faces using these points.
Definition: removePoints.H:59
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
int main(int argc, char *argv[])
Definition: financialFoam.C:44
Functions for checking mesh topology and geometry.
bool checkMesh(const bool report, const polyMesh &mesh, const dictionary &dict, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet &wrongFaces)
Check (subset of mesh including baffles) with mesh settings in dict.
Definition: checkMesh.C:33
bool checkFacePyramids(const bool report, const scalar minPyrVol, const polyMesh &, const vectorField &cellCentres, const pointField &p, const labelList &checkFaces, const List< labelPair > &baffles, labelHashSet *)
Check face pyramid volumes.
Namespace for OpenFOAM.
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
Foam::argList args(argc, argv)