autoPatch.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  autoPatch
28 
29 Group
30  grpMeshManipulationUtilities
31 
32 Description
33  Divides external faces into patches based on (user supplied) feature
34  angle.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "argList.H"
39 #include "polyMesh.H"
40 #include "Time.H"
41 #include "boundaryMesh.H"
42 #include "repatchPolyTopoChanger.H"
43 #include "unitConversion.H"
44 #include "OFstream.H"
45 #include "ListOps.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 // Get all feature edges.
52 void collectFeatureEdges(const boundaryMesh& bMesh, labelList& markedEdges)
53 {
54  markedEdges.setSize(bMesh.mesh().nEdges());
55 
56  label markedI = 0;
57 
58  forAll(bMesh.featureSegments(), i)
59  {
60  const labelList& segment = bMesh.featureSegments()[i];
61 
62  forAll(segment, j)
63  {
64  label featEdgeI = segment[j];
65 
66  label meshEdgeI = bMesh.featureToEdge()[featEdgeI];
67 
68  markedEdges[markedI++] = meshEdgeI;
69  }
70  }
71  markedEdges.setSize(markedI);
72 }
73 
74 
75 
76 int main(int argc, char *argv[])
77 {
79  (
80  "Divides external faces into patches based on feature angle"
81  );
82 
83  #include "addOverwriteOption.H"
84 
86  argList::noFunctionObjects(); // Never use function objects
87 
88  argList::addArgument("featureAngle", "in degrees [0-180]");
89 
90  #include "setRootCase.H"
91  #include "createTime.H"
92  #include "createPolyMesh.H"
93 
94  const word oldInstance = mesh.pointsInstance();
95 
96  Info<< "Mesh read in = "
98  << " s\n" << endl << endl;
99 
100 
101  const scalar featureAngle = args.get<scalar>(1);
102  const bool overwrite = args.found("overwrite");
103 
104  const scalar minCos = Foam::cos(degToRad(featureAngle));
105 
106  Info<< "Feature:" << featureAngle << endl
107  << "minCos :" << minCos << endl
108  << endl;
109 
110  //
111  // Use boundaryMesh to reuse all the featureEdge stuff in there.
112  //
113 
115  bMesh.read(mesh);
116 
117  // Set feature angle (calculate feature edges)
118  bMesh.setFeatureEdges(minCos);
119 
120  // Collect all feature edges as edge labels
121  labelList markedEdges;
122 
123  collectFeatureEdges(bMesh, markedEdges);
124 
125 
126 
127  // (new) patch ID for every face in mesh.
128  labelList patchIDs(bMesh.mesh().size(), -1);
129 
130  //
131  // Fill patchIDs with values for every face by floodfilling without
132  // crossing feature edge.
133  //
134 
135  // Current patch number.
136  label newPatchi = bMesh.patches().size();
137 
138  label suffix = 0;
139 
140  while (true)
141  {
142  // Find first unset face.
143  label unsetFacei = patchIDs.find(-1);
144 
145  if (unsetFacei == -1)
146  {
147  // All faces have patchID set. Exit.
148  break;
149  }
150 
151  // Found unset face. Create patch for it.
152  word patchName;
153  do
154  {
155  patchName = "auto" + name(suffix++);
156  }
157  while (bMesh.findPatchID(patchName) != -1);
158 
159  bMesh.addPatch(patchName);
160 
161  bMesh.changePatchType(patchName, "patch");
162 
163 
164  // Fill visited with all faces reachable from unsetFacei.
165  boolList visited(bMesh.mesh().size());
166 
167  bMesh.markFaces(markedEdges, unsetFacei, visited);
168 
169 
170  // Assign all visited faces to current patch
171  label nVisited = 0;
172 
173  forAll(visited, facei)
174  {
175  if (visited[facei])
176  {
177  nVisited++;
178 
179  patchIDs[facei] = newPatchi;
180  }
181  }
182 
183  Info<< "Assigned " << nVisited << " faces to patch " << patchName
184  << endl << endl;
185 
186  newPatchi++;
187  }
188 
189 
190 
191  const PtrList<boundaryPatch>& patches = bMesh.patches();
192 
193  // Create new list of patches with old ones first
194  List<polyPatch*> newPatchPtrList(patches.size());
195 
196  newPatchi = 0;
197 
198  // Copy old patches
199  forAll(mesh.boundaryMesh(), patchi)
200  {
201  const polyPatch& patch = mesh.boundaryMesh()[patchi];
202 
203  newPatchPtrList[newPatchi] =
204  patch.clone
205  (
206  mesh.boundaryMesh(),
207  newPatchi,
208  patch.size(),
209  patch.start()
210  ).ptr();
211 
212  newPatchi++;
213  }
214 
215  // Add new ones with empty size.
216  for (label patchi = newPatchi; patchi < patches.size(); patchi++)
217  {
218  const boundaryPatch& bp = patches[patchi];
219 
220  newPatchPtrList[newPatchi] = polyPatch::New
221  (
222  polyPatch::typeName,
223  bp.name(),
224  0,
225  mesh.nFaces(),
226  newPatchi,
228  ).ptr();
229 
230  newPatchi++;
231  }
232 
233  if (!overwrite)
234  {
235  ++runTime;
236  }
237 
238 
239  // Change patches
240  repatchPolyTopoChanger polyMeshRepatcher(mesh);
241  polyMeshRepatcher.changePatches(newPatchPtrList);
242 
243 
244  // Change face ordering
245 
246  // Since bMesh read from mesh there is one to one mapping so we don't
247  // have to do the geometric stuff.
248  const labelList& meshFace = bMesh.meshFace();
249 
250  forAll(patchIDs, facei)
251  {
252  label meshFacei = meshFace[facei];
253 
254  polyMeshRepatcher.changePatchID(meshFacei, patchIDs[facei]);
255  }
256 
257  polyMeshRepatcher.repatch();
258 
259  // Write resulting mesh
260  if (overwrite)
261  {
262  mesh.setInstance(oldInstance);
263  }
264 
265  // Set the precision of the points data to 10
267 
268  mesh.write();
269 
270  Info<< "End\n" << endl;
271 
272  return 0;
273 }
274 
275 
276 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fvMesh::write
virtual bool write(const bool valid=true) const
Definition: fvMesh.C:1034
addOverwriteOption.H
Foam::PrimitivePatch::nEdges
label nEdges() const
Definition: PrimitivePatch.H:318
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
unitConversion.H
Unit conversion functions.
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Definition: polyMesh.H:440
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
polyMesh.H
Foam::cpuTimePosix::cpuTimeIncrement
double cpuTimeIncrement() const
Definition: cpuTimePosix.C:80
Foam::argList::get
T get(const label index) const
Definition: argListI.H:271
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
OFstream.H
Foam::boundaryPatch
Like polyPatch but without reference to mesh. Used in boundaryMesh to hold data on patches....
Definition: boundaryPatch.H:50
Foam::argList::addArgument
static void addArgument(const string &argName, const string &usage="")
Definition: argList.C:294
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Definition: polyMesh.C:839
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
Foam::repatchPolyTopoChanger
A mesh which allows changes in the patch distribution of the boundary faces. The change in patching i...
Definition: repatchPolyTopoChanger.H:50
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
Foam::bMesh
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points)
Definition: bMesh.H:42
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
argList.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
Foam::polyPatch::New
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Definition: polyPatchNew.C:28
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Definition: atmBoundaryLayer.C:26
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Definition: unitConversion.H:44
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Definition: IOstream.H:338
Time.H
setRootCase.H
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:50
Foam::foamVersion::patch
const std::string patch
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
createTime.H
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:59
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
boundaryMesh.H
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:58
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Definition: primitiveMeshI.H:83
Foam::patchIdentifier::name
const word & name() const noexcept
Definition: patchIdentifier.H:131
ListOps.H
Various functions to operate on Lists.
Foam::argList::noParallel
static void noParallel()
Definition: argList.C:503
args
Foam::argList args(argc, argv)
repatchPolyTopoChanger.H
createPolyMesh.H
Required Variables.
Foam::polyMesh::setInstance
void setInstance(const fileName &instance, const IOobject::writeOption wOpt=IOobject::AUTO_WRITE)
Definition: polyMeshIO.C:29
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:75
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258