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