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-2023 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 inline Foam::face::face(const label sz)
32 :
33  labelList(sz, -1)
34 {}
35 
36 
37 inline Foam::face::face(const labelUList& list)
38 :
39  labelList(list)
40 {}
41 
42 
43 inline Foam::face::face(labelList&& list)
44 :
45  labelList(std::move(list))
46 {}
47 
48 
49 inline Foam::face::face(std::initializer_list<label> list)
50 :
51  labelList(list)
52 {}
53 
54 
55 template<unsigned N>
56 inline Foam::face::face(const FixedList<label, N>& list)
57 :
58  labelList(list)
59 {}
60 
61 
62 inline Foam::face::face(const labelUList& list, const labelUList& indices)
63 :
64  labelList(list, indices)
65 {}
66 
67 
68 template<unsigned N>
69 inline Foam::face::face
70 (
71  const labelUList& list,
72  const FixedList<label, N>& indices
73 )
74 :
75  labelList(list, indices)
76 {}
77 
78 
79 inline Foam::face::face(Istream& is)
80 :
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 
104 inline Foam::vector Foam::face::unitNormal(const UList<point>& p) const
105 {
106  vector n(areaNormal(p));
107  (void) n.normalise(ROOTVSMALL);
108  return n;
109 }
110 
111 
112 inline Foam::scalar Foam::face::mag(const UList<point>& p) const
113 {
114  return ::Foam::mag(areaNormal(p));
115 }
116 
117 
119 Foam::face::box(const UList<point>& pts) const
120 {
121  Pair<point> bb(point::rootMax, point::rootMin);
122 
123  for (const label pointi : *this)
124  {
125  bb.first() = min(bb.first(), pts[pointi]);
126  bb.second() = max(bb.second(), pts[pointi]);
127  }
128 
129  return bb;
130 }
131 
132 
133 inline Foam::label Foam::face::nEdges() const noexcept
134 {
135  // for a closed polygon a number of edges is the same as number of points
136  return size();
137 }
138 
139 
140 inline Foam::edge Foam::face::edge(const label edgei) const
141 {
142  return Foam::edge(thisLabel(edgei), nextLabel(edgei));
143 }
144 
145 
147 (
148  const label edgei,
150 ) const
151 {
152  return vector(pts[nextLabel(edgei)] - pts[thisLabel(edgei)]);
153 }
154 
155 
156 inline Foam::edge Foam::face::rcEdge(const label edgei) const
157 {
158  // Edge 0 (forward and reverse) always starts at [0]
159  // for consistency with face flipping
160  const label pointi = edgei ? (nEdges() - edgei) : 0;
161  return Foam::edge(thisLabel(pointi), prevLabel(pointi));
162 }
163 
164 
166 (
167  const label edgei,
168  const UList<point>& pts
169 ) const
170 {
171  // Edge 0 (forward and reverse) always starts at [0]
172  // for consistency with face flipping
173  const label pointi = edgei ? (nEdges() - edgei) : 0;
174  return vector(pts[prevLabel(pointi)] - pts[thisLabel(pointi)]);
175 }
176 
178 inline Foam::label Foam::face::which(const label vertex) const
179 {
180  return labelList::find(vertex);
181 }
182 
184 inline Foam::label Foam::face::thisLabel(const label i) const
185 {
186  return labelList::operator[](i);
187 }
188 
190 inline Foam::label Foam::face::nextLabel(const label i) const
191 {
192  return labelList::fcValue(i);
193 }
194 
196 inline Foam::label Foam::face::prevLabel(const label i) const
197 {
198  return labelList::rcValue(i);
199 }
200 
202 inline Foam::label Foam::face::nTriangles() const
203 {
204  return size() - 2;
205 }
206 
207 
208 inline bool Foam::face::connected(const labelUList& other) const
209 {
210  for (const label pointi : *this)
211  {
212  if (other.contains(pointi))
213  {
214  return true;
215  }
216  }
217  return false;
218 }
219 
220 
221 template<unsigned N>
222 inline bool Foam::face::connected(const FixedList<label, N>& other) const
223 {
224  for (const label pointi : *this)
225  {
226  if (other.contains(pointi))
227  {
228  return true;
229  }
230  }
231  return false;
232 }
233 
234 
235 inline bool Foam::face::contains(const Foam::edge& e) const
236 {
237  // or (find(e) >= 0)
238  return (edgeDirection(e) != 0);
239 }
240 
241 
242 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
243 
244 inline void Foam::face::operator+=(const label vertexOffset)
245 {
246  if (vertexOffset)
247  {
248  for (label& vrt : static_cast<labelList&>(*this))
249  {
250  vrt += vertexOffset;
251  }
252  }
253 }
254 
255 
256 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
257 
258 inline bool Foam::operator==(const face& a, const face& b)
259 {
260  return face::compare(a,b) != 0;
261 }
262 
263 inline bool Foam::operator!=(const face& a, const face& b)
264 {
265  return face::compare(a,b) == 0;
266 }
267 
268 
269 // ************************************************************************* //
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:105
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:107
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr face() noexcept=default
Default construct.
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:183
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition: faceI.H:149
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition: faceI.H:97
Pair< point > box(const UList< point > &pts) const
The enclosing (bounding) box for the face.
Definition: faceI.H:112
bool contains(const T &val) const
True if the value is contained in the list.
Definition: UListI.H:300
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:189
label & operator[](const label i)
Return element of UList.
Definition: UListI.H:361
pointField points(const UList< point > &pts) const
Return the points corresponding to this face.
Definition: faceI.H:80
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
const label & fcValue(const label i) const
Return forward circular value (ie, next value in the list)
Definition: UListI.H:111
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:59
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
label which(const label vertex) const
Find local vertex on face for the vertex label, same as find()
Definition: faceI.H:171
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:391
label find(const label &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
void operator+=(const label vertexOffset)
Increment (offset) vertices by given amount.
Definition: faceI.H:237
const direction noexcept
Definition: Scalar.H:258
label nEdges() const noexcept
Return number of edges.
Definition: faceI.H:126
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Foam::edge edge(const label edgei) const
Return i-th face edge (forward walk order).
Definition: faceI.H:133
bool connected(const labelUList &other) const
True if the face has at least one vertex in common with other.
Definition: faceI.H:201
bool contains(const Foam::edge &e) const
True if face contains(edge)
Definition: faceI.H:228
label thisLabel(const label i) const
The vertex on face - identical to operator[], but with naming similar to nextLabel(), prevLabel()
Definition: faceI.H:177
label n
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:195
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:273
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
const label & rcValue(const label i) const
Return reverse circular value (ie, previous value in the list)
Definition: UListI.H:125
const pointField & pts