Loading...
Searching...
No Matches
faceI.H
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 OpenFOAM Foundation
9 Copyright (C) 2017-2025 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// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31inline Foam::face::face(const label sz)
32:
33 labelList(sz, -1)
34{}
35
37inline Foam::face::face(const labelUList& list)
38:
39 labelList(list)
40{}
41
43inline Foam::face::face(labelList&& list)
44:
45 labelList(std::move(list))
46{}
47
48
49inline Foam::face::face(std::initializer_list<label> list)
50:
51 labelList(list)
52{}
53
54
55template<unsigned N>
56inline Foam::face::face(const FixedList<label, N>& list)
57:
58 labelList(list)
59{}
60
61
62inline Foam::face::face(const labelUList& list, const labelUList& indices)
63:
64 labelList(list, indices)
65{}
66
67
68template<unsigned N>
70(
71 const labelUList& list,
72 const FixedList<label, N>& indices
73)
74:
75 labelList(list, indices)
76{}
77
78
79inline Foam::face::face(Istream& is)
81 labelList(is)
82{}
83
84
85// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86
88{
89 // There are as many points as there are labels for them
90 pointField p(size());
91
92 auto iter = p.begin();
93
94 for (const label pointi : *this)
95 {
96 *iter = pts[pointi];
97 ++iter;
98 }
99
100 return p;
101}
102
103
104inline Foam::vector Foam::face::unitNormal(const UList<point>& p) const
107 (void) n.normalise(ROOTVSMALL);
108 return n;
109}
110
112inline Foam::scalar Foam::face::mag(const UList<point>& p) const
113{
114 return ::Foam::mag(areaNormal(p));
115}
116
117
118inline Foam::scalar Foam::face::magSqr(const UList<point>& p) const
119{
120 return ::Foam::magSqr(areaNormal(p));
121}
122
123
125Foam::face::box(const UList<point>& pts) const
126{
128
129 for (const label pointi : *this)
130 {
131 bb.first() = min(bb.first(), pts[pointi]);
132 bb.second() = max(bb.second(), pts[pointi]);
134
135 return bb;
136}
137
138
139inline Foam::point
140Foam::face::average(const UList<point>& points) const
141{
142 // Note: use double precision to avoid overflows when summing
143 doubleVector avg(Zero);
144
145 if (const auto n = size(); n)
146 {
147 for (const auto pointi : *this)
148 {
149 avg += points[pointi];
150 }
151 avg /= n;
152 }
153
154 return avg;
155}
156
157
158inline Foam::label Foam::face::nEdges() const noexcept
159{
160 // for a closed polygon a number of edges is the same as number of points
161 return size();
162}
163
165inline Foam::edge Foam::face::edge(const label edgei) const
166{
167 return Foam::edge(thisLabel(edgei), nextLabel(edgei));
168}
169
170
172(
173 const label edgei,
175) const
176{
177 return vector(pts[nextLabel(edgei)] - pts[thisLabel(edgei)]);
178}
179
180
181inline Foam::edge Foam::face::rcEdge(const label edgei) const
182{
183 // Edge 0 (forward and reverse) always starts at [0]
184 // for consistency with face flipping
185 const label pointi = edgei ? (nEdges() - edgei) : 0;
186 return Foam::edge(thisLabel(pointi), prevLabel(pointi));
187}
188
189
191(
192 const label edgei,
193 const UList<point>& pts
194) const
195{
196 // Edge 0 (forward and reverse) always starts at [0]
197 // for consistency with face flipping
198 const label pointi = edgei ? (nEdges() - edgei) : 0;
199 return vector(pts[prevLabel(pointi)] - pts[thisLabel(pointi)]);
200}
201
203inline Foam::label Foam::face::which(const label vertex) const
204{
205 return labelList::find(vertex);
206}
207
209inline Foam::label Foam::face::thisLabel(const label i) const
210{
211 return labelList::operator[](i);
212}
213
215inline Foam::label Foam::face::nextLabel(const label i) const
216{
217 return labelList::fcValue(i);
218}
219
221inline Foam::label Foam::face::prevLabel(const label i) const
222{
223 return labelList::rcValue(i);
224}
225
227inline Foam::label Foam::face::nTriangles() const noexcept
228{
229 return labelList::size() - 2;
230}
231
232
233inline bool Foam::face::connected(const labelUList& other) const
234{
235 for (const label pointi : *this)
236 {
237 if (other.contains(pointi))
238 {
239 return true;
241 }
242 return false;
243}
244
245
246template<unsigned N>
247inline bool Foam::face::connected(const FixedList<label, N>& other) const
248{
249 for (const label pointi : *this)
250 {
251 if (other.contains(pointi))
252 {
253 return true;
254 }
255 }
256 return false;
257}
258
259
260inline bool Foam::face::contains(const Foam::edge& e) const
261{
262 // or (find(e) >= 0)
263 return (edgeDirection(e) != 0);
264}
265
266
267// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
268
269inline void Foam::face::operator+=(const label vertexOffset)
270{
271 if (vertexOffset)
272 {
273 for (auto& vrt : static_cast<labelList&>(*this))
274 {
275 vrt += vertexOffset;
277 }
278}
279
280
281// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
282
283inline bool Foam::operator==(const face& a, const face& b)
284{
285 return face::compare(a,b) != 0;
286}
287
288inline bool Foam::operator!=(const face& a, const face& b)
289{
290 return face::compare(a,b) == 0;
291}
292
293
294// ************************************************************************* //
label n
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
bool contains(const T &val) const
True if the value is contained in the list.
Definition FixedListI.H:294
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool contains(const T &val) const
True if the value is contained in the list.
Definition UListI.H:302
const label & rcValue(const label i) const
const label & fcValue(const label i) const
label size() const noexcept
Definition UList.H:706
void size(const label n)
Definition UList.H:118
label & operator[](const label i)
label find(const label &val) const
static const Form rootMin
static const Form rootMax
A Vector of values with double precision.
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
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition faceI.H:158
Pair< point > box(const UList< point > &pts) const
The enclosing (bounding) box for the face.
Definition faceI.H:118
bool contains(const Foam::edge &e) const
True if face contains(edge).
Definition faceI.H:253
label nEdges() const noexcept
Return number of edges.
Definition faceI.H:151
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition faceI.H:105
static int compare(const face &a, const face &b)
Compare faces.
Definition face.C:273
bool connected(const labelUList &other) const
True if the face has at least one vertex in common with other.
Definition faceI.H:226
int edgeDirection(const Foam::edge &e) const
Test the edge direction on the face.
Definition face.C:816
void operator+=(const label vertexOffset)
Increment (offset) vertices by given amount.
Definition faceI.H:262
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition faceI.H:174
point average(const UList< point > &points) const
Average of face points: a quick estimate of the face centre.
Definition faceI.H:133
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition faceI.H:80
label which(const label vertex) const
Find local vertex on face for the vertex label, same as find().
Definition faceI.H:196
scalar magSqr(const UList< point > &p) const
Magnitude squared of face area.
Definition faceI.H:111
label nextLabel(const label i) const
Next vertex on face.
Definition faceI.H:208
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition faceI.H:97
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(),...
Definition faceI.H:202
vector areaNormal(const UList< point > &p) const
The area normal - with magnitude equal to area of face.
Definition face.C:569
label nTriangles() const noexcept
Number of triangles after splitting.
Definition faceI.H:220
constexpr face() noexcept=default
Default construct.
label prevLabel(const label i) const
Previous vertex on face.
Definition faceI.H:214
volScalarField & p
const pointField & points
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Vector< double > doubleVector
Definition vector.H:54
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
const pointField & pts
volScalarField & b
volScalarField & e