PrimitivePatchBdryPoints.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) 2020-2022 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 "PrimitivePatch.H"
30 #include "HashSet.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class FaceList, class PointField>
35 void
37 {
38  if (boundaryPointsPtr_)
39  {
40  // Error to recalculate if already allocated
42  << "boundaryPoints already calculated"
43  << abort(FatalError);
44  }
45 
46  labelHashSet bp(0);
47 
48  if (hasEdges())
49  {
51  << "Calculating boundary points from existing addressing"
52  << nl;
53 
54  bp.reserve(2*nBoundaryEdges());
55 
56  for (const edge& e : boundaryEdges())
57  {
58  bp.insert(e.first());
59  bp.insert(e.second());
60  }
61  }
62  else
63  {
65  << "Calculating boundary points with manual edge addressing"
66  << nl;
67 
68 
69  // Calculate manually.
70  // Needs localFaces, but uses local hashes of the edges here
71  // instead of forcing a full faceFaces/edgeFaces/faceEdges calculation
72 
73  // Get reference to localFaces
74  const List<face_type>& locFcs = localFaces();
75 
76  // Guess the max number of edges/neighbours for a face
77  label edgeCount = 0;
78  for (const auto& f : locFcs)
79  {
80  edgeCount += f.nEdges();
81  }
82 
83  // ie, EdgeMap<label> to keep counts
84  HashTable<label, edge, Hash<edge>> knownEdges(2*edgeCount);
85 
86  for (const auto& f : locFcs)
87  {
88  const label numEdges = f.nEdges();
89 
90  for (label edgei = 0; edgei < numEdges; ++edgei)
91  {
92  ++ knownEdges(f.edge(edgei));
93  }
94  }
95 
96  edgeCount = 0;
97 
98  forAllConstIters(knownEdges, iter)
99  {
100  if (1 == iter.val()) // Singly connected edge
101  {
102  ++edgeCount;
103  }
104  }
105 
106  bp.reserve(2*edgeCount);
107 
108  forAllConstIters(knownEdges, iter)
109  {
110  const edge& e = iter.key();
111 
112  if (1 == iter.val()) // Singly connected edge
113  {
114  bp.insert(e.first());
115  bp.insert(e.second());
116  }
117  }
118  }
119 
120  boundaryPointsPtr_.reset(new labelList(bp.sortedToc()));
121  DebugInfo << " Finished." << nl;
122 }
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
127 template<class FaceList, class PointField>
130 {
132  (
133  point_type::uniform(ROOTVGREAT),
134  point_type::uniform(-ROOTVGREAT)
135  );
136 
137  if (hasMeshPoints())
138  {
139  // Less looping if meshPoints() are already available
140  for (const label pointi : meshPoints())
141  {
142  bb.first() = min(bb.first(), points_[pointi]);
143  bb.second() = max(bb.second(), points_[pointi]);
144  }
145  }
146  else
147  {
148  // Walk the points on each face
149  for (const face_type& f : *this)
150  {
151  for (const label pointi : f)
152  {
153  bb.first() = min(bb.first(), points_[pointi]);
154  bb.second() = max(bb.second(), points_[pointi]);
155  }
156  }
157  }
159  return bb;
160 }
161 
162 
163 template<class FaceList, class PointField>
164 Foam::scalar
166 {
167  scalar radiusSqr = 0;
168 
169  const point_type& fc = this->faceCentres()[facei];
170 
171  for (const label fp : this->operator[](facei))
172  {
173  const scalar sqrDist = magSqr(fc - points_[fp]);
174  if (radiusSqr < sqrDist)
175  {
176  radiusSqr = sqrDist;
177  }
178  }
179 
180  return radiusSqr;
181 }
182 
183 
184 // ************************************************************************* //
const T & first() const noexcept
Access the first element.
Definition: Pair.H:137
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
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Pair< point_type > box() const
The enclosing (bounding) box for the patch points.
A list of faces which address into the list of points.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
#define DebugInFunction
Report an information message using Foam::Info.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
#define DebugInfo
Report an information message using Foam::Info.
labelList f(nPoints)
scalar sphere(const label facei) const
The enclosing (bounding) sphere radius^2 for specified face.
List< label > labelList
A List of labels.
Definition: List.H:62
const T & second() const noexcept
Access the second element.
Definition: Pair.H:147
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28