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-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  mergeBaffles
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 "indirectPrimitivePatch.H"
41 #include "processorPolyPatch.H"
42 #include "localPointRegion.H"
43 #include "duplicatePoints.H"
44 #include "ReadFields.H"
45 #include "volFields.H"
46 #include "surfaceFields.H"
47 #include "pointFields.H"
48 
49 using namespace Foam;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 void mergeDuplicateBoundaryFaces
54 (
55  const polyMesh& mesh,
56  polyTopoChange& meshMod
57 )
58 {
59  // Get all duplicate face labels in the boundary
61  (
62  mesh,
63  identityMap(mesh.nFaces() - mesh.nInternalFaces())
64  + mesh.nInternalFaces()
65  );
66 
67  // Check that none are on processor patches
68  const polyBoundaryMesh& patches = mesh.boundaryMesh();
69  forAll(duplicates, bFacei)
70  {
71  if (duplicates[bFacei] != -1)
72  {
73  label facei = mesh.nInternalFaces() + bFacei;
74  label patchi = patches.whichPatch(facei);
75 
76  if (isA<processorPolyPatch>(patches[patchi]))
77  {
79  << "Duplicate face " << facei
80  << " is on a processorPolyPatch."
81  << "This is not allowed." << nl
82  << "Face:" << facei
83  << " is on patch:" << patches[patchi].name()
84  << abort(FatalError);
85  }
86  }
87  }
88 
89  const faceList& faces = mesh.faces();
90  const labelList& faceOwner = mesh.faceOwner();
91 
92  forAll(duplicates, bFacei)
93  {
94  label otherFacei = duplicates[bFacei];
95 
96  if (otherFacei != -1 && otherFacei > bFacei)
97  {
98  // Two duplicate faces. Merge.
99 
100  label face0 = mesh.nInternalFaces() + bFacei;
101  label face1 = mesh.nInternalFaces() + otherFacei;
102 
103  label own0 = faceOwner[face0];
104  label own1 = faceOwner[face1];
105 
106  if (own0 < own1)
107  {
108  // Use face0 as the new internal face.
109 
110  meshMod.removeFace(face1, -1);
111  meshMod.modifyFace
112  (
113  faces[face0], // modified face
114  face0, // label of face being modified
115  own0, // owner
116  own1, // neighbour
117  false, // face flip
118  -1 // patch for face
119  );
120  }
121  else
122  {
123  // Use face1 as the new internal face.
124 
125  meshMod.removeFace(face0, -1);
126  meshMod.modifyFace
127  (
128  faces[face1], // modified face
129  face1, // label of face being modified
130  own1, // owner
131  own0, // neighbour
132  false, // face flip
133  -1 // patch for face
134  );
135  }
136  }
137  }
138 }
139 
140 
141 int main(int argc, char *argv[])
142 {
144  (
145  "Detect faces that share points (baffles)\n"
146  "and merge them into internal faces."
147  );
148 
149  #include "addOverwriteOption.H"
150  #include "addRegionOption.H"
152  (
153  "fields",
154  "update fields"
155  );
156 
157  #include "setRootCase.H"
160 
161  const bool overwrite = args.optionFound("overwrite");
162  const bool fields = args.optionFound("fields");
163 
164  const word oldInstance = mesh.pointsInstance();
165 
166  // Read objects in time directory
167  IOobjectList objects(mesh, runTime.name());
168 
169  if (fields) Info<< "Reading geometric fields" << nl << endl;
170 
171  #include "readVolFields.H"
172  #include "readSurfaceFields.H"
173  #include "readPointFields.H"
174 
175  // Mesh change engine
176  polyTopoChange meshMod(mesh);
177 
178  // Merge duplicate boundary faces into internal faces.
179  mergeDuplicateBoundaryFaces(mesh, meshMod);
180 
181  if (!overwrite)
182  {
183  runTime++;
184  }
185 
186  // Change the mesh
187  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh);
188 
189  // Update fields
190  mesh.topoChange(map);
191 
192  if (overwrite)
193  {
194  mesh.setInstance(oldInstance);
195  }
196 
197  Info<< "Writing mesh to time " << runTime.name() << endl;
198  mesh.write();
199 
200  Info<< "End\n" << endl;
201 
202  return 0;
203 }
204 
205 
206 // ************************************************************************* //
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
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
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
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.
void removeFace(const label, const label)
Remove face / merge faces.
autoPtr< polyTopoChangeMap > changeMesh(polyMesh &mesh, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
void modifyFace(const face &f, const label facei, const label own, const label nei, const bool flipFaceFlux, const label patchID)
Modify vertices or cell of face.
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:334
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:228
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:257
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:266
objects
Foam::argList args(argc, argv)
Foam::surfaceFields.