cellFeatures.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) 2019-2021 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 \*---------------------------------------------------------------------------*/
28 
29 #include "cellFeatures.H"
30 #include "primitiveMesh.H"
31 #include "HashSet.H"
32 #include "Map.H"
33 #include "ListOps.H"
34 #include "meshTools.H"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 
39 // Return true if edge start and end are on increasing face vertices. (edge is
40 // guaranteed to be on face)
41 bool Foam::cellFeatures::faceAlignedEdge(const label facei, const label edgeI)
42  const
43 {
44  const edge& e = mesh_.edges()[edgeI];
45 
46  const face& f = mesh_.faces()[facei];
47 
48  forAll(f, fp)
49  {
50  if (f[fp] == e.start())
51  {
52  label fp1 = f.fcIndex(fp);
53 
54  return f[fp1] == e.end();
55  }
56  }
57 
59  << "Can not find edge " << mesh_.edges()[edgeI]
60  << " on face " << facei << abort(FatalError);
61 
62  return false;
63 }
64 
65 
66 // Return edge in featureEdge that uses vertI and is on same superface
67 // but is not edgeI
68 Foam::label Foam::cellFeatures::nextEdge
69 (
70  const Map<label>& toSuperFace,
71  const label superFacei,
72  const label thisEdgeI,
73  const label thisVertI
74 ) const
75 {
76  const labelList& pEdges = mesh_.pointEdges()[thisVertI];
77 
78  forAll(pEdges, pEdgeI)
79  {
80  label edgeI = pEdges[pEdgeI];
81 
82  if ((edgeI != thisEdgeI) && featureEdge_.found(edgeI))
83  {
84  // Check that edge is used by a face on same superFace
85 
86  const labelList& eFaces = mesh_.edgeFaces()[edgeI];
87 
88  forAll(eFaces, eFacei)
89  {
90  label facei = eFaces[eFacei];
91 
92  if
93  (
94  meshTools::faceOnCell(mesh_, celli_, facei)
95  && (toSuperFace[facei] == superFacei)
96  )
97  {
98  return edgeI;
99  }
100  }
101  }
102  }
103 
105  << "Can not find edge in " << featureEdge_ << " connected to edge "
106  << thisEdgeI << " at vertex " << thisVertI << endl
107  << "This might mean that the externalEdges do not form a closed loop"
108  << abort(FatalError);
109 
110  return -1;
111 }
112 
113 
114 // Return true if angle between faces using it is larger than certain value.
115 bool Foam::cellFeatures::isCellFeatureEdge
116 (
117  const scalar minCos,
118  const label edgeI
119 ) const
120 {
121  // Get the two faces using this edge
122 
123  label face0;
124  label face1;
125  meshTools::getEdgeFaces(mesh_, celli_, edgeI, face0, face1);
126 
127  // Check the angle between them by comparing the face normals.
128 
129  const vector n0 = normalised(mesh_.faceAreas()[face0]);
130  const vector n1 = normalised(mesh_.faceAreas()[face1]);
131 
132  scalar cosAngle = n0 & n1;
133 
134  const edge& e = mesh_.edges()[edgeI];
135 
136  const face& f0 = mesh_.faces()[face0];
137 
138  label face0Start = f0.find(e.start());
139  label face0End = f0.fcIndex(face0Start);
140 
141  const face& f1 = mesh_.faces()[face1];
142 
143  label face1Start = f1.find(e.start());
144  label face1End = f1.fcIndex(face1Start);
145 
146  if
147  (
148  (
149  (f0[face0End] == e.end())
150  && (f1[face1End] != e.end())
151  )
152  || (
153  (f0[face0End] != e.end())
154  && (f1[face1End] == e.end())
155  )
156  )
157  {
158  }
159  else
160  {
161  cosAngle = -cosAngle;
162  }
163 
164  if (cosAngle < minCos)
165  {
166  return true;
167  }
168 
169  return false;
170 }
171 
172 
173 // Recursively mark (on toSuperFace) all face reachable across non-feature
174 // edges.
175 void Foam::cellFeatures::walkSuperFace
176 (
177  const label facei,
178  const label superFacei,
179  Map<label>& toSuperFace
180 ) const
181 {
182  if (toSuperFace.insert(facei, superFacei))
183  {
184  const labelList& fEdges = mesh_.faceEdges()[facei];
185 
186  for (const label edgeI : fEdges)
187  {
188  if (!featureEdge_.found(edgeI))
189  {
190  label face0;
191  label face1;
192  meshTools::getEdgeFaces(mesh_, celli_, edgeI, face0, face1);
193 
194  if (face0 == facei)
195  {
196  face0 = face1;
197  }
198 
199  walkSuperFace
200  (
201  face0,
202  superFacei,
203  toSuperFace
204  );
205  }
206  }
207  }
208 }
209 
210 
211 void Foam::cellFeatures::calcSuperFaces() const
212 {
213  // Determine superfaces by edge walking across non-feature edges
214 
215  const labelList& cFaces = mesh_.cells()[celli_];
216 
217  // Mapping from old to super face:
218  // <not found> : not visited
219  // >=0 : superFace
220  Map<label> toSuperFace(10*cFaces.size());
221 
222  label superFacei = 0;
223 
224  forAll(cFaces, cFacei)
225  {
226  label facei = cFaces[cFacei];
227 
228  if (!toSuperFace.found(facei))
229  {
230  walkSuperFace
231  (
232  facei,
233  superFacei,
234  toSuperFace
235  );
236  superFacei++;
237  }
238  }
239 
240  // Construct superFace-to-oldface mapping.
241 
242  faceMap_.setSize(superFacei);
243 
244  forAll(cFaces, cFacei)
245  {
246  label facei = cFaces[cFacei];
247 
248  faceMap_[toSuperFace[facei]].append(facei);
249  }
250 
251  forAll(faceMap_, superI)
252  {
253  faceMap_[superI].shrink();
254  }
255 
256 
257  // Construct superFaces
258 
259  facesPtr_.reset(new faceList(superFacei));
260  auto& faces = *facesPtr_;
261 
262  forAll(cFaces, cFacei)
263  {
264  label facei = cFaces[cFacei];
265 
266  label superFacei = toSuperFace[facei];
267 
268  if (faces[superFacei].empty())
269  {
270  // Superface not yet constructed.
271 
272  // Find starting feature edge on face.
273  label startEdgeI = -1;
274 
275  const labelList& fEdges = mesh_.faceEdges()[facei];
276 
277  forAll(fEdges, fEdgeI)
278  {
279  label edgeI = fEdges[fEdgeI];
280 
281  if (featureEdge_.found(edgeI))
282  {
283  startEdgeI = edgeI;
284 
285  break;
286  }
287  }
288 
289 
290  if (startEdgeI != -1)
291  {
292  // Walk point-edge-point along feature edges
293 
294  DynamicList<label> superFace(10*mesh_.faces()[facei].size());
295 
296  const edge& e = mesh_.edges()[startEdgeI];
297 
298  // Walk either start-end or end-start depending on orientation
299  // of face. SuperFace will have celli as owner.
300  bool flipOrientation =
301  (mesh_.faceOwner()[facei] == celli_)
302  ^ (faceAlignedEdge(facei, startEdgeI));
303 
304  label startVertI = -1;
305 
306  if (flipOrientation)
307  {
308  startVertI = e.end();
309  }
310  else
311  {
312  startVertI = e.start();
313  }
314 
315  label edgeI = startEdgeI;
316 
317  label vertI = e.otherVertex(startVertI);
318 
319  do
320  {
321  label newEdgeI = nextEdge
322  (
323  toSuperFace,
324  superFacei,
325  edgeI,
326  vertI
327  );
328 
329  // Determine angle between edges.
330  if (isFeaturePoint(edgeI, newEdgeI))
331  {
332  superFace.append(vertI);
333  }
334 
335  edgeI = newEdgeI;
336 
337  if (vertI == startVertI)
338  {
339  break;
340  }
341 
342  vertI = mesh_.edges()[edgeI].otherVertex(vertI);
343  }
344  while (true);
345 
346  if (superFace.size() <= 2)
347  {
349  << " Can not collapse faces " << faceMap_[superFacei]
350  << " into one big face on cell " << celli_ << endl
351  << "Try decreasing minCos:" << minCos_ << endl;
352  }
353  else
354  {
355  faces[superFacei].transfer(superFace);
356  }
357  }
358  }
359  }
360 }
361 
362 
363 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
364 
365 // Construct from components
366 Foam::cellFeatures::cellFeatures
367 (
368  const primitiveMesh& mesh,
369  const scalar minCos,
370  const label celli
371 )
372 :
373  mesh_(mesh),
374  minCos_(minCos),
375  celli_(celli),
376  featureEdge_(10*mesh.cellEdges()[celli].size()),
377  facesPtr_(nullptr),
378  faceMap_(0)
379 {
380  const labelList& cEdges = mesh_.cellEdges()[celli_];
381 
382  forAll(cEdges, cEdgeI)
383  {
384  label edgeI = cEdges[cEdgeI];
385 
386  if (isCellFeatureEdge(minCos_, edgeI))
387  {
388  featureEdge_.insert(edgeI);
389  }
390  }
391 }
392 
393 
394 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
397 {}
398 
399 
400 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
401 
402 bool Foam::cellFeatures::isFeaturePoint(const label edge0, const label edge1)
403  const
404 {
405  if
406  (
407  (edge0 < 0)
408  || (edge0 >= mesh_.nEdges())
409  || (edge1 < 0)
410  || (edge1 >= mesh_.nEdges())
411  )
412  {
414  << "Illegal edge labels : edge0:" << edge0 << " edge1:" << edge1
415  << abort(FatalError);
416  }
417 
418  const edge& e0 = mesh_.edges()[edge0];
419 
420  const vector e0Vec = e0.unitVec(mesh_.points());
421 
422  const edge& e1 = mesh_.edges()[edge1];
423 
424  const vector e1Vec = e1.unitVec(mesh_.points());
425 
426  scalar cosAngle;
427 
428  if
429  (
430  (e0.start() == e1.end())
431  || (e0.end() == e1.start())
432  )
433  {
434  // Same direction
435  cosAngle = e0Vec & e1Vec;
436  }
437  else if
438  (
439  (e0.start() == e1.start())
440  || (e0.end() == e1.end())
441  )
442  {
443  // back on back
444  cosAngle = - e0Vec & e1Vec;
445  }
446  else
447  {
448  cosAngle = GREAT; // satisfy compiler
449 
451  << "Edges do not share common vertex. e0:" << e0
452  << " e1:" << e1 << abort(FatalError);
453  }
454 
455  if (cosAngle < minCos_)
456  {
457  // Angle larger than criterion
458  return true;
459  }
460 
461  return false;
462 }
463 
464 
466 (
467  const label facei,
468  const label vertI
469 ) const
470 {
471  if
472  (
473  (facei < 0)
474  || (facei >= mesh_.nFaces())
475  || (vertI < 0)
476  || (vertI >= mesh_.nPoints())
477  )
478  {
480  << "Illegal face " << facei << " or vertex " << vertI
481  << abort(FatalError);
482  }
483 
484  const labelList& pEdges = mesh_.pointEdges()[vertI];
485 
486  label edge0 = -1;
487  label edge1 = -1;
488 
489  forAll(pEdges, pEdgeI)
490  {
491  label edgeI = pEdges[pEdgeI];
492 
493  if (meshTools::edgeOnFace(mesh_, facei, edgeI))
494  {
495  if (edge0 == -1)
496  {
497  edge0 = edgeI;
498  }
499  else
500  {
501  edge1 = edgeI;
502 
503  // Found the two edges.
504  break;
505  }
506  }
507  }
508 
509  if (edge1 == -1)
510  {
512  << "Did not find two edges sharing vertex " << vertI
513  << " on face " << facei << " vertices:" << mesh_.faces()[facei]
514  << abort(FatalError);
515  }
516 
517  return isFeaturePoint(edge0, edge1);
518 }
519 
520 
521 // ************************************************************************* //
const labelListList & cellEdges() const
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
bool isFeatureVertex(const label facei, const label vertI) const
Is vertexI on facei used by two edges that form feature.
Definition: cellFeatures.C:459
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool edgeOnFace(const primitiveMesh &mesh, const label facei, const label edgeI)
Is edge used by face.
Definition: meshTools.C:312
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
~cellFeatures()
Destructor.
Definition: cellFeatures.C:389
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
bool isFeaturePoint(const label edge0, const label edge1) const
Are two edges connected at feature point?
Definition: cellFeatures.C:395
labelList f(nPoints)
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
#define WarningInFunction
Report a warning using Foam::Warning.
virtual const faceList & faces() const =0
Return faces.
List< label > labelList
A List of labels.
Definition: List.H:62
bool faceOnCell(const primitiveMesh &mesh, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:323
void getEdgeFaces(const primitiveMesh &mesh, 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:472