regionSide.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-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "regionSide.H"
27 #include "meshTools.H"
28 #include "primitiveMesh.H"
29 
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 defineTypeNameAndDebug(regionSide, 0);
37 
38 }
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
42 // Step across edge onto other face on cell
44 (
45  const primitiveMesh& mesh,
46  const label celli,
47  const label facei,
48  const label edgeI
49 )
50 {
51  label f0I, f1I;
52 
53  meshTools::getEdgeFaces(mesh, celli, edgeI, f0I, f1I);
54 
55  if (f0I == facei)
56  {
57  return f1I;
58  }
59  else
60  {
61  return f0I;
62  }
63 }
64 
65 
66 // Step across point to other edge on face
67 Foam::label Foam::regionSide::otherEdge
68 (
69  const primitiveMesh& mesh,
70  const label facei,
71  const label edgeI,
72  const label pointi
73 )
74 {
75  const edge& e = mesh.edges()[edgeI];
76 
77  // Get other point on edge.
78  label freePointi = e.otherVertex(pointi);
79 
80  const labelList& fEdges = mesh.faceEdges()[facei];
81 
82  forAll(fEdges, fEdgeI)
83  {
84  const label otherEdgeI = fEdges[fEdgeI];
85  const edge& otherE = mesh.edges()[otherEdgeI];
86 
87  if
88  (
89  (
90  otherE.start() == pointi
91  && otherE.end() != freePointi
92  )
93  || (
94  otherE.end() == pointi
95  && otherE.start() != freePointi
96  )
97  )
98  {
99  // otherE shares one (but not two) points with e.
100  return otherEdgeI;
101  }
102  }
103 
105  << "Cannot find other edge on face " << facei << " that uses point "
106  << pointi << " but not point " << freePointi << endl
107  << "Edges on face:" << fEdges
108  << " verts:" << UIndirectList<edge>(mesh.edges(), fEdges)()
109  << " Vertices on face:"
110  << mesh.faces()[facei]
111  << " Vertices on original edge:" << e << abort(FatalError);
112 
113  return -1;
114 }
115 
116 
117 // Step from facei (on side celli) to connected face & cell without crossing
118 // fenceEdges.
119 void Foam::regionSide::visitConnectedFaces
120 (
121  const primitiveMesh& mesh,
122  const labelHashSet& region,
123  const labelHashSet& fenceEdges,
124  const label celli,
125  const label facei,
126  labelHashSet& visitedFace
127 )
128 {
129  if (!visitedFace.found(facei))
130  {
131  if (debug)
132  {
133  Info<< "visitConnectedFaces : celli:" << celli << " facei:"
134  << facei << " isOwner:" << (celli == mesh.faceOwner()[facei])
135  << endl;
136  }
137 
138  // Mark as visited
139  visitedFace.insert(facei);
140 
141  // Mark which side of face was visited.
142  if (celli == mesh.faceOwner()[facei])
143  {
144  sideOwner_.insert(facei);
145  }
146 
147 
148  // Visit all neighbouring faces on faceSet. Stay on this 'side' of
149  // face by doing edge-face-cell walk.
150  const labelList& fEdges = mesh.faceEdges()[facei];
151 
152  forAll(fEdges, fEdgeI)
153  {
154  label edgeI = fEdges[fEdgeI];
155 
156  if (!fenceEdges.found(edgeI))
157  {
158  // Step along faces on edge from cell to cell until
159  // we hit face on faceSet.
160 
161  // Find face reachable from edge
162  label otherFacei = otherFace(mesh, celli, facei, edgeI);
163 
164  if (mesh.isInternalFace(otherFacei))
165  {
166  label otherCelli = celli;
167 
168  // Keep on crossing faces/cells until back on face on
169  // surface
170  while (!region.found(otherFacei))
171  {
172  visitedFace.insert(otherFacei);
173 
174  if (debug)
175  {
176  Info<< "visitConnectedFaces : celli:" << celli
177  << " found insideEdgeFace:" << otherFacei
178  << endl;
179  }
180 
181 
182  // Cross otherFacei into neighbouring cell
183  otherCelli =
185  (
186  mesh,
187  otherCelli,
188  otherFacei
189  );
190 
191  otherFacei =
192  otherFace
193  (
194  mesh,
195  otherCelli,
196  otherFacei,
197  edgeI
198  );
199  }
200 
201  visitConnectedFaces
202  (
203  mesh,
204  region,
205  fenceEdges,
206  otherCelli,
207  otherFacei,
208  visitedFace
209  );
210  }
211  }
212  }
213  }
214 }
215 
216 
217 // From edge on face connected to point on region (regionPointi) cross
218 // to all other edges using this point by walking across faces
219 // Does not cross regionEdges so stays on one side
220 // of region
221 void Foam::regionSide::walkPointConnectedFaces
222 (
223  const primitiveMesh& mesh,
224  const labelHashSet& regionEdges,
225  const label regionPointi,
226  const label startFacei,
227  const label startEdgeI,
228  labelHashSet& visitedEdges
229 )
230 {
231  // Mark as visited
232  insidePointFaces_.insert(startFacei);
233 
234  if (debug)
235  {
236  Info<< "walkPointConnectedFaces : regionPointi:" << regionPointi
237  << " facei:" << startFacei
238  << " edgeI:" << startEdgeI << " verts:"
239  << mesh.edges()[startEdgeI]
240  << endl;
241  }
242 
243  // Cross facei i.e. get edge not startEdgeI which uses regionPointi
244  label edgeI = otherEdge(mesh, startFacei, startEdgeI, regionPointi);
245 
246  if (!regionEdges.found(edgeI))
247  {
248  if (!visitedEdges.found(edgeI))
249  {
250  visitedEdges.insert(edgeI);
251 
252  if (debug)
253  {
254  Info<< "Crossed face from "
255  << " edgeI:" << startEdgeI << " verts:"
256  << mesh.edges()[startEdgeI]
257  << " to edge:" << edgeI << " verts:"
258  << mesh.edges()[edgeI]
259  << endl;
260  }
261 
262  // Cross edge to all faces connected to it.
263 
264  const labelList& eFaces = mesh.edgeFaces()[edgeI];
265 
266  forAll(eFaces, eFacei)
267  {
268  label facei = eFaces[eFacei];
269 
270  walkPointConnectedFaces
271  (
272  mesh,
273  regionEdges,
274  regionPointi,
275  facei,
276  edgeI,
277  visitedEdges
278  );
279  }
280  }
281  }
282 }
283 
284 
285 // Find all faces reachable from all non-fence points and staying on
286 // regionFaces side.
287 void Foam::regionSide::walkAllPointConnectedFaces
288 (
289  const primitiveMesh& mesh,
290  const labelHashSet& regionFaces,
291  const labelHashSet& fencePoints
292 )
293 {
294  //
295  // Get all (internal and external) edges on region.
296  //
297  labelHashSet regionEdges(4*regionFaces.size());
298 
299  forAllConstIter(labelHashSet, regionFaces, iter)
300  {
301  const label facei = iter.key();
302  const labelList& fEdges = mesh.faceEdges()[facei];
303 
304  forAll(fEdges, fEdgeI)
305  {
306  regionEdges.insert(fEdges[fEdgeI]);
307  }
308  }
309 
310 
311  //
312  // Visit all internal points on surface.
313  //
314 
315  // Storage for visited points
316  labelHashSet visitedPoint(4*regionFaces.size());
317 
318  // Insert fence points so we don't visit them
319  forAllConstIter(labelHashSet, fencePoints, iter)
320  {
321  visitedPoint.insert(iter.key());
322  }
323 
324  labelHashSet visitedEdges(2*fencePoints.size());
325 
326 
327  if (debug)
328  {
329  Info<< "Excluding visit of points:" << visitedPoint << endl;
330  }
331 
332  forAllConstIter(labelHashSet, regionFaces, iter)
333  {
334  const label facei = iter.key();
335 
336  // Get side of face.
337  label celli;
338 
339  if (sideOwner_.found(facei))
340  {
341  celli = mesh.faceOwner()[facei];
342  }
343  else
344  {
345  celli = mesh.faceNeighbour()[facei];
346  }
347 
348  // Find starting point and edge on face.
349  const labelList& fEdges = mesh.faceEdges()[facei];
350 
351  forAll(fEdges, fEdgeI)
352  {
353  label edgeI = fEdges[fEdgeI];
354 
355  // Get the face 'perpendicular' to facei on region.
356  label otherFacei = otherFace(mesh, celli, facei, edgeI);
357 
358  // Edge
359  const edge& e = mesh.edges()[edgeI];
360 
361  if (!visitedPoint.found(e.start()))
362  {
363  Info<< "Determining visibility from point " << e.start()
364  << endl;
365 
366  visitedPoint.insert(e.start());
367 
368  // edgeI = otherEdge(mesh, otherFacei, edgeI, e.start());
369 
370  walkPointConnectedFaces
371  (
372  mesh,
373  regionEdges,
374  e.start(),
375  otherFacei,
376  edgeI,
377  visitedEdges
378  );
379  }
380  if (!visitedPoint.found(e.end()))
381  {
382  Info<< "Determining visibility from point " << e.end()
383  << endl;
384 
385  visitedPoint.insert(e.end());
386 
387  // edgeI = otherEdge(mesh, otherFacei, edgeI, e.end());
388 
389  walkPointConnectedFaces
390  (
391  mesh,
392  regionEdges,
393  e.end(),
394  otherFacei,
395  edgeI,
396  visitedEdges
397  );
398  }
399  }
400  }
401 }
402 
403 
404 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
405 
406 // Construct from components
408 (
409  const primitiveMesh& mesh,
410  const labelHashSet& region, // faces of region
411  const labelHashSet& fenceEdges, // outside edges
412  const label startCelli,
413  const label startFacei
414 )
415 :
416  sideOwner_(region.size()),
417  insidePointFaces_(region.size())
418 {
419  // Storage for visited faces
420  labelHashSet visitedFace(region.size());
421 
422  // Visit all faces on this side of region.
423  // Mark for each face whether owner (or neighbour) has been visited.
424  // Sets sideOwner_
425  visitConnectedFaces
426  (
427  mesh,
428  region,
429  fenceEdges,
430  startCelli,
431  startFacei,
432  visitedFace
433  );
434 
435  //
436  // Visit all internal points on region and mark faces visitable from these.
437  // Sets insidePointFaces_.
438  //
439 
440  labelHashSet fencePoints(fenceEdges.size());
441 
442  forAllConstIter(labelHashSet, fenceEdges, iter)
443  {
444  const edge& e = mesh.edges()[iter.key()];
445 
446  fencePoints.insert(e.start());
447  fencePoints.insert(e.end());
448  }
449 
450  walkAllPointConnectedFaces(mesh, region, fencePoints);
451 }
452 
453 
454 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
label otherEdge(const primitiveMesh &, const labelList &edgeLabels, const label edgeI, const label vertI)
Return label of other edge (out of candidates edgeLabels)
Definition: meshTools.C:501
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label otherCell(const primitiveMesh &, const label celli, const label facei)
Return cell on other side of face. Throws error.
Definition: meshTools.C:561
static label otherFace(const primitiveMesh &mesh, const label celli, const label excludeFacei, const label edgeI)
Step across edge onto other face on cell.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
regionSide(const primitiveMesh &mesh, const labelHashSet &region, const labelHashSet &fenceEdges, const label startCell, const label startFace)
Construct from components.
defineTypeNameAndDebug(combustionModel, 0)
messageStream Info
void getEdgeFaces(const primitiveMesh &, const label celli, const label edgeI, label &face0, label &face1)
Get faces on cell using edgeI. Throws error if no two found.
Definition: meshTools.C:458
label otherFace(const primitiveMesh &, const label celli, const label facei, const label edgeI)
Return face on cell using edgeI but not facei. Throws error.
Definition: meshTools.C:536
Namespace for OpenFOAM.