surfaceToPatch.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  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  surfaceToPatch
29 
30 Group
31  grpSurfaceUtilities
32 
33 Description
34  Reads surface and applies surface regioning to a mesh. Uses boundaryMesh
35  to do the hard work.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "argList.H"
40 #include "Time.H"
41 #include "boundaryMesh.H"
42 #include "polyMesh.H"
43 #include "faceSet.H"
44 #include "polyTopoChange.H"
45 #include "polyModifyFace.H"
46 #include "globalMeshData.H"
47 
48 using namespace Foam;
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 // Adds empty patch if not yet there. Returns patchID.
53 label addPatch(polyMesh& mesh, const word& patchName)
54 {
55  label patchi = mesh.boundaryMesh().findPatchID(patchName);
56 
57  if (patchi == -1)
58  {
60 
61  polyPatchList newPatches(patches.size() + 1);
62 
63  patchi = 0;
64 
65  // Copy all old patches
66  forAll(patches, i)
67  {
68  const polyPatch& pp = patches[i];
69 
70  newPatches.set
71  (
72  patchi,
73  pp.clone
74  (
75  patches,
76  patchi,
77  pp.size(),
78  pp.start()
79  )
80  );
81 
82  patchi++;
83  }
84 
85  // Add zero-sized patch
86  newPatches.set
87  (
88  patchi,
89  new polyPatch
90  (
91  patchName,
92  0,
93  mesh.nFaces(),
94  patchi,
95  patches,
96  polyPatch::typeName
97  )
98  );
99 
101  mesh.addPatches(newPatches);
102 
103  Pout<< "Created patch " << patchName << " at " << patchi << endl;
104  }
105  else
106  {
107  Pout<< "Reusing patch " << patchName << " at " << patchi << endl;
108  }
109 
110  return patchi;
111 }
112 
113 
114 // Repatch single face. Return true if patch changed.
115 bool repatchFace
116 (
117  const polyMesh& mesh,
118  const boundaryMesh& bMesh,
119  const labelList& nearest,
120  const labelList& surfToMeshPatch,
121  const label facei,
122  polyTopoChange& meshMod
123 )
124 {
125  bool changed = false;
126 
127  label bFacei = facei - mesh.nInternalFaces();
128 
129  if (nearest[bFacei] != -1)
130  {
131  // Use boundary mesh one.
132  label bMeshPatchID = bMesh.whichPatch(nearest[bFacei]);
133 
134  label patchID = surfToMeshPatch[bMeshPatchID];
135 
136  if (patchID != mesh.boundaryMesh().whichPatch(facei))
137  {
138  label own = mesh.faceOwner()[facei];
139 
140  label zoneID = mesh.faceZones().whichZone(facei);
141 
142  bool zoneFlip = false;
143 
144  if (zoneID >= 0)
145  {
146  const faceZone& fZone = mesh.faceZones()[zoneID];
147 
148  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
149  }
150 
151  meshMod.setAction
152  (
154  (
155  mesh.faces()[facei],// modified face
156  facei, // label of face being modified
157  own, // owner
158  -1, // neighbour
159  false, // face flip
160  patchID, // patch for face
161  false, // remove from zone
162  zoneID, // zone for face
163  zoneFlip // face flip in zone
164  )
165  );
166 
167  changed = true;
168  }
169  }
170  else
171  {
172  changed = false;
173  }
174  return changed;
175 }
176 
177 
178 
179 int main(int argc, char *argv[])
180 {
182  (
183  "Reads surface and applies surface regioning to a mesh"
184  );
185 
187  argList::addArgument("surfaceFile");
189  (
190  "faceSet",
191  "name",
192  "Only repatch the faces in specified faceSet"
193  );
195  (
196  "tol",
197  "scalar",
198  "Search tolerance as fraction of mesh size (default 1e-3)"
199  );
200 
201  #include "setRootCase.H"
202  #include "createTime.H"
203  #include "createPolyMesh.H"
204 
205  const auto surfName = args.get<fileName>(1);
206 
207  Info<< "Reading surface from " << surfName << " ..." << endl;
208 
209  word setName;
210  const bool readSet = args.readIfPresent("faceSet", setName);
211 
212  if (readSet)
213  {
214  Info<< "Repatching only the faces in faceSet " << setName
215  << " according to nearest surface triangle ..." << endl;
216  }
217  else
218  {
219  Info<< "Patching all boundary faces according to nearest surface"
220  << " triangle ..." << endl;
221  }
222 
223  const scalar searchTol = args.getOrDefault<scalar>("tol", 1e-3);
224 
225  // Get search box. Anything not within this box will not be considered.
226  const boundBox& meshBb = mesh.bounds();
227  const vector searchSpan = searchTol * meshBb.span();
228 
229  Info<< "All boundary faces further away than " << searchTol
230  << " of mesh bounding box " << meshBb
231  << " will keep their patch label ..." << endl;
232 
233 
234  Info<< "Before patching:" << nl
235  << " patch\tsize" << endl;
236 
237  forAll(mesh.boundaryMesh(), patchi)
238  {
239  Info<< " " << mesh.boundaryMesh()[patchi].name() << '\t'
240  << mesh.boundaryMesh()[patchi].size() << nl;
241  }
242  Info<< endl;
243 
244 
246 
247  // Load in the surface.
248  bMesh.readTriSurface(surfName);
249 
250  // Add all the boundaryMesh patches to the mesh.
251  const PtrList<boundaryPatch>& bPatches = bMesh.patches();
252 
253  // Map from surface patch ( = boundaryMesh patch) to polyMesh patch
254  labelList patchMap(bPatches.size());
255 
256  forAll(bPatches, i)
257  {
258  patchMap[i] = addPatch(mesh, bPatches[i].name());
259  }
260 
261  // Obtain nearest face in bMesh for each boundary face in mesh that
262  // is within search span.
263  // Note: should only determine for faceSet if working with that.
264  labelList nearest(bMesh.getNearest(mesh, searchSpan));
265 
266  {
267  // Dump unmatched faces to faceSet for debugging.
268  faceSet unmatchedFaces(mesh, "unmatchedFaces", nearest.size()/100);
269 
270  forAll(nearest, bFacei)
271  {
272  if (nearest[bFacei] == -1)
273  {
274  unmatchedFaces.insert(mesh.nInternalFaces() + bFacei);
275  }
276  }
277 
278  Pout<< "Writing all " << unmatchedFaces.size()
279  << " unmatched faces to faceSet "
280  << unmatchedFaces.name()
281  << endl;
282 
283  unmatchedFaces.write();
284  }
285 
286 
287  polyTopoChange meshMod(mesh);
288 
289  label nChanged = 0;
290 
291  if (readSet)
292  {
293  faceSet faceLabels(mesh, setName);
294  Info<< "Read " << faceLabels.size() << " faces to repatch ..." << endl;
295 
296  for (const label facei : faceLabels)
297  {
298  if (repatchFace(mesh, bMesh, nearest, patchMap, facei, meshMod))
299  {
300  nChanged++;
301  }
302  }
303  }
304  else
305  {
306  forAll(nearest, bFacei)
307  {
308  const label facei = mesh.nInternalFaces() + bFacei;
309 
310  if (repatchFace(mesh, bMesh, nearest, patchMap, facei, meshMod))
311  {
312  ++nChanged;
313  }
314  }
315  }
316 
317  Pout<< "Changed " << nChanged << " boundary faces." << nl << endl;
318 
319  if (nChanged > 0)
320  {
321  meshMod.changeMesh(mesh, false);
322 
323  Info<< "After patching:" << nl
324  << " patch\tsize" << endl;
325 
326  forAll(mesh.boundaryMesh(), patchi)
327  {
328  Info<< " " << mesh.boundaryMesh()[patchi].name() << '\t'
329  << mesh.boundaryMesh()[patchi].size() << endl;
330  }
331  Info<< endl;
332 
333 
334  ++runTime;
335 
336  // Write resulting mesh
337  Info<< "Writing modified mesh to time " << runTime.value() << endl;
338  mesh.write();
339  }
340 
341 
342  Info<< "End\n" << endl;
343 
344  return 0;
345 }
346 
347 
348 // ************************************************************************* //
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
const Type & value() const noexcept
Return const reference to value.
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Required Variables.
A class for handling file names.
Definition: fileName.H:72
A list of face labels.
Definition: faceSet.H:47
Class describing modification of a face.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void noParallel()
Remove the parallel options.
Definition: argList.C:584
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:397
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
label nFaces() const noexcept
Number of mesh faces.
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const labelUList &patchMap, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
labelList faceLabels(nFaceLabels)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
A list of faces which address into the list of points.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:385
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nInternalFaces() const noexcept
Number of internal faces.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label whichPatch(const label meshFacei) const
Return patch index for a given mesh face index. Uses binary search.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1102
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:32
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:346
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Direct mesh changes based on v1.3 polyTopoChange syntax.
const polyBoundaryMesh & patches
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:351
const boolList & flipMap() const noexcept
Return face flip map.
Definition: faceZone.H:367
messageStream Info
Information stream (stdout output on master, null elsewhere)
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:60
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:971
Foam::argList args(argc, argv)
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points) ...
Definition: bMesh.H:39
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.