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-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// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30
31inline Foam::face::face(const label sz)
32:
33 labelList(sz, -1)
34{}
35
36
37inline Foam::face::face(const labelUList& list)
38:
39 labelList(list)
40{}
41
42
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>
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
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
105{
106 const vector n(areaNormal(p));
107 const scalar s(Foam::mag(n));
108 return s < ROOTVSMALL ? Zero : n/s;
109}
110
111
112inline Foam::scalar Foam::face::mag(const UList<point>& p) const
113{
114 return ::Foam::mag(areaNormal(p));
115}
116
117
118inline Foam::label Foam::face::nEdges() const noexcept
119{
120 // for a closed polygon a number of edges is the same as number of points
121 return size();
122}
123
124
125inline Foam::edge Foam::face::edge(const label edgei) const
126{
127 return Foam::edge(thisLabel(edgei), nextLabel(edgei));
128}
129
130
132(
133 const label edgei,
134 const UList<point>& pts
135) const
136{
137 return vector(pts[nextLabel(edgei)] - pts[thisLabel(edgei)]);
138}
139
140
141inline Foam::edge Foam::face::rcEdge(const label edgei) const
142{
143 // Edge 0 (forward and reverse) always starts at [0]
144 // for consistency with face flipping
145 const label pointi = edgei ? (nEdges() - edgei) : 0;
146 return Foam::edge(thisLabel(pointi), prevLabel(pointi));
147}
148
149
151(
152 const label edgei,
153 const UList<point>& pts
154) const
155{
156 // Edge 0 (forward and reverse) always starts at [0]
157 // for consistency with face flipping
158 const label pointi = edgei ? (nEdges() - edgei) : 0;
159 return vector(pts[prevLabel(pointi)] - pts[thisLabel(pointi)]);
160}
161
162
163inline Foam::label Foam::face::which(const label pointLabel) const
164{
165 return labelList::find(pointLabel);
166}
167
168
169inline Foam::label Foam::face::thisLabel(const label i) const
170{
171 return labelList::operator[](i);
172}
173
174
175inline Foam::label Foam::face::nextLabel(const label i) const
176{
177 return labelList::fcValue(i);
178}
179
180
181inline Foam::label Foam::face::prevLabel(const label i) const
182{
183 return labelList::rcValue(i);
184}
185
186
187inline Foam::label Foam::face::nTriangles() const
188{
189 return size() - 2;
190}
191
192
193// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
194
195inline bool Foam::operator==(const face& a, const face& b)
196{
197 return face::compare(a,b) != 0;
198}
199
200inline bool Foam::operator!=(const face& a, const face& b)
201{
202 return face::compare(a,b) == 0;
203}
204
205
206// ************************************************************************* //
label n
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
const T & rcValue(const label i) const
Return reverse circular value (ie, previous value in the list)
const T & fcValue(const label i) const
Return forward circular value (ie, next value in the list)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
iterator find(const word &key)
Find a list element has a name matching key.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
label operator[](const label i) const
Processor-local element id from linear-list of addresses.
Definition: ensightPart.H:189
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
label nEdges() const noexcept
Return number of edges.
Definition: faceI.H:118
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:281
Foam::edge rcEdge(const label edgei) const
Return i-th face edge in reverse walk order.
Definition: faceI.H:141
label which(const label pointLabel) const
Find local index on face for the point label, same as find()
Definition: faceI.H:163
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:175
label thisLabel(const label i) const
Definition: faceI.H:169
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:187
constexpr face() noexcept=default
Default construct.
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:181
Computes the magnitude of an input field.
Definition: mag.H:153
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:119
volScalarField & p
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:239
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const direction noexcept
Definition: Scalar.H:223
Vector< scalar > vector
Definition: vector.H:61
volScalarField & b
Definition: createFields.H:27