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-2018 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 // Edge to the right of face vertex i
32 inline Foam::label Foam::face::right(const label i) const
33 {
34  return i;
35 }
36 
37 
38 // Edge to the left of face vertex i
39 inline Foam::label Foam::face::left(const label i) const
40 {
41  return rcIndex(i);
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 {}
49 
50 
51 inline Foam::face::face(const label sz)
52 :
53  labelList(sz, -1)
54 {}
55 
56 
57 inline Foam::face::face(const labelUList& list)
58 :
59  labelList(list)
60 {}
61 
62 
63 template<unsigned N>
65 :
66  labelList(list)
67 {}
68 
69 
70 inline Foam::face::face(std::initializer_list<label> list)
71 :
72  labelList(list)
73 {}
74 
75 
77 :
78  labelList(std::move(list))
79 {}
80 
81 
83 {
84  is >> *this;
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 (
92  const UList<point>& meshPoints
93 ) const
94 {
95  // There are as many points as there are labels for them
96  pointField p(size());
97 
98  // For each point in list, set it to the point in 'pnts' addressed
99  // by 'labs'
100  label i = 0;
101  for (const label pointi : *this)
102  {
103  p[i++] = meshPoints[pointi];
104  }
105 
106  // Return list
107  return p;
108 }
109 
110 
112 {
113  const vector n(areaNormal(p));
114  const scalar s(Foam::mag(n));
115  return s < ROOTVSMALL ? Zero : n/s;
116 }
117 
118 
119 inline Foam::scalar Foam::face::mag(const UList<point>& p) const
120 {
121  return ::Foam::mag(areaNormal(p));
122 }
123 
124 
126 {
127  // for a closed polygon a number of edges is the same as number of points
128  return size();
129 }
130 
131 
133 {
134  return edge(operator[](n), operator[](fcIndex(n)));
135 }
136 
137 
138 inline bool Foam::face::found(const label pointLabel) const
139 {
140  return labelList::found(pointLabel);
141 }
142 
143 
144 inline Foam::label Foam::face::which(const label pointLabel) const
145 {
146  return labelList::find(pointLabel);
147 }
148 
149 
151 {
152  return labelList::fcValue(i);
153 }
154 
155 
157 {
158  return labelList::rcValue(i);
159 }
160 
161 
163 {
164  return size() - 2;
165 }
166 
167 
168 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
169 
170 inline bool Foam::operator==(const face& a, const face& b)
171 {
172  return face::compare(a,b) != 0;
173 }
174 
175 inline bool Foam::operator!=(const face& a, const face& b)
176 {
177  return face::compare(a,b) == 0;
178 }
179 
180 
181 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
182 
184 {
186  {
187  // Read starting (
188  is.readBegin("face");
189 
190  // Read the 'name' token for the face
191  token t(is);
192 
193  // Read labels
194  is >> static_cast<labelList&>(f);
195 
196  // Read end)
197  is.readEnd("face");
198  }
199  else
200  {
201  is >> static_cast<labelList&>(f);
202  }
203 
204  is.check(FUNCTION_NAME);
205  return is;
206 }
207 
208 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::face::points
pointField points(const UList< point > &points) const
Return the points corresponding to this face.
Definition: faceI.H:91
s
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))
Definition: gmvOutputSpray.H:25
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::face::mag
scalar mag(const UList< point > &p) const
Magnitude of face area.
Definition: faceI.H:119
Foam::face::nEdges
label nEdges() const
Return number of edges.
Definition: faceI.H:125
Foam::face::faceEdge
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:132
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::face::nTriangles
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:162
Foam::operator>>
Istream & operator>>(Istream &, directionInfo &)
Definition: directionInfo.C:228
Foam::token
A token holds an item read from Istream.
Definition: token.H:69
Foam::Istream::readEnd
bool readEnd(const char *funcName)
End read of data chunk, ends with ')'.
Definition: Istream.C:127
Foam::Istream::readBegin
bool readBegin(const char *funcName)
Begin read of data chunk, starts with '('.
Definition: Istream.C:109
Foam::face::compare
static int compare(const face &a, const face &b)
Compare faces.
Definition: face.C:300
Foam::operator!=
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:235
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::face::unitNormal
vector unitNormal(const UList< point > &p) const
The unit normal.
Definition: faceI.H:111
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::face::which
label which(const label pointLabel) const
Find local index on face for the point label,.
Definition: faceI.H:144
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::IOstreamOption::version
versionNumber version() const noexcept
Get the stream version.
Definition: IOstreamOption.H:321
Foam::IOstreamOption::originalVersion
static const versionNumber originalVersion
The original version number.
Definition: IOstreamOption.H:190
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
Foam::ListOps::find
label find(const ListType &input, const UnaryPredicate &pred, const label start=0)
Find index of the first occurrence that satisfies the predicate.
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::face::prevLabel
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:156
Foam::face::nextLabel
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:150
f
labelList f(nPoints)
Foam::Vector< scalar >
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::UList< label >
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:261
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::face::face
face()
Construct null.
Definition: faceI.H:47
Foam::face::found
bool found(const label pointLabel) const
Return true if the point label is found in face.
Definition: faceI.H:138