mergeBaffles.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-2023 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) and merge them into internal
29  faces.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "Time.H"
35 #include "syncTools.H"
36 #include "faceSet.H"
37 #include "pointSet.H"
38 #include "meshTools.H"
39 #include "polyTopoChange.H"
40 #include "polyRemoveFace.H"
41 #include "polyModifyFace.H"
42 #include "indirectPrimitivePatch.H"
43 #include "processorPolyPatch.H"
44 #include "localPointRegion.H"
45 #include "duplicatePoints.H"
46 #include "ReadFields.H"
47 #include "volFields.H"
48 #include "surfaceFields.H"
49 #include "pointFields.H"
50 
51 using namespace Foam;
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 void mergeDuplicateBoundaryFaces
56 (
57  const polyMesh& mesh,
58  polyTopoChange& meshMod
59 )
60 {
61  // Get all duplicate face labels in the boundary
63  (
64  mesh,
65  identityMap(mesh.nFaces() - mesh.nInternalFaces())
66  + mesh.nInternalFaces()
67  );
68 
69  // Check that none are on processor patches
70  const polyBoundaryMesh& patches = mesh.boundaryMesh();
71  forAll(duplicates, bFacei)
72  {
73  if (duplicates[bFacei] != -1)
74  {
75  label facei = mesh.nInternalFaces() + bFacei;
76  label patchi = patches.whichPatch(facei);
77 
78  if (isA<processorPolyPatch>(patches[patchi]))
79  {
81  << "Duplicate face " << facei
82  << " is on a processorPolyPatch."
83  << "This is not allowed." << nl
84  << "Face:" << facei
85  << " is on patch:" << patches[patchi].name()
86  << abort(FatalError);
87  }
88  }
89  }
90 
91  const faceList& faces = mesh.faces();
92  const labelList& faceOwner = mesh.faceOwner();
93  const meshFaceZones& faceZones = mesh.faceZones();
94 
95  forAll(duplicates, bFacei)
96  {
97  label otherFacei = duplicates[bFacei];
98 
99  if (otherFacei != -1 && otherFacei > bFacei)
100  {
101  // Two duplicate faces. Merge.
102 
103  label face0 = mesh.nInternalFaces() + bFacei;
104  label face1 = mesh.nInternalFaces() + otherFacei;
105 
106  label own0 = faceOwner[face0];
107  label own1 = faceOwner[face1];
108 
109  if (own0 < own1)
110  {
111  // Use face0 as the new internal face.
112  label zoneID = faceZones.whichZone(face0);
113  bool zoneFlip = false;
114 
115  if (zoneID >= 0)
116  {
117  const faceZone& fZone = faceZones[zoneID];
118  zoneFlip = fZone.flipMap()[fZone.whichFace(face0)];
119  }
120 
121  meshMod.setAction(polyRemoveFace(face1));
122  meshMod.setAction
123  (
125  (
126  faces[face0], // modified face
127  face0, // label of face being modified
128  own0, // owner
129  own1, // neighbour
130  false, // face flip
131  -1, // patch for face
132  false, // remove from zone
133  zoneID, // zone for face
134  zoneFlip // face flip in zone
135  )
136  );
137  }
138  else
139  {
140  // Use face1 as the new internal face.
141  label zoneID = faceZones.whichZone(face1);
142  bool zoneFlip = false;
143 
144  if (zoneID >= 0)
145  {
146  const faceZone& fZone = faceZones[zoneID];
147  zoneFlip = fZone.flipMap()[fZone.whichFace(face1)];
148  }
149 
150  meshMod.setAction(polyRemoveFace(face0));
151  meshMod.setAction
152  (
154  (
155  faces[face1], // modified face
156  face1, // label of face being modified
157  own1, // owner
158  own0, // neighbour
159  false, // face flip
160  -1, // patch for face
161  false, // remove from zone
162  zoneID, // zone for face
163  zoneFlip // face flip in zone
164  )
165  );
166  }
167  }
168  }
169 }
170 
171 
172 int main(int argc, char *argv[])
173 {
175  (
176  "Detect faces that share points (baffles)\n"
177  "and merge them into internal faces."
178  );
179 
180  #include "addOverwriteOption.H"
181  #include "addRegionOption.H"
183  (
184  "fields",
185  "update fields"
186  );
187 
188  #include "setRootCase.H"
190  #include "createNamedMesh.H"
191 
192  const bool overwrite = args.optionFound("overwrite");
193  const bool fields = args.optionFound("fields");
194 
195  const word oldInstance = mesh.pointsInstance();
196 
197  // Read objects in time directory
198  IOobjectList objects(mesh, runTime.name());
199 
200  if (fields) Info<< "Reading geometric fields" << nl << endl;
201 
202  #include "readVolFields.H"
203  #include "readSurfaceFields.H"
204  #include "readPointFields.H"
205 
206  // Mesh change engine
207  polyTopoChange meshMod(mesh);
208 
209  // Merge duplicate boundary faces into internal faces.
210  mergeDuplicateBoundaryFaces(mesh, meshMod);
211 
212  if (!overwrite)
213  {
214  runTime++;
215  }
216 
217  // Change the mesh. No inflation.
218  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh, false);
219 
220  // Update fields
221  mesh.topoChange(map);
222 
223  // Move mesh (since morphing does not do this)
224  if (map().hasMotionPoints())
225  {
226  mesh.setPoints(map().preMotionPoints());
227  }
228 
229  if (overwrite)
230  {
231  mesh.setInstance(oldInstance);
232  }
233 
234  Info<< "Writing mesh to time " << runTime.name() << endl;
235  mesh.write();
236 
237  Info<< "End\n" << endl;
238 
239  return 0;
240 }
241 
242 
243 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: MeshZones.C:221
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:68
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:307
static labelList findDuplicateFaces(const primitiveMesh &, const labelList &)
Helper routine to find baffles (two boundary faces using the.
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:447
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
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:1019
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1448
Class describing modification of a face.
Class containing data for face removal.
Direct mesh changes based on v1.3 polyTopoChange syntax.
autoPtr< polyTopoChangeMap > 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.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
label nInternalFaces() const
label nFaces() const
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
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:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static const char nl
Definition: Ostream.H:260
objects
Foam::argList args(argc, argv)
Foam::surfaceFields.