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-------------------------------------------------------------------------------
11License
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)
42bool 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
69Foam::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.
116bool 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.
176void Foam::cellFeatures::walkSuperFace
177(
178 const label facei,
179 const label superFacei,
180 Map<label>& toSuperFace
181) const
182{
183 if (toSuperFace.insert(facei, superFacei))
184 {
185 const labelList& fEdges = mesh_.faceEdges()[facei];
186
187 for (const label edgeI : fEdges)
188 {
189 if (!featureEdge_.found(edgeI))
190 {
191 label face0;
192 label face1;
193 meshTools::getEdgeFaces(mesh_, celli_, edgeI, face0, face1);
194
195 if (face0 == facei)
196 {
197 face0 = face1;
198 }
199
200 walkSuperFace
201 (
202 face0,
203 superFacei,
204 toSuperFace
205 );
206 }
207 }
208 }
209}
210
211
212void Foam::cellFeatures::calcSuperFaces() const
213{
214 // Determine superfaces by edge walking across non-feature edges
215
216 const labelList& cFaces = mesh_.cells()[celli_];
217
218 // Mapping from old to super face:
219 // <not found> : not visited
220 // >=0 : superFace
221 Map<label> toSuperFace(10*cFaces.size());
222
223 label superFacei = 0;
224
225 forAll(cFaces, cFacei)
226 {
227 label facei = cFaces[cFacei];
228
229 if (!toSuperFace.found(facei))
230 {
231 walkSuperFace
232 (
233 facei,
234 superFacei,
235 toSuperFace
236 );
237 superFacei++;
238 }
239 }
240
241 // Construct superFace-to-oldface mapping.
242
243 faceMap_.setSize(superFacei);
244
245 forAll(cFaces, cFacei)
246 {
247 label facei = cFaces[cFacei];
248
249 faceMap_[toSuperFace[facei]].append(facei);
250 }
251
252 forAll(faceMap_, superI)
253 {
254 faceMap_[superI].shrink();
255 }
256
257
258 // Construct superFaces
259
260 facesPtr_ = new faceList(superFacei);
261
262 faceList& faces = *facesPtr_;
263
264 forAll(cFaces, cFacei)
265 {
266 label facei = cFaces[cFacei];
267
268 label superFacei = toSuperFace[facei];
269
270 if (faces[superFacei].empty())
271 {
272 // Superface not yet constructed.
273
274 // Find starting feature edge on face.
275 label startEdgeI = -1;
276
277 const labelList& fEdges = mesh_.faceEdges()[facei];
278
279 forAll(fEdges, fEdgeI)
280 {
281 label edgeI = fEdges[fEdgeI];
282
283 if (featureEdge_.found(edgeI))
284 {
285 startEdgeI = edgeI;
286
287 break;
288 }
289 }
290
291
292 if (startEdgeI != -1)
293 {
294 // Walk point-edge-point along feature edges
295
296 DynamicList<label> superFace(10*mesh_.faces()[facei].size());
297
298 const edge& e = mesh_.edges()[startEdgeI];
299
300 // Walk either start-end or end-start depending on orientation
301 // of face. SuperFace will have celli as owner.
302 bool flipOrientation =
303 (mesh_.faceOwner()[facei] == celli_)
304 ^ (faceAlignedEdge(facei, startEdgeI));
305
306 label startVertI = -1;
307
308 if (flipOrientation)
309 {
310 startVertI = e.end();
311 }
312 else
313 {
314 startVertI = e.start();
315 }
316
317 label edgeI = startEdgeI;
318
319 label vertI = e.otherVertex(startVertI);
320
321 do
322 {
323 label newEdgeI = nextEdge
324 (
325 toSuperFace,
326 superFacei,
327 edgeI,
328 vertI
329 );
330
331 // Determine angle between edges.
332 if (isFeaturePoint(edgeI, newEdgeI))
333 {
334 superFace.append(vertI);
335 }
336
337 edgeI = newEdgeI;
338
339 if (vertI == startVertI)
340 {
341 break;
342 }
343
344 vertI = mesh_.edges()[edgeI].otherVertex(vertI);
345 }
346 while (true);
347
348 if (superFace.size() <= 2)
349 {
351 << " Can not collapse faces " << faceMap_[superFacei]
352 << " into one big face on cell " << celli_ << endl
353 << "Try decreasing minCos:" << minCos_ << endl;
354 }
355 else
356 {
357 faces[superFacei].transfer(superFace);
358 }
359 }
360 }
361 }
362}
363
364
365// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
366
367// Construct from components
369(
370 const primitiveMesh& mesh,
371 const scalar minCos,
372 const label celli
373)
374:
375 mesh_(mesh),
376 minCos_(minCos),
377 celli_(celli),
378 featureEdge_(10*mesh.cellEdges()[celli].size()),
379 facesPtr_(nullptr),
380 faceMap_(0)
381{
382 const labelList& cEdges = mesh_.cellEdges()[celli_];
383
384 forAll(cEdges, cEdgeI)
385 {
386 label edgeI = cEdges[cEdgeI];
387
388 if (isCellFeatureEdge(minCos_, edgeI))
389 {
390 featureEdge_.insert(edgeI);
391 }
392 }
393
394}
395
396
397// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
398
400{
401 deleteDemandDrivenData(facesPtr_);
402}
403
404
405// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
406
407bool Foam::cellFeatures::isFeaturePoint(const label edge0, const label edge1)
408 const
409{
410 if
411 (
412 (edge0 < 0)
413 || (edge0 >= mesh_.nEdges())
414 || (edge1 < 0)
415 || (edge1 >= mesh_.nEdges())
416 )
417 {
419 << "Illegal edge labels : edge0:" << edge0 << " edge1:" << edge1
420 << abort(FatalError);
421 }
422
423 const edge& e0 = mesh_.edges()[edge0];
424
425 const vector e0Vec = e0.unitVec(mesh_.points());
426
427 const edge& e1 = mesh_.edges()[edge1];
428
429 const vector e1Vec = e1.unitVec(mesh_.points());
430
431 scalar cosAngle;
432
433 if
434 (
435 (e0.start() == e1.end())
436 || (e0.end() == e1.start())
437 )
438 {
439 // Same direction
440 cosAngle = e0Vec & e1Vec;
441 }
442 else if
443 (
444 (e0.start() == e1.start())
445 || (e0.end() == e1.end())
446 )
447 {
448 // back on back
449 cosAngle = - e0Vec & e1Vec;
450 }
451 else
452 {
453 cosAngle = GREAT; // satisfy compiler
454
456 << "Edges do not share common vertex. e0:" << e0
457 << " e1:" << e1 << abort(FatalError);
458 }
459
460 if (cosAngle < minCos_)
461 {
462 // Angle larger than criterion
463 return true;
464 }
465
466 return false;
467}
468
469
471(
472 const label facei,
473 const label vertI
474) const
475{
476 if
477 (
478 (facei < 0)
479 || (facei >= mesh_.nFaces())
480 || (vertI < 0)
481 || (vertI >= mesh_.nPoints())
482 )
483 {
485 << "Illegal face " << facei << " or vertex " << vertI
486 << abort(FatalError);
487 }
488
489 const labelList& pEdges = mesh_.pointEdges()[vertI];
490
491 label edge0 = -1;
492 label edge1 = -1;
493
494 forAll(pEdges, pEdgeI)
495 {
496 label edgeI = pEdges[pEdgeI];
497
498 if (meshTools::edgeOnFace(mesh_, facei, edgeI))
499 {
500 if (edge0 == -1)
501 {
502 edge0 = edgeI;
503 }
504 else
505 {
506 edge1 = edgeI;
507
508 // Found the two edges.
509 break;
510 }
511 }
512 }
513
514 if (edge1 == -1)
515 {
517 << "Did not find two edges sharing vertex " << vertI
518 << " on face " << facei << " vertices:" << mesh_.faces()[facei]
519 << abort(FatalError);
520 }
521
522 return isFeaturePoint(edge0, edge1);
523}
524
525
526// ************************************************************************* //
Various functions to operate on Lists.
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
label fcIndex(const label i) const noexcept
Definition: UListI.H:60
Cell analysis class.
Definition: cellFeatures.H:64
bool isFeaturePoint(const label edge0, const label edge1) const
Are two edges connected at feature point?
Definition: cellFeatures.C:407
bool isFeatureVertex(const label facei, const label vertI) const
Is vertexI on facei used by two edges that form feature.
Definition: cellFeatures.C:471
~cellFeatures()
Destructor.
Definition: cellFeatures.C:399
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
label end() const
Return end (last/second) vertex label.
Definition: edgeI.H:106
vector unitVec(const UList< point > &pts) const
Return the unit vector (end - start)
Definition: edgeI.H:432
label start() const
Return start (first) vertex label.
Definition: edgeI.H:95
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const labelListList & cellEdges() const
virtual const faceList & faces() const =0
Return faces.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
bool faceOnCell(const primitiveMesh &mesh, const label celli, const label facei)
Is face used by cell.
Definition: meshTools.C:330
bool edgeOnFace(const primitiveMesh &mesh, const label facei, const label edgeI)
Is edge used by face.
Definition: meshTools.C:319
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
List< label > labelList
A List of labels.
Definition: List.H:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
void deleteDemandDrivenData(DataPtr &dataPtr)
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333