Loading...
Searching...
No Matches
VectorI.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-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2023 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
31template<class Cmpt>
35{}
36
37
38template<class Cmpt>
39template<class Cmpt2>
41(
42 const VectorSpace<Vector<Cmpt2>, Cmpt2, 3>& vs
44:
45 Vector::vsType(vs)
46{}
47
48
49template<class Cmpt>
51(
52 const Cmpt& vx,
53 const Cmpt& vy,
54 const Cmpt& vz
55)
56{
57 this->x() = vx;
58 this->y() = vy;
59 this->z() = vz;
60}
61
62
63template<class Cmpt>
65:
67{}
68
69
70// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71
72template<class Cmpt>
74(
76) const noexcept
77{
78 return *this;
79}
80
81
82template<class Cmpt>
83inline Foam::scalar Foam::Vector<Cmpt>::magSqr() const
84{
85 return
86 (
87 Foam::magSqr(this->x())
88 + Foam::magSqr(this->y())
89 + Foam::magSqr(this->z())
90 );
91}
92
93
94template<class Cmpt>
95inline Foam::scalar Foam::Vector<Cmpt>::mag() const
96{
97 return ::sqrt(this->magSqr());
98}
99
100
101template<class Cmpt>
102inline Foam::scalar Foam::Vector<Cmpt>::distSqr(const Vector<Cmpt>& v2) const
103{
104 return
105 (
106 Foam::magSqr(v2.x() - this->x())
107 + Foam::magSqr(v2.y() - this->y())
108 + Foam::magSqr(v2.z() - this->z())
109 );
110}
111
112
113template<class Cmpt>
114inline Foam::scalar Foam::Vector<Cmpt>::dist(const Vector<Cmpt>& v2) const
115{
116 return ::sqrt(this->distSqr(v2));
117}
118
119
120template<class Cmpt>
122{
123 #ifdef __clang__
124 volatile // Use volatile to avoid aggressive branch optimization
125 #endif
126 const scalar s = this->mag();
127
128 if (s < tol)
129 {
130 *this = Zero;
131 }
132 else
133 {
134 *this /= s;
135 }
137 return *this;
138}
139
140
141template<class Cmpt>
142inline Foam::Vector<Cmpt>&
144{
145 *this -= (*this & unitVec) * unitVec;
146 return *this;
147}
148
149
150// * * * * * * * * * * * * * * * Member Operations * * * * * * * * * * * * * //
151
152template<class Cmpt>
153inline Cmpt Foam::Vector<Cmpt>::inner(const Vector<Cmpt>& v2) const
154{
155 const Vector<Cmpt>& v1 = *this;
157 return (v1.x()*v2.x() + v1.y()*v2.y() + v1.z()*v2.z());
158}
159
160
161template<class Cmpt>
164{
165 const Vector<Cmpt>& v1 = *this;
166
167 return Vector<Cmpt>
169 (v1.y()*v2.z() - v1.z()*v2.y()),
170 (v1.z()*v2.x() - v1.x()*v2.z()),
171 (v1.x()*v2.y() - v1.y()*v2.x())
172 );
174
175
176// * * * * * * * * * * * * * Comparison Operations * * * * * * * * * * * * * //
177
178template<class Cmpt>
179inline bool
182 return
183 (
184 (a.x() < b.x()) // Component is less
185 ||
186 (
187 !(b.x() < a.x()) // Equal? Check next component
188 &&
189 (
190 (a.y() < b.y()) // Component is less
191 ||
192 (
193 !(b.y() < a.y()) // Equal? Check next component
194 && (a.z() < b.z())
195 )
196 )
198 );
199}
200
202template<class Cmpt>
203inline bool
205{
206 return
208 (a.y() < b.y()) // Component is less
209 ||
210 (
211 !(b.y() < a.y()) // Equal? Check next component
212 &&
213 (
214 (a.z() < b.z()) // Component is less
215 ||
216 (
217 !(b.z() < a.z()) // Equal? Check next component
218 && (a.x() < b.x())
219 )
220 )
222 );
223}
224
226template<class Cmpt>
227inline bool
229{
230 return
231 (
232 (a.z() < b.z()) // Component is less
233 ||
235 !(b.z() < a.z()) // Equal? Check next component
236 &&
237 (
238 (a.x() < b.x()) // Component is less
239 ||
240 (
241 !(b.x() < a.x()) // Equal? Check next component
242 && (a.y() < b.y())
244 )
245 )
246 );
247}
248
249
250// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251
252namespace Foam
253{
255// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
256
257//- Linear interpolation of vectors a and b by factor t
258template<class Cmpt>
259inline Vector<Cmpt> lerp
260(
261 const Vector<Cmpt>& a,
262 const Vector<Cmpt>& b,
263 const scalar t
264)
265{
266 const scalar onet = (1-t);
267
268 return Vector<Cmpt>
269 (
270 onet*a.x() + t*b.x(),
271 onet*a.y() + t*b.y(),
272 onet*a.z() + t*b.z()
273 );
274}
275
276
277// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
278
279//- Dummy innerProduct for scalar
280// Allows the construction of vtables for virtual member functions
281// involving the inner-products of fields
282// for which a "NotImplemented" specialization for scalar is provided.
283template<class Cmpt>
284class innerProduct<Vector<Cmpt>, scalar>
285{
286public:
287
288 typedef scalar type;
289};
290
292template<class Cmpt>
293inline Cmpt operator&(const Vector<Cmpt>& v1, const Vector<Cmpt>& v2)
294{
295 return v1.inner(v2);
296}
297
299template<class Cmpt>
300inline Vector<Cmpt> operator^(const Vector<Cmpt>& v1, const Vector<Cmpt>& v2)
301{
302 return v1.cross(v2);
303}
304
306template<class Cmpt>
307inline Vector<Cmpt> operator*(const Cmpt& s, const Vector<Cmpt>& v)
308{
309 return Vector<Cmpt>(s*v.x(), s*v.y(), s*v.z());
310}
311
313template<class Cmpt>
314inline Vector<Cmpt> operator*(const Vector<Cmpt>& v, const Cmpt& s)
315{
316 return s*v;
317}
318
319
320// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321
322} // End namespace Foam
323
324// ************************************************************************* //
scalar y
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
Templated vector space.
Definition VectorSpace.H:75
VectorSpace< Vector< Cmpt >, Cmpt, Ncmpts > vsType
Definition VectorSpace.H:86
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
static bool less_xyz(const Vector< Cmpt > &a, const Vector< Cmpt > &b)
Lexicographically compare a and b with order (x:y:z).
Definition VectorI.H:173
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
scalar dist(const Vector< Cmpt > &v2) const
The L2-norm distance from another vector. The mag() of the difference.
Definition VectorI.H:107
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
Cmpt inner(const Vector< Cmpt > &v2) const
Scalar-product of this with another Vector.
Definition VectorI.H:146
Vector< Cmpt > cross(const Vector< Cmpt > &v2) const
Cross-product of this with another Vector.
Definition VectorI.H:156
Vector< Cmpt > & normalise(const scalar tol=ROOTVSMALL)
Inplace normalise the vector by its magnitude.
Definition VectorI.H:114
const Vector< Cmpt > & centre(const Foam::UList< Vector< Cmpt > > &) const noexcept
Return this (for point which is a typedef to Vector<scalar>).
Definition VectorI.H:67
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
Definition VectorI.H:95
static bool less_zxy(const Vector< Cmpt > &a, const Vector< Cmpt > &b)
Lexicographically compare a and b with order (z:x:y).
Definition VectorI.H:221
scalar mag() const
The length (L2-norm) of the vector.
Definition VectorI.H:88
static bool less_yzx(const Vector< Cmpt > &a, const Vector< Cmpt > &b)
Lexicographically compare a and b with order (y:z:x).
Definition VectorI.H:197
Vector()=default
Default construct.
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
scalar magSqr() const
The length (L2-norm) squared of the vector.
Definition VectorI.H:76
Vector< Cmpt > & removeCollinear(const Vector< Cmpt > &unitVec)
Inplace removal of components that are collinear to the given unit vector.
Definition VectorI.H:136
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
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))
Namespace for OpenFOAM.
tmp< faMatrix< Type > > operator*(const areaScalarField::Internal &, const faMatrix< Type > &)
bitSet operator^(const bitSet &a, const bitSet &b)
Bitwise-XOR of two bitsets to form a unique bit-set.
Definition bitSetI.H:702
tmp< GeometricField< Type, faPatchField, areaMesh > > operator&(const faMatrix< Type > &, const DimensionedField< Type, areaMesh > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const direction noexcept
Definition scalarImpl.H:265
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
volScalarField & b