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 | Copyright (C) 2011-2016 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 boundaryMesh
29  to do the hard work.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "Time.H"
35 #include "boundaryMesh.H"
36 #include "polyMesh.H"
37 #include "faceSet.H"
38 #include "polyTopoChange.H"
39 #include "polyModifyFace.H"
40 #include "globalMeshData.H"
41 
42 using namespace Foam;
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 // Adds empty patch if not yet there. Returns patchID.
47 label addPatch(polyMesh& mesh, const word& patchName)
48 {
49  label patchi = mesh.boundaryMesh().findPatchID(patchName);
50 
51  if (patchi == -1)
52  {
53  const polyBoundaryMesh& patches = mesh.boundaryMesh();
54 
55  List<polyPatch*> newPatches(patches.size() + 1);
56 
57  patchi = 0;
58 
59  // Copy all old patches
60  forAll(patches, i)
61  {
62  const polyPatch& pp = patches[i];
63 
64  newPatches[patchi] =
65  pp.clone
66  (
67  patches,
68  patchi,
69  pp.size(),
70  pp.start()
71  ).ptr();
72 
73  patchi++;
74  }
75 
76  // Add zero-sized patch
77  newPatches[patchi] =
78  new polyPatch
79  (
80  patchName,
81  0,
82  mesh.nFaces(),
83  patchi,
84  patches,
85  polyPatch::typeName
86  );
87 
88  mesh.removeBoundary();
89  mesh.addPatches(newPatches);
90 
91  Pout<< "Created patch " << patchName << " at " << patchi << endl;
92  }
93  else
94  {
95  Pout<< "Reusing patch " << patchName << " at " << patchi << endl;
96  }
97 
98  return patchi;
99 }
100 
101 
102 // Repatch single face. Return true if patch changed.
103 bool repatchFace
104 (
105  const polyMesh& mesh,
106  const boundaryMesh& bMesh,
107  const labelList& nearest,
108  const labelList& surfToMeshPatch,
109  const label facei,
110  polyTopoChange& meshMod
111 )
112 {
113  bool changed = false;
114 
115  label bFacei = facei - mesh.nInternalFaces();
116 
117  if (nearest[bFacei] != -1)
118  {
119  // Use boundary mesh one.
120  label bMeshPatchID = bMesh.whichPatch(nearest[bFacei]);
121 
122  label patchID = surfToMeshPatch[bMeshPatchID];
123 
124  if (patchID != mesh.boundaryMesh().whichPatch(facei))
125  {
126  label own = mesh.faceOwner()[facei];
127 
128  label zoneID = mesh.faceZones().whichZone(facei);
129 
130  bool zoneFlip = false;
131 
132  if (zoneID >= 0)
133  {
134  const faceZone& fZone = mesh.faceZones()[zoneID];
135 
136  zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
137  }
138 
139  meshMod.setAction
140  (
142  (
143  mesh.faces()[facei],// modified face
144  facei, // label of face being modified
145  own, // owner
146  -1, // neighbour
147  false, // face flip
148  patchID, // patch for face
149  false, // remove from zone
150  zoneID, // zone for face
151  zoneFlip // face flip in zone
152  )
153  );
154 
155  changed = true;
156  }
157  }
158  else
159  {
160  changed = false;
161  }
162  return changed;
163 }
164 
165 
166 
167 int main(int argc, char *argv[])
168 {
170  (
171  "reads surface and applies surface regioning to a mesh"
172  );
173 
175  argList::validArgs.append("surfaceFile");
177  (
178  "faceSet",
179  "name",
180  "only repatch the faces in specified faceSet"
181  );
183  (
184  "tol",
185  "scalar",
186  "search tolerance as fraction of mesh size (default 1e-3)"
187  );
188 
189  #include "setRootCase.H"
190  #include "createTime.H"
191  #include "createPolyMesh.H"
192 
193  const fileName surfName = args[1];
194 
195  Info<< "Reading surface from " << surfName << " ..." << endl;
196 
197  word setName;
198  const bool readSet = args.optionReadIfPresent("faceSet", setName);
199 
200  if (readSet)
201  {
202  Info<< "Repatching only the faces in faceSet " << setName
203  << " according to nearest surface triangle ..." << endl;
204  }
205  else
206  {
207  Info<< "Patching all boundary faces according to nearest surface"
208  << " triangle ..." << endl;
209  }
210 
211  const scalar searchTol = args.optionLookupOrDefault("tol", 1e-3);
212 
213  // Get search box. Anything not within this box will not be considered.
214  const boundBox& meshBb = mesh.bounds();
215  const vector searchSpan = searchTol * meshBb.span();
216 
217  Info<< "All boundary faces further away than " << searchTol
218  << " of mesh bounding box " << meshBb
219  << " will keep their patch label ..." << endl;
220 
221 
222  Info<< "Before patching:" << nl
223  << " patch\tsize" << endl;
224 
225  forAll(mesh.boundaryMesh(), patchi)
226  {
227  Info<< " " << mesh.boundaryMesh()[patchi].name() << '\t'
228  << mesh.boundaryMesh()[patchi].size() << nl;
229  }
230  Info<< endl;
231 
232 
234 
235  // Load in the surface.
236  bMesh.readTriSurface(surfName);
237 
238  // Add all the boundaryMesh patches to the mesh.
239  const PtrList<boundaryPatch>& bPatches = bMesh.patches();
240 
241  // Map from surface patch ( = boundaryMesh patch) to polyMesh patch
242  labelList patchMap(bPatches.size());
243 
244  forAll(bPatches, i)
245  {
246  patchMap[i] = addPatch(mesh, bPatches[i].name());
247  }
248 
249  // Obtain nearest face in bMesh for each boundary face in mesh that
250  // is within search span.
251  // Note: should only determine for faceSet if working with that.
252  labelList nearest(bMesh.getNearest(mesh, searchSpan));
253 
254  {
255  // Dump unmatched faces to faceSet for debugging.
256  faceSet unmatchedFaces(mesh, "unmatchedFaces", nearest.size()/100);
257 
258  forAll(nearest, bFacei)
259  {
260  if (nearest[bFacei] == -1)
261  {
262  unmatchedFaces.insert(mesh.nInternalFaces() + bFacei);
263  }
264  }
265 
266  Pout<< "Writing all " << unmatchedFaces.size()
267  << " unmatched faces to faceSet "
268  << unmatchedFaces.name()
269  << endl;
270 
271  unmatchedFaces.write();
272  }
273 
274 
275  polyTopoChange meshMod(mesh);
276 
277  label nChanged = 0;
278 
279  if (readSet)
280  {
281  faceSet faceLabels(mesh, setName);
282  Info<< "Read " << faceLabels.size() << " faces to repatch ..." << endl;
283 
284  forAllConstIter(faceSet, faceLabels, iter)
285  {
286  label facei = iter.key();
287 
288  if (repatchFace(mesh, bMesh, nearest, patchMap, facei, meshMod))
289  {
290  nChanged++;
291  }
292  }
293  }
294  else
295  {
296  forAll(nearest, bFacei)
297  {
298  label facei = mesh.nInternalFaces() + bFacei;
299 
300  if (repatchFace(mesh, bMesh, nearest, patchMap, facei, meshMod))
301  {
302  nChanged++;
303  }
304  }
305  }
306 
307  Pout<< "Changed " << nChanged << " boundary faces." << nl << endl;
308 
309  if (nChanged > 0)
310  {
311  meshMod.changeMesh(mesh, false);
312 
313  Info<< "After patching:" << nl
314  << " patch\tsize" << endl;
315 
316  forAll(mesh.boundaryMesh(), patchi)
317  {
318  Info<< " " << mesh.boundaryMesh()[patchi].name() << '\t'
319  << mesh.boundaryMesh()[patchi].size() << endl;
320  }
321  Info<< endl;
322 
323 
324  runTime++;
325 
326  // Write resulting mesh
327  Info<< "Writing modified mesh to time " << runTime.value() << endl;
328  mesh.write();
329  }
330 
331 
332  Info<< "End\n" << endl;
333 
334  return 0;
335 }
336 
337 
338 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A class for handling file names.
Definition: fileName.H:69
A list of face labels.
Definition: faceSet.H:48
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:219
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:253
Class describing modification of a face.
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
const double e
Elementary charge.
Definition: doubleFloat.H:78
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
patches[0]
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
T optionLookupOrDefault(const word &opt, const T &deflt) const
Read a value from the named option if present.
Definition: argListI.H:237
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:305
labelList getNearest(const primitiveMesh &pMesh, const vector &searchSpan) const
Get bMesh index of nearest face for every boundary face in.
Definition: boundaryMesh.C:861
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
A class for handling words, derived from string.
Definition: word.H:59
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:93
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Foam::polyBoundaryMesh.
void removeBoundary()
Remove boundary patches.
Definition: polyMeshClear.C:36
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const char nl
Definition: Ostream.H:262
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:878
label nFaces() const
label patchi
void readTriSurface(const fileName &)
Read from triSurface.
Definition: boundaryMesh.C:602
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
label whichZone(const label objectIndex) const
Given a global object index, return the zone it is in.
Definition: ZoneMesh.C:221
Direct mesh changes based on v1.3 polyTopoChange syntax.
virtual bool write() const
Write using setting from DB.
messageStream Info
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
const PtrList< boundaryPatch > & patches() const
Definition: boundaryMesh.H:213
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
autoPtr< mapPolyMesh > 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.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::argList args(argc, argv)
label findPatchID(const word &patchName) const
Find patch index given a name.
PrimitivePatch< face, List, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points) ...
Definition: bMesh.H:45
label nInternalFaces() const
const word & name() const
Return name.
Definition: IOobject.H:260
Namespace for OpenFOAM.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.