edgeI.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-2018 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 #include "IOstreams.H"
27 #include "Swap.H"
28 
29 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
30 
31 
32 // return
33 // - 0: different
34 // - +1: identical
35 // - -1: same edge, but different orientation
36 inline int Foam::edge::compare(const edge& a, const edge& b)
37 {
38  if (a[0] == b[0] && a[1] == b[1])
39  {
40  return 1;
41  }
42  else if (a[0] == b[1] && a[1] == b[0])
43  {
44  return -1;
45  }
46  else
47  {
48  return 0;
49  }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 {}
57 
58 
59 inline Foam::edge::edge(const label a, const label b)
60 {
61  start() = a;
62  end() = b;
63 }
64 
65 
67 {
68  start() = a[0];
69  end() = a[1];
70 }
71 
72 
74 :
75  FixedList<label, 2>(is)
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  return operator[](0);
84 }
85 
87 {
88  return operator[](0);
89 }
90 
91 
93 {
94  return operator[](1);
95 }
96 
98 {
99  return operator[](1);
100 }
101 
102 
103 inline bool Foam::edge::connected(const edge& a) const
104 {
105  if
106  (
107  start() == a.start()
108  || start() == a.end()
109  || end() == a.start()
110  || end() == a.end()
111  )
112  {
113  return true;
114  }
115  else
116  {
117  return false;
118  }
119 }
120 
121 
123 {
124  if (start() == a.start() || start() == a.end())
125  {
126  return start();
127  }
128  else if (end() == a.start() || end() == a.end())
129  {
130  return end();
131  }
132  else
133  {
134  // No shared vertex.
135  return -1;
136  }
137 }
138 
139 
141 {
142  if (a == start())
143  {
144  return end();
145  }
146  else if (a == end())
147  {
148  return start();
149  }
150  else
151  {
152  // The given vertex is not on the edge in the first place.
153  return -1;
154  }
155 }
156 
157 inline void Foam::edge::flip()
158 {
159  Swap(operator[](0), operator[](1));
160 }
161 
162 
164 {
165  return edge(end(), start());
166 }
167 
168 
170 {
171  return 0.5*(p[start()] + p[end()]);
172 }
173 
174 
176 {
177  return p[end()] - p[start()];
178 }
179 
180 
181 inline Foam::scalar Foam::edge::mag(const pointField& p) const
182 {
184 }
185 
186 
188 {
189  return linePointRef(p[start()], p[end()]);
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * //
194 
195 inline bool Foam::operator==(const edge& a, const edge& b)
196 {
197  return edge::compare(a,b) != 0;
198 }
199 
200 
201 inline bool Foam::operator!=(const edge& a, const edge& b)
202 {
203  return edge::compare(a,b) == 0;
204 }
205 
206 
207 // ************************************************************************* //
A line primitive.
Definition: line.H:56
edge reverseEdge() const
Return reverse edge.
Definition: edgeI.H:163
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static int compare(const edge &, const edge &)
Compare edges.
Definition: edgeI.H:36
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
label otherVertex(const label a) const
Given one vertex, return the other.
Definition: edgeI.H:140
label commonVertex(const edge &a) const
Return common vertex.
Definition: edgeI.H:122
point centre(const pointField &) const
Return centre (centroid)
Definition: edgeI.H:169
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
linePointRef line(const pointField &) const
Return edge line.
Definition: edgeI.H:187
void Swap(T &a, T &b)
Definition: Swap.H:43
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
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 > &)
label & operator[](const label)
Return element of FixedList.
Definition: FixedListI.H:247
scalar mag(const pointField &) const
Return scalar magnitude.
Definition: edgeI.H:181
bool connected(const edge &a) const
Return true if connected to given edge.
Definition: edgeI.H:103
vector vec(const pointField &) const
Return the vector (end - start)
Definition: edgeI.H:175
label end() const
Return end vertex label.
Definition: edgeI.H:92
Swap its arguments.
edge()
Null constructor for lists.
Definition: edgeI.H:55
dimensioned< scalar > mag(const dimensioned< Type > &)
void flip()
Flip the edge in-place.
Definition: edgeI.H:157
volScalarField & p
bool operator!=(const particle &, const particle &)
Definition: particle.C:1282
label start() const
Return start vertex label.
Definition: edgeI.H:81