surfaceToPatch.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  surfaceToPatch
26 
27 Description
28  Reads surface and applies surface regioning to a mesh. Uses repatchMesh
29  to do the hard work.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "Time.H"
35 #include "repatchMesh.H"
36 #include "polyMesh.H"
37 #include "faceSet.H"
38 #include "polyTopoChange.H"
39 #include "globalMeshData.H"
40 
41 using namespace Foam;
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 // Adds empty patch if not yet there. Returns patchID.
46 label addPatch(polyMesh& mesh, const word& patchName)
47 {
48  label patchi = mesh.boundaryMesh().findIndex(patchName);
49 
50  if (patchi == -1)
51  {
52  const polyBoundaryMesh& patches = mesh.boundaryMesh();
53 
54  List<polyPatch*> newPatches(patches.size() + 1);
55 
56  patchi = 0;
57 
58  // Copy all old patches
59  forAll(patches, i)
60  {
61  const polyPatch& pp = patches[i];
62 
63  newPatches[patchi] =
64  pp.clone
65  (
66  patches,
67  patchi,
68  pp.size(),
69  pp.start()
70  ).ptr();
71 
72  patchi++;
73  }
74 
75  // Add zero-sized patch
76  newPatches[patchi] =
77  new polyPatch
78  (
79  patchName,
80  0,
81  mesh.nFaces(),
82  patchi,
83  patches,
84  polyPatch::typeName
85  );
86 
87  mesh.removeBoundary();
88  mesh.addPatches(newPatches);
89 
90  Pout<< "Created patch " << patchName << " at " << patchi << endl;
91  }
92  else
93  {
94  Pout<< "Reusing patch " << patchName << " at " << patchi << endl;
95  }
96 
97  return patchi;
98 }
99 
100 
101 // Repatch single face. Return true if patch changed.
102 bool repatchFace
103 (
104  const polyMesh& mesh,
105  const repatchMesh& rMesh,
106  const labelList& nearest,
107  const labelList& surfToMeshPatch,
108  const label facei,
109  polyTopoChange& meshMod
110 )
111 {
112  bool changed = false;
113 
114  label bFacei = facei - mesh.nInternalFaces();
115 
116  if (nearest[bFacei] != -1)
117  {
118  // Use boundary mesh one.
119  const label rMeshPatchID = rMesh.whichPatch(nearest[bFacei]);
120 
121  const label patchID = surfToMeshPatch[rMeshPatchID];
122 
123  if (patchID != mesh.boundaryMesh().whichPatch(facei))
124  {
125  const label own = mesh.faceOwner()[facei];
126 
127  meshMod.modifyFace
128  (
129  mesh.faces()[facei],// modified face
130  facei, // label of face being modified
131  own, // owner
132  -1, // neighbour
133  false, // face flip
134  patchID // patch for face
135  );
136 
137  changed = true;
138  }
139  }
140  else
141  {
142  changed = false;
143  }
144 
145  return changed;
146 }
147 
148 
149 
150 int main(int argc, char *argv[])
151 {
153  (
154  "reads surface and applies surface regioning to a mesh"
155  );
156 
158  argList::validArgs.append("surface file");
160  (
161  "faceSet",
162  "name",
163  "only repatch the faces in specified faceSet"
164  );
166  (
167  "tol",
168  "scalar",
169  "search tolerance as fraction of mesh size (default 1e-3)"
170  );
171 
172  #include "setRootCase.H"
173  #include "createTime.H"
174  #include "createPolyMesh.H"
175 
176  const fileName surfName = args[1];
177 
178  Info<< "Reading surface from " << surfName << " ..." << endl;
179 
180  word setName;
181  const bool readSet = args.optionReadIfPresent("faceSet", setName);
182 
183  if (readSet)
184  {
185  Info<< "Repatching only the faces in faceSet " << setName
186  << " according to nearest surface triangle ..." << endl;
187  }
188  else
189  {
190  Info<< "Patching all boundary faces according to nearest surface"
191  << " triangle ..." << endl;
192  }
193 
194  const scalar searchTol = args.optionLookupOrDefault("tol", 1e-3);
195 
196  // Get search box. Anything not within this box will not be considered.
197  const boundBox& meshBb = mesh.bounds();
198  const vector searchSpan = searchTol * meshBb.span();
199 
200  Info<< "All boundary faces further away than " << searchTol
201  << " of mesh bounding box " << meshBb
202  << " will keep their patch label ..." << endl;
203 
204 
205  Info<< "Before patching:" << nl
206  << " patch\tsize" << endl;
207 
208  forAll(mesh.boundaryMesh(), patchi)
209  {
210  Info<< " " << mesh.boundaryMesh()[patchi].name() << '\t'
211  << mesh.boundaryMesh()[patchi].size() << nl;
212  }
213  Info<< endl;
214 
215 
216  repatchMesh rMesh;
217 
218  // Load in the surface.
219  rMesh.readTriSurface(surfName);
220 
221  // Add all the boundaryMesh patches to the mesh.
222  const PtrList<repatchPatch>& bPatches = rMesh.patches();
223 
224  // Map from surface patch ( = boundaryMesh patch) to polyMesh patch
225  labelList patchMap(bPatches.size());
226 
227  forAll(bPatches, i)
228  {
229  patchMap[i] = addPatch(mesh, bPatches[i].name());
230  }
231 
232  // Obtain nearest face in rMesh for each boundary face in mesh that
233  // is within search span.
234  // Note: should only determine for faceSet if working with that.
235  labelList nearest(rMesh.getNearest(mesh, searchSpan));
236 
237  {
238  // Dump unmatched faces to faceSet for debugging.
239  faceSet unmatchedFaces(mesh, "unmatchedFaces", nearest.size()/100);
240 
241  forAll(nearest, bFacei)
242  {
243  if (nearest[bFacei] == -1)
244  {
245  unmatchedFaces.insert(mesh.nInternalFaces() + bFacei);
246  }
247  }
248 
249  Pout<< "Writing all " << unmatchedFaces.size()
250  << " unmatched faces to faceSet "
251  << unmatchedFaces.name()
252  << endl;
253 
254  unmatchedFaces.write();
255  }
256 
257 
258  polyTopoChange meshMod(mesh);
259 
260  label nChanged = 0;
261 
262  if (readSet)
263  {
264  faceSet faceLabels(mesh, setName);
265  Info<< "Read " << faceLabels.size() << " faces to repatch ..." << endl;
266 
267  forAllConstIter(faceSet, faceLabels, iter)
268  {
269  label facei = iter.key();
270 
271  if (repatchFace(mesh, rMesh, nearest, patchMap, facei, meshMod))
272  {
273  nChanged++;
274  }
275  }
276  }
277  else
278  {
279  forAll(nearest, bFacei)
280  {
281  label facei = mesh.nInternalFaces() + bFacei;
282 
283  if (repatchFace(mesh, rMesh, nearest, patchMap, facei, meshMod))
284  {
285  nChanged++;
286  }
287  }
288  }
289 
290  Pout<< "Changed " << nChanged << " boundary faces." << nl << endl;
291 
292  if (nChanged > 0)
293  {
294  meshMod.changeMesh(mesh);
295 
296  Info<< "After patching:" << nl
297  << " patch\tsize" << endl;
298 
299  forAll(mesh.boundaryMesh(), patchi)
300  {
301  Info<< " " << mesh.boundaryMesh()[patchi].name() << '\t'
302  << mesh.boundaryMesh()[patchi].size() << endl;
303  }
304  Info<< endl;
305 
306 
307  runTime++;
308 
309  // Write resulting mesh
310  Info<< "Writing modified mesh to time "
311  << runTime.userTimeName() << endl;
312  mesh.write();
313  }
314 
315 
316  Info<< "End\n" << endl;
317 
318  return 0;
319 }
320 
321 
322 // ************************************************************************* //
#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
const word & name() const
Return name.
Definition: IOobject.H:310
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:85
virtual Ostream & write(const char)=0
Write character.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
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 addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:243
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
A list of face labels.
Definition: faceSet.H:51
A class for handling file names.
Definition: fileName.H:82
Foam::polyBoundaryMesh.
label findIndex(const word &patchName) const
Find patch index given a name.
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
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1084
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:410
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:215
Direct mesh changes based on v1.3 polyTopoChange syntax.
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.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: repatchMesh.H:58
void readTriSurface(const fileName &)
Read from triSurface.
Definition: repatchMesh.C:474
labelList getNearest(const primitiveMesh &pMesh, const vector &searchSpan) const
Get rMesh index of nearest face for every boundary face in.
Definition: repatchMesh.C:651
label whichPatch(const label facei) const
Get index of patch face is in.
Definition: repatchMesh.C:1106
const PtrList< repatchPatch > & patches() const
Access the patches.
Definition: repatchMesh.H:192
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:266
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(1e-3))
Foam::argList args(argc, argv)