faceI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
27 
29 {}
30 
31 
33 :
34  labelList(s, -1)
35 {}
36 
37 
38 inline Foam::face::face(const labelUList& lst)
39 :
40  labelList(lst)
41 {}
42 
43 
44 inline Foam::face::face(const labelList& lst)
45 :
46  labelList(lst)
47 {}
48 
49 
51 :
52  labelList(move(lst))
53 {}
54 
55 
57 {
58  is >> *this;
59 }
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
64 inline Foam::pointField Foam::face::points(const pointField& meshPoints) const
65 {
66  // There are as many points as there labels for them
67  pointField p(size());
68 
69  // For each point in list, set it to the point in 'pnts' addressed
70  // by 'labs'
71  forAll(p, i)
72  {
73  p[i] = meshPoints[operator[](i)];
74  }
75 
76  // Return list
77  return p;
78 }
79 
80 
81 inline Foam::scalar Foam::face::mag(const pointField& points) const
82 {
83  return ::Foam::mag(area(points));
84 }
85 
86 
88 {
89  // for a closed polygon a number of edges is the same as number of points
90  return size();
91 }
92 
93 
95 {
96  return edge(operator[](n), operator[](fcIndex(n)));
97 }
98 
99 
101 {
102  return operator[](fcIndex(i));
103 }
104 
105 
107 {
108  return operator[](rcIndex(i));
109 }
110 
111 
113 {
114  return size() - 2;
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
119 
121 {
122  labelList::operator=(move(l));
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
127 
128 inline bool Foam::operator==(const face& a, const face& b)
129 {
130  return face::compare(a,b) != 0;
131 }
132 
133 
134 inline bool Foam::operator!=(const face& a, const face& b)
135 {
136  return face::compare(a,b) == 0;
137 }
138 
139 
140 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
141 
143 {
144  is >> static_cast<labelList&>(f);
145 
146  // Check state of Ostream
147  is.check("Istream& operator>>(Istream&, face&)");
148 
149  return is;
150 }
151 
152 
153 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
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:28
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:87
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:100
static vector area(const PointField &ps)
Return vector area given face points.
label prevLabel(const label i) const
Previous vertex on face.
Definition: faceI.H:106
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
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:48
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: faceI.H:81
pointField points(const pointField &) const
Return the points corresponding to this face.
Definition: faceI.H:64
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:60
labelList f(nPoints)
void operator=(const UList< label > &)
Assignment to UList operator. Takes linear time.
Definition: List.C:376
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
label nTriangles() const
Size of the face&#39;s triangulation.
Definition: faceI.H:112
volScalarField & p
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:94
bool operator!=(const particle &, const particle &)
Definition: particle.C:1282
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
void operator=(labelList &&)
Move assignment labelList.
Definition: faceI.H:120