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 | Copyright (C) 2011 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
27 
28 // Edge to the right of face vertex i
29 inline Foam::label Foam::face::right(const label i) const
30 {
31  return i;
32 }
33 
34 
35 // Edge to the left of face vertex i
36 inline Foam::label Foam::face::left(const label i) const
37 {
38  return rcIndex(i);
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 {}
46 
47 
49 :
50  labelList(s, -1)
51 {}
52 
53 
54 inline Foam::face::face(const labelUList& lst)
55 :
56  labelList(lst)
57 {}
58 
59 
60 inline Foam::face::face(const labelList& lst)
61 :
62  labelList(lst)
63 {}
64 
65 
67 :
68  labelList(lst)
69 {}
70 
71 
73 {
74  is >> *this;
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
80 inline Foam::pointField Foam::face::points(const pointField& meshPoints) const
81 {
82  // There are as many points as there labels for them
83  pointField p(size());
84 
85  // For each point in list, set it to the point in 'pnts' addressed
86  // by 'labs'
87  forAll(p, i)
88  {
89  p[i] = meshPoints[operator[](i)];
90  }
91 
92  // Return list
93  return p;
94 }
95 
96 
97 inline Foam::scalar Foam::face::mag(const pointField& p) const
98 {
100 }
101 
102 
104 {
105  // for a closed polygon a number of edges is the same as number of points
106  return size();
107 }
108 
109 
111 {
112  return edge(operator[](n), operator[](fcIndex(n)));
113 }
114 
115 
116 // Next vertex on face
118 {
119  return operator[](fcIndex(i));
120 }
121 
122 
123 // Previous vertex on face
125 {
126  return operator[](rcIndex(i));
127 }
128 
129 // Number of triangles directly known from number of vertices
131 {
132  return size() - 2;
133 }
134 
135 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
136 
137 inline bool Foam::operator==(const face& a, const face& b)
138 {
139  return face::compare(a,b) != 0;
140 }
141 
142 
143 inline bool Foam::operator!=(const face& a, const face& b)
144 {
145  return face::compare(a,b) == 0;
146 }
147 
148 
149 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
150 
152 {
154  {
155  // Read starting (
156  is.readBegin("face");
157 
158  // Read the 'name' token for the face
159  token t(is);
160 
161  // Read labels
162  is >> static_cast<labelList&>(f);
163 
164  // Read end)
165  is.readEnd("face");
166  }
167  else
168  {
169  is >> static_cast<labelList&>(f);
170  }
171 
172  // Check state of Ostream
173  is.check("Istream& operator>>(Istream&, face&)");
174 
175  return is;
176 }
177 
178 // ************************************************************************* //
Istream & readBegin(const char *funcName)
Definition: Istream.C:86
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static const versionNumber originalVersion
Original version number.
Definition: IOstream.H:203
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
face()
Construct null.
Definition: faceI.H:44
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
label nEdges() const
Return number of edges.
Definition: faceI.H:103
A token holds items read from Istream.
Definition: token.H:69
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
versionNumber version() const
Return the stream version.
Definition: IOstream.H:399
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:124
vector normal(const pointField &) const
Vector normal; magnitude is equal to area of face.
Definition: face.C:552
Istream & readEnd(const char *funcName)
Definition: Istream.C:103
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Istream & operator>>(Istream &, directionInfo &)
static int compare(const face &, const face &)
Compare faces.
Definition: face.C:298
scalar mag(const pointField &) const
Magnitude of face area.
Definition: faceI.H:97
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: faceI.H:80
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
labelList f(nPoints)
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
label nTriangles() const
Number of triangles after splitting.
Definition: faceI.H:130
volScalarField & p
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:110
bool operator!=(const particle &, const particle &)
Definition: particle.C:1106
label size() const
Return the number of elements in the UList.
Definition: ListI.H:170