mergeOrSplitBaffles.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-2016 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  mergeOrSplitBaffles
26 
27 Description
28  Detects faces that share points (baffles). Either merge them or
29  duplicate the points.
30 
31  Notes:
32  - can only handle pairwise boundary faces. So three faces using
33  the same points is not handled (is illegal mesh anyway)
34 
35  - there is no option to only split/merge some baffles.
36 
37  - surfaces consisting of duplicate faces can be topologically split
38  if the points on the interior of the surface cannot walk to all the
39  cells that use them in one go.
40 
41  - Parallel operation (where duplicate face is perpendicular to a coupled
42  boundary) is supported but not really tested.
43  (Note that coupled faces themselves are not seen as duplicate faces)
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "argList.H"
48 #include "Time.H"
49 #include "syncTools.H"
50 #include "faceSet.H"
51 #include "pointSet.H"
52 #include "meshTools.H"
53 #include "polyTopoChange.H"
54 #include "polyRemoveFace.H"
55 #include "polyModifyFace.H"
56 #include "indirectPrimitivePatch.H"
57 #include "processorPolyPatch.H"
58 #include "localPointRegion.H"
59 #include "duplicatePoints.H"
60 #include "ReadFields.H"
61 #include "volFields.H"
62 #include "surfaceFields.H"
63 
64 using namespace Foam;
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 void insertDuplicateMerge
69 (
70  const polyMesh& mesh,
71  const labelList& duplicates,
72  polyTopoChange& meshMod
73 )
74 {
75  const faceList& faces = mesh.faces();
76  const labelList& faceOwner = mesh.faceOwner();
77  const faceZoneMesh& faceZones = mesh.faceZones();
78 
79  forAll(duplicates, bFacei)
80  {
81  label otherFacei = duplicates[bFacei];
82 
83  if (otherFacei != -1 && otherFacei > bFacei)
84  {
85  // Two duplicate faces. Merge.
86 
87  label face0 = mesh.nInternalFaces() + bFacei;
88  label face1 = mesh.nInternalFaces() + otherFacei;
89 
90  label own0 = faceOwner[face0];
91  label own1 = faceOwner[face1];
92 
93  if (own0 < own1)
94  {
95  // Use face0 as the new internal face.
96  label zoneID = faceZones.whichZone(face0);
97  bool zoneFlip = false;
98 
99  if (zoneID >= 0)
100  {
101  const faceZone& fZone = faceZones[zoneID];
102  zoneFlip = fZone.flipMap()[fZone.whichFace(face0)];
103  }
104 
105  meshMod.setAction(polyRemoveFace(face1));
106  meshMod.setAction
107  (
109  (
110  faces[face0], // modified face
111  face0, // label of face being modified
112  own0, // owner
113  own1, // neighbour
114  false, // face flip
115  -1, // patch for face
116  false, // remove from zone
117  zoneID, // zone for face
118  zoneFlip // face flip in zone
119  )
120  );
121  }
122  else
123  {
124  // Use face1 as the new internal face.
125  label zoneID = faceZones.whichZone(face1);
126  bool zoneFlip = false;
127 
128  if (zoneID >= 0)
129  {
130  const faceZone& fZone = faceZones[zoneID];
131  zoneFlip = fZone.flipMap()[fZone.whichFace(face1)];
132  }
133 
134  meshMod.setAction(polyRemoveFace(face0));
135  meshMod.setAction
136  (
138  (
139  faces[face1], // modified face
140  face1, // label of face being modified
141  own1, // owner
142  own0, // neighbour
143  false, // face flip
144  -1, // patch for face
145  false, // remove from zone
146  zoneID, // zone for face
147  zoneFlip // face flip in zone
148  )
149  );
150  }
151  }
152  }
153 }
154 
155 
156 labelList findBaffles(const polyMesh& mesh, const labelList& boundaryFaces)
157 {
158  // Get all duplicate face labels (in boundaryFaces indices!).
160  (
161  mesh,
162  boundaryFaces
163  );
164 
165 
166  // Check that none are on processor patches
167  const polyBoundaryMesh& patches = mesh.boundaryMesh();
168 
169  forAll(duplicates, bFacei)
170  {
171  if (duplicates[bFacei] != -1)
172  {
173  label facei = mesh.nInternalFaces() + bFacei;
174  label patchi = patches.whichPatch(facei);
175 
176  if (isA<processorPolyPatch>(patches[patchi]))
177  {
179  << "Duplicate face " << facei
180  << " is on a processorPolyPatch."
181  << "This is not allowed." << nl
182  << "Face:" << facei
183  << " is on patch:" << patches[patchi].name()
184  << abort(FatalError);
185  }
186  }
187  }
188 
189 
190  // Write to faceSet for ease of postprocessing.
191  {
192  faceSet duplicateSet
193  (
194  mesh,
195  "duplicateFaces",
196  (mesh.nFaces() - mesh.nInternalFaces())/256
197  );
198 
199  forAll(duplicates, bFacei)
200  {
201  label otherFacei = duplicates[bFacei];
202 
203  if (otherFacei != -1 && otherFacei > bFacei)
204  {
205  duplicateSet.insert(mesh.nInternalFaces() + bFacei);
206  duplicateSet.insert(mesh.nInternalFaces() + otherFacei);
207  }
208  }
209 
210  Pout<< "Writing " << duplicateSet.size()
211  << " duplicate faces to faceSet " << duplicateSet.objectPath()
212  << nl << endl;
213  duplicateSet.write();
214  }
215 
216  return duplicates;
217 }
218 
219 
220 
221 
222 int main(int argc, char *argv[])
223 {
225  (
226  "Detect faces that share points (baffles).\n"
227  "Merge them or duplicate the points."
228  );
229 
230  #include "addOverwriteOption.H"
231  #include "addRegionOption.H"
233  (
234  "detectOnly",
235  "find baffles only, but do not merge or split them"
236  );
238  (
239  "split",
240  "topologically split duplicate surfaces"
241  );
242 
243  #include "setRootCase.H"
244  #include "createTime.H"
245  runTime.functionObjects().off();
246  #include "createNamedMesh.H"
247 
248  const word oldInstance = mesh.pointsInstance();
249 
250  const bool split = args.optionFound("split");
251  const bool overwrite = args.optionFound("overwrite");
252  const bool detectOnly = args.optionFound("detectOnly");
253 
254  // Collect all boundary faces
255  labelList boundaryFaces(mesh.nFaces() - mesh.nInternalFaces());
256 
257  forAll(boundaryFaces, i)
258  {
259  boundaryFaces[i] = i+mesh.nInternalFaces();
260  }
261 
262 
263  if (detectOnly)
264  {
265  findBaffles(mesh, boundaryFaces);
266  return 0;
267  }
268 
269 
270  // Read objects in time directory
271  IOobjectList objects(mesh, runTime.timeName());
272 
273  // Read vol fields.
274 
276  ReadFields(mesh, objects, vsFlds);
277 
279  ReadFields(mesh, objects, vvFlds);
280 
282  ReadFields(mesh, objects, vstFlds);
283 
284  PtrList<volSymmTensorField> vsymtFlds;
285  ReadFields(mesh, objects, vsymtFlds);
286 
288  ReadFields(mesh, objects, vtFlds);
289 
290  // Read surface fields.
291 
293  ReadFields(mesh, objects, ssFlds);
294 
296  ReadFields(mesh, objects, svFlds);
297 
299  ReadFields(mesh, objects, sstFlds);
300 
302  ReadFields(mesh, objects, ssymtFlds);
303 
305  ReadFields(mesh, objects, stFlds);
306 
307 
308  // Mesh change engine
309  polyTopoChange meshMod(mesh);
310 
311 
312  if (split)
313  {
314  Pout<< "Topologically splitting duplicate surfaces"
315  << ", i.e. duplicating points internal to duplicate surfaces."
316  << nl << endl;
317 
318  // Analyse which points need to be duplicated
320 
321  // Point duplication engine
322  duplicatePoints pointDuplicator(mesh);
323 
324  // Insert topo changes
325  pointDuplicator.setRefinement(regionSide, meshMod);
326  }
327  else
328  {
329  Pout<< "Merging duplicate faces."
330  << nl << endl;
331 
332  // Get all duplicate face labels (in boundaryFaces indices!).
333  labelList duplicates(findBaffles(mesh, boundaryFaces));
334 
335  // Merge into internal faces.
336  insertDuplicateMerge(mesh, duplicates, meshMod);
337  }
338 
339  if (!overwrite)
340  {
341  runTime++;
342  }
343 
344  // Change the mesh. No inflation.
345  autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);
346 
347  // Update fields
348  mesh.updateMesh(map);
349 
350  // Move mesh (since morphing does not do this)
351  if (map().hasMotionPoints())
352  {
353  mesh.movePoints(map().preMotionPoints());
354  }
355 
356  if (overwrite)
357  {
358  mesh.setInstance(oldInstance);
359  }
360  Pout<< "Writing mesh to time " << runTime.timeName() << endl;
361  mesh.write();
362 
363  // Dump duplicated points (if any)
364  if (split)
365  {
366  const labelList& pointMap = map().pointMap();
367 
368  labelList nDupPerPoint(map().nOldPoints(), 0);
369 
370  pointSet dupPoints(mesh, "duplicatedPoints", 100);
371 
372  forAll(pointMap, pointi)
373  {
374  label oldPointi = pointMap[pointi];
375 
376  nDupPerPoint[oldPointi]++;
377 
378  if (nDupPerPoint[oldPointi] > 1)
379  {
380  dupPoints.insert(map().reversePointMap()[oldPointi]);
381  dupPoints.insert(pointi);
382  }
383  }
384 
385  Pout<< "Writing " << dupPoints.size()
386  << " duplicated points to pointSet "
387  << dupPoints.objectPath() << nl << endl;
388 
389  dupPoints.write();
390  }
391 
392  Info<< "End\n" << endl;
393 
394  return 0;
395 }
396 
397 
398 // ************************************************************************* //
Foam::surfaceFields.
Class containing data for face removal.
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1048
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A list of face labels.
Definition: faceSet.H:48
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:253
Class describing modification of a face.
A set of point labels.
Definition: pointSet.H:48
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Field reading functions for post-processing utilities.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Takes mesh with &#39;baffles&#39; (= boundary faces sharing points). Determines for selected points on bounda...
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:305
Determines the &#39;side&#39; for every face and connected to a singly-connected (through edges) region of fa...
Definition: regionSide.H:61
Duplicate points.
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
virtual void updateMesh(const mapPolyMesh &mpm)
Update the mesh corresponding to given map.
A class for handling words, derived from string.
Definition: word.H:59
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::polyBoundaryMesh.
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:262
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
label nFaces() const
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:221
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual bool write() const
Write using setting from DB.
static labelList findDuplicateFaces(const primitiveMesh &, const labelList &)
Helper routine to find baffles (two boundary faces using the.
messageStream Info
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::argList args(argc, argv)
label nInternalFaces() const
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.