Loading...
Searching...
No Matches
pairPatchAgglomeration.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) 2016-2023 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
30#include "meshTools.H"
31#include "edgeHashes.H"
33#include "OBJstream.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
40}
41
42
43// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44
45void Foam::pairPatchAgglomeration::compactLevels(const label nCreatedLevels)
46{
47 nFaces_.setSize(nCreatedLevels);
48 restrictAddressing_.setSize(nCreatedLevels);
49 patchLevels_.setSize(nCreatedLevels);
50}
51
52
53bool Foam::pairPatchAgglomeration::continueAgglomerating
54(
55 const label nLocal,
56 const label nLocalOld,
57 const label nMarkedEdges
58)
59{
60 // Keep agglomerating
61 // - if global number of faces is still changing
62 // - and if local number of faces still too large (on any processor)
63 // or if global number of faces still too large
64
65 label nGlobal = returnReduce(nLocal, sumOp<label>());
66 label nGlobalOld = returnReduce(nLocalOld, sumOp<label>());
67 label nGlobalMarked = returnReduce(nMarkedEdges, sumOp<label>());
68
69 return
70 (
71 returnReduceOr(nLocal > nFacesInCoarsestLevel_)
72 || nGlobal > nGlobalFacesInCoarsestLevel_
73 )
74 && (nGlobal != nGlobalOld || nGlobalMarked > 0);
75}
76
77
78void Foam::pairPatchAgglomeration::setLevel0EdgeWeights()
79{
80 const bPatch& coarsePatch = patchLevels_[0];
81 const auto& coarseEdges = coarsePatch.edges();
82
83 // Statistics on edges
84 label nNonManif = 0;
85 label nFeat = 0;
86
87 for (label i = 0; i < coarsePatch.nInternalEdges(); i++)
88 {
89 scalar edgeLength = coarseEdges[i].mag(coarsePatch.localPoints());
90
91 const labelList& eFaces = coarsePatch.edgeFaces()[i];
92
93 if (eFaces.size() == 2)
94 {
95 scalar cosI =
96 coarsePatch.faceNormals()[eFaces[0]]
97 & coarsePatch.faceNormals()[eFaces[1]];
98
99 const edge edgeCommon = edge(eFaces[0], eFaces[1]);
100
101 if (facePairWeight_.found(edgeCommon))
102 {
103 facePairWeight_[edgeCommon] += edgeLength;
104 }
105 else
106 {
107 facePairWeight_.insert(edgeCommon, edgeLength);
108 }
109
110 if (cosI < Foam::cos(degToRad(featureAngle_)))
111 {
112 facePairWeight_[edgeCommon] = -1.0;
113 nFeat++;
114 }
115 }
116 else
117 {
118 forAll(eFaces, j)
119 {
120 for (label k = j+1; k<eFaces.size(); k++)
121 {
122 facePairWeight_.insert
123 (
124 edge(eFaces[j], eFaces[k]),
125 -1.0
126 );
127 }
128 }
129 nNonManif++;
130 }
131 }
132
133 if (debug)
134 {
135 Pout<< "Level:" << 0
136 << " nEdges:" << coarsePatch.nEdges() << " of which:" << nl
137 << " boundary:" << coarsePatch.nBoundaryEdges() << nl
138 << " non-manifold:" << nNonManif << nl
139 << " feature (angle < " << featureAngle_ << "):" << nFeat << nl
140 << endl;
141 }
142}
143
144
145void Foam::pairPatchAgglomeration::setEdgeWeights
146(
147 const label fineLevelIndex
148)
149{
150 const bPatch& coarsePatch = patchLevels_[fineLevelIndex];
151 const auto& coarseEdges = coarsePatch.edges();
152 const labelList& fineToCoarse = restrictAddressing_[fineLevelIndex];
153
154 edgeHashSet fineFeaturedFaces(coarsePatch.nEdges()/10);
155
156 // Map fine faces with featured edge into coarse faces
157 forAllConstIters(facePairWeight_, iter)
158 {
159 if (iter() == -1.0)
160 {
161 const edge e = iter.key();
162 const edge edgeFeatured
163 (
164 fineToCoarse[e[0]],
165 fineToCoarse[e[1]]
166 );
167 fineFeaturedFaces.insert(edgeFeatured);
168 }
169 }
170
171
172 // Statistics on edges
173 label nNonManif = 0;
174 label nFeat = 0;
175
176 // Clean old weights
177 facePairWeight_.clear();
178 facePairWeight_.reserve(coarsePatch.nEdges());
179
180 for (label i = 0; i < coarsePatch.nInternalEdges(); i++)
181 {
182 scalar edgeLength = coarseEdges[i].mag(coarsePatch.localPoints());
183
184 const labelList& eFaces = coarsePatch.edgeFaces()[i];
185
186 if (eFaces.size() == 2)
187 {
188 const edge edgeCommon(eFaces[0], eFaces[1]);
189 // If the fine 'pair' faces was featured edge so it is
190 // the coarse 'pair'
191 if (fineFeaturedFaces.found(edgeCommon))
192 {
193 auto w = facePairWeight_.find(edgeCommon);
194 if (!w.good())
195 {
196 facePairWeight_.insert(edgeCommon, -1.0);
197 nFeat++;
198 }
199 else if (w() != -1.0)
200 {
201 // Mark as feature edge
202 w() = -1.0;
203 nFeat++;
204 }
205 }
206 else
207 {
208 auto w = facePairWeight_.find(edgeCommon);
209 if (w)
210 {
211 if (w() != -1.0)
212 {
213 w() += edgeLength;
214 }
215 }
216 else
217 {
218 facePairWeight_.insert(edgeCommon, edgeLength);
219 }
220 }
221 }
222 else
223 {
224 // Set edge as barrier by setting weight to -1
225 forAll(eFaces, j)
226 {
227 for (label k = j+1; k<eFaces.size(); k++)
228 {
229 facePairWeight_.insert
230 (
231 edge(eFaces[j], eFaces[k]),
232 -1.0
233 );
234 }
235 }
236 nNonManif++;
237 }
238 }
239
240 if (debug)
241 {
242 Pout<< "Level:" << fineLevelIndex
243 << " nEdges:" << coarsePatch.nEdges() << " of which:" << nl
244 << " boundary:" << coarsePatch.nBoundaryEdges() << nl
245 << " non-manifold:" << nNonManif << nl
246 << " feature (angle < " << featureAngle_ << "):" << nFeat << nl
247 << endl;
248 }
249}
250
251
252bool Foam::pairPatchAgglomeration::isSingleEdgeLoop
253(
254 const bPatch& patch,
255 const labelList& faceIDs,
256 const label facei
257) const
258{
259 // Does combining facei with faceIDs produce a valid single face?
260
261 labelList allFaces(faceIDs.size()+1);
262 SubList<label>(allFaces, faceIDs.size()) = faceIDs;
263 allFaces.last() = facei;
264
265 // Construct single face
266 const indirectPrimitivePatch upp
267 (
268 IndirectList<face>(patch, allFaces),
269 patch.points()
270 );
271
272 return (upp.edgeLoops().size() == 1);
273}
274
275
276Foam::label Foam::pairPatchAgglomeration::maxValidNeighbour
277(
278 const bool addToCluster,
279 const bPatch& patch,
280 const label facei,
281 const labelList& fineToCoarse
282 //const labelListList& coarseToFine
283) const
284{
285 // Return index of neighbour face with max edge weight. Either looks
286 // at clustered faces (addToCluster=true) or at unclustered faces
287
288 const auto& fFaces = patch.faceFaces()[facei];
289
290 label matchFaceNeibNo = -1;
291 scalar maxFaceWeight = -0.5; // negative but larger than -1 (= feature)
292
293 if (addToCluster)
294 {
295 // Check faces to find grouped neighbour with largest face weight
296 // and that forms a single edge cut
297 for (const label faceNeig : fFaces)
298 {
299 const label coarsei = fineToCoarse[faceNeig];
300
301 if (coarsei >= 0)
302 {
303 const edge edgeCommon = edge(facei, faceNeig);
304 const auto& weight = facePairWeight_[edgeCommon];
305 if
306 (
307 (weight > maxFaceWeight)
308 //&& (isSingleEdgeLoop(patch, coarseToFine[coarsei], facei))
309 )
310 {
311 maxFaceWeight = weight;
312 matchFaceNeibNo = faceNeig;
313 }
314 }
315 }
316 }
317 else
318 {
319 // Check faces to find ungrouped neighbour with largest face weight
320 for (const label faceNeig : fFaces)
321 {
322 const label coarsei = fineToCoarse[faceNeig];
323
324 if (coarsei < 0) // ungrouped
325 {
326 const edge edgeCommon = edge(facei, faceNeig);
327 const auto& weight = facePairWeight_[edgeCommon];
328 if (weight > maxFaceWeight)
329 {
330 maxFaceWeight = weight;
331 matchFaceNeibNo = faceNeig;
332 }
333 }
334 }
335 }
337 return matchFaceNeibNo;
338}
339
340
341// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
342
343Foam::pairPatchAgglomeration::pairPatchAgglomeration
344(
345 const faceList& faces,
346 const pointField& points,
348)
349:
350 mergeLevels_
351 (
352 controlDict.getOrDefault<label>("mergeLevels", 2)
353 ),
354 maxLevels_(50),
355 nFacesInCoarsestLevel_
356 (
357 controlDict.get<label>("nFacesInCoarsestLevel")
358 ),
359 nGlobalFacesInCoarsestLevel_(labelMax),
360 //(
361 // controlDict.get<label>("nGlobalFacesInCoarsestLevel")
362 //),
363 featureAngle_
364 (
365 controlDict.getOrDefault<scalar>("featureAngle", 0)
366 ),
367 nFaces_(maxLevels_),
368 restrictAddressing_(maxLevels_),
369 restrictTopBottomAddressing_(identity(faces.size())),
370 patchLevels_(maxLevels_),
371 facePairWeight_(faces.size())
372{
373 // Set base fine patch
374 patchLevels_.set(0, new bPatch(faces, points));
375
376 // Set number of faces for the base patch
377 nFaces_[0] = faces.size();
378
379 // Set edge weights for level 0
380 setLevel0EdgeWeights();
381}
382
383
384Foam::pairPatchAgglomeration::pairPatchAgglomeration
385(
386 const faceList& faces,
387 const pointField& points,
388 const label mergeLevels,
389 const label maxLevels,
390 const label nFacesInCoarsestLevel, // local number of cells
391 const label nGlobalFacesInCoarsestLevel, // global number of cells
392 const scalar featureAngle
393)
394:
395 mergeLevels_(mergeLevels),
396 maxLevels_(maxLevels),
397 nFacesInCoarsestLevel_(nFacesInCoarsestLevel),
398 nGlobalFacesInCoarsestLevel_(nGlobalFacesInCoarsestLevel),
399 featureAngle_(featureAngle),
400 nFaces_(maxLevels_),
401 restrictAddressing_(maxLevels_),
402 restrictTopBottomAddressing_(identity(faces.size())),
403 patchLevels_(maxLevels_),
404 facePairWeight_(faces.size())
405{
406 // Set base fine patch
407 patchLevels_.set(0, new bPatch(faces, points));
408
409 // Set number of faces for the base patch
410 nFaces_[0] = faces.size();
411
412 // Set edge weights for level 0
413 setLevel0EdgeWeights();
414}
415
416
417// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
418
420{}
421
422
423// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
424
427(
428 const label i
429) const
430{
431 return patchLevels_[i];
432}
433
434
435void Foam::pairPatchAgglomeration::mapBaseToTopAgglom
436(
437 const label fineLevelIndex
438)
439{
440 const labelList& fineToCoarse = restrictAddressing_[fineLevelIndex];
441 forAll(restrictTopBottomAddressing_, i)
442 {
443 restrictTopBottomAddressing_[i] =
444 fineToCoarse[restrictTopBottomAddressing_[i]];
445 }
446}
447
448
449bool Foam::pairPatchAgglomeration::agglomeratePatch
450(
451 const bPatch& patch,
452 const labelList& fineToCoarse,
453 const label fineLevelIndex,
454 label& nMarkedEdges
455)
456{
457 if (min(fineToCoarse) == -1)
458 {
460 << "min(fineToCoarse) == -1" << exit(FatalError);
461 }
462
463 if (fineToCoarse.size() == 0)
464 {
465 return false;
466 }
467
468 if (fineToCoarse.size() != patch.size())
469 {
471 << "restrict map does not correspond to fine level. " << endl
472 << " Sizes: restrictMap: " << fineToCoarse.size()
473 << " nEqns: " << patch.size()
474 << abort(FatalError);
475 }
476
477 const label nCoarseI = max(fineToCoarse) + 1;
478 List<face> patchFaces(nCoarseI);
479
480 // Patch faces per agglomeration
481 labelListList coarseToFine(invertOneToMany(nCoarseI, fineToCoarse));
482
483 // Additional feature edges created
484 nMarkedEdges = 0;
485
486
487 for (label coarseI = 0; coarseI < nCoarseI; coarseI++)
488 {
489 const labelList& fineFaces = coarseToFine[coarseI];
490
491 // Construct single face
493 (
494 IndirectList<face>(patch, fineFaces),
495 patch.points()
496 );
497
498 if (upp.edgeLoops().size() != 1)
499 {
500 // More than one outside loop. Possible because pair-wise
501 // agglomeration has e.g. walked a path leaving a hole in the
502 // middle of a coarse face.
503 // - mark as feature edges to avoid locally further agglomeration
504 // - but further agglomeration might e.g. 'fill in the hole'
505 // - or immediately leave all agglomeration
506 // Currently we have no choice but to leave agglomeration since
507 // we cannot store a face-with-hole.
508
509
510 //{
511 // OBJstream os
512 // (
513 // "error_agglomeration_"+Foam::name(fineLevelIndex)+".obj"
514 // );
515 // Pout<< "Writing error patch at level:" << fineLevelIndex
516 // << " to:" << os.name() << endl;
517 // os.write(upp.localFaces(), upp.localPoints(), true);
518 //}
519
520 if (fineFaces.size() >= 2)
521 {
522 forAll(fineFaces, j)
523 {
524 for (label k = j+1; k<fineFaces.size(); k++)
525 {
526 const edge e(fineFaces[j], fineFaces[k]);
527
528 auto w = facePairWeight_.find(e);
529 if (!w.good())
530 {
531 facePairWeight_.insert(e, -1.0);
532 nMarkedEdges++;
533 }
534 else if (w() != -1.0)
535 {
536 // Mark as feature edge
537 w() = -1.0;
538 nMarkedEdges++;
539 }
540 }
541 }
542 }
543
544 return false;
545 }
546
547 // In-place override face
548 patchFaces[coarseI] = face
549 (
551 (
552 upp.meshPoints(),
553 upp.edgeLoops()[0]
554 )
555 );
556 }
557
558 patchLevels_.set
559 (
560 fineLevelIndex,
561 new bPatch
562 (
563 SubList<face>(patchFaces, nCoarseI, 0),
564 patch.points()
566 );
567
568 return true;
569}
570
571
573{
574 label nPairLevels = 0;
575 label nCreatedLevels = 1; // 0 level is the base patch
576
577 label nCoarseFaces = 0;
578 label nCoarseFacesOld = 0;
579 label nMarkedEdges = 0;
580
581 while (nCreatedLevels < maxLevels_)
582 {
583 const bPatch& patch = patchLevels_[nCreatedLevels - 1];
584
585 // Agglomerate locally
586 tmp<labelField> tfinalAgglom;
587
588 bool createdLevel = false;
589 while (!createdLevel)
590 {
591 // Agglomerate locally using edge weights
592 // - calculates nCoarseFaces; returns fine to coarse addressing
593 tfinalAgglom = agglomerateOneLevel(nCoarseFaces, patch);
594
595 if (nCoarseFaces == 0)
596 {
597 break;
598 }
599 else
600 {
601 // Attempt to create coarse face addressing
602 // - returns true if successful; otherwise resets edge weights
603 // and tries again...
604 createdLevel = agglomeratePatch
605 (
606 patch,
607 tfinalAgglom,
608 nCreatedLevels,
609 nMarkedEdges
610 );
611 }
612 }
613
614 if (createdLevel)
615 {
616 if (debug)
617 {
618 const auto& agglomPatch = patchLevels_[nCreatedLevels];
619 OBJstream os("agglomPatch"+Foam::name(nCreatedLevels)+".obj");
620 Pout<< "Writing new patch at level:" << nCreatedLevels
621 << " to:" << os.name() << endl;
622 os.write(agglomPatch, agglomPatch.points(), true);
623 }
624
625 restrictAddressing_.set(nCreatedLevels, tfinalAgglom);
626
627 mapBaseToTopAgglom(nCreatedLevels);
628
629 setEdgeWeights(nCreatedLevels);
630
631 if (nPairLevels % mergeLevels_)
632 {
633 combineLevels(nCreatedLevels);
634 }
635 else
636 {
637 nCreatedLevels++;
638 }
639
640 nPairLevels++;
641
642 nFaces_[nCreatedLevels] = nCoarseFaces;
643 }
644
645 // Check to see if we need to continue agglomerating
646 // - Note: performs parallel reductions
647 if (!continueAgglomerating(nCoarseFaces, nCoarseFacesOld, nMarkedEdges))
648 {
649 break;
650 }
651
652 nCoarseFacesOld = nCoarseFaces;
653 }
654
655 compactLevels(nCreatedLevels);
656}
657
658
659Foam::tmp<Foam::labelField> Foam::pairPatchAgglomeration::agglomerateOneLevel
660(
661 label& nCoarseFaces,
662 const bPatch& patch
663)
664{
665 const label nFineFaces = patch.size();
666
667 auto tcoarseCellMap = tmp<labelField>::New(nFineFaces, -1);
668 auto& coarseCellMap = tcoarseCellMap.ref();
669
670 const labelListList& faceFaces = patch.faceFaces();
671
672 //labelListList coarseToFine(nFineFaces);
673 nCoarseFaces = 0;
674
675 forAll(faceFaces, facei)
676 {
677 if (coarseCellMap[facei] < 0)
678 {
679 const label matchFaceNeibNo = maxValidNeighbour
680 (
681 false, // ungrouped neighbours only
682 patch,
683 facei,
684 coarseCellMap
685 //coarseToFine
686 );
687
688 if (matchFaceNeibNo >= 0)
689 {
690 // Make a new group
691 coarseCellMap[facei] = nCoarseFaces;
692 coarseCellMap[matchFaceNeibNo] = nCoarseFaces;
693 //coarseToFine[nCoarseFaces] =
694 // labelList({facei, matchFaceNeibNo});
695 nCoarseFaces++;
696 }
697 else
698 {
699 // No match. Find the best neighbouring cluster and
700 // put the cell there
701 const label clusterMatchFaceNo = maxValidNeighbour
702 (
703 true, // grouped neighbours only
704 patch,
705 facei,
706 coarseCellMap
707 //coarseToFine
708 );
709
710 if (clusterMatchFaceNo >= 0)
711 {
712 // Add the cell to the best cluster
713 const label coarsei = coarseCellMap[clusterMatchFaceNo];
714 coarseCellMap[facei] = coarsei;
715 //coarseToFine[coarsei].append(facei);
716 }
717 else
718 {
719 // If not create single-cell "clusters" for each
720 coarseCellMap[facei] = nCoarseFaces;
721 //coarseToFine[nCoarseFaces] = labelList({facei});
722 nCoarseFaces++;
723 }
724 }
725 }
726 }
727
728 //coarseToFine.setSize(nCoarseFaces);
729
730 // Check that all faces are part of clusters,
731 for (label facei=0; facei<nFineFaces; facei++)
732 {
733 if (coarseCellMap[facei] < 0)
734 {
736 << " face " << facei
737 << " is not part of a cluster"
738 << exit(FatalError);
739 }
740 }
741
742 return tcoarseCellMap;
743}
744
745
746void Foam::pairPatchAgglomeration::combineLevels(const label curLevel)
747{
748 label prevLevel = curLevel - 1;
749
750 // Set the previous level nCells to the current
751 nFaces_[prevLevel] = nFaces_[curLevel];
752
753 // Map the restrictAddressing from the coarser level into the previous
754 // finer level
755
756 const labelList& curResAddr = restrictAddressing_[curLevel];
757 labelList& prevResAddr = restrictAddressing_[prevLevel];
758
759 forAll(prevResAddr, i)
760 {
761 prevResAddr[i] = curResAddr[prevResAddr[i]];
762 }
763
764 // Delete the restrictAddressing for the coarser level
765 restrictAddressing_.set(curLevel, nullptr);
766
767 patchLevels_.set(prevLevel, patchLevels_.set(curLevel, nullptr));
768}
769
770
771// ************************************************************************* //
label k
A List with indirect addressing.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
T & last()
Access last element of the list, position [size()-1].
Definition UList.H:971
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Primitive patch pair agglomerate method.
PtrList< bPatch > patchLevels_
Hierarchy of patch addressing.
labelList nFaces_
The number of faces in each level.
const bPatch & patchLevel(const label leveli) const
Return primitivePatch of given level.
void agglomerate()
Agglomerate patch.
EdgeMap< scalar > facePairWeight_
Edge weights.
label mergeLevels_
Number of levels to merge, 1 = don't merge, 2 = merge pairs etc.
label maxLevels_
Max number of levels.
label nFacesInCoarsestLevel_
Number of faces in coarsest level.
PtrList< labelField > restrictAddressing_
Cell restriction addressing array.
PrimitivePatch< List< face >, const pointField > bPatch
label nGlobalFacesInCoarsestLevel_
Global number of faces in coarsest level.
labelList restrictTopBottomAddressing_
Maps from finest to coarsest.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
runTime controlDict().readEntry("adjustTimeStep"
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
HashSet< edge, Hash< edge > > edgeHashSet
A HashSet with edge for its key. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:48
List< label > labelList
A List of labels.
Definition List.H:62
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
List< face > faceList
List of faces.
Definition faceListFwd.H:41
constexpr label labelMax
Definition label.H:55
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
PrimitivePatch<::Foam::List< face >, const pointField > bPatch
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition ListOps.C:125
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensionedScalar cos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Unit conversion functions.