Scalar.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-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2021 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
27Typedef
28 Foam::Scalar
29
30Description
31 Floating-point number (float or double)
32
33Note
34 The floor/ceil/round operations could easily be extended to handle
35 VectorSpace when the need arises. For example,
36
37 \code
38 template<class T>
39 struct floorOp
40 {
41 T operator()(const T& x) const WARNRETURN
42 {
43 T ret;
44 for (direction cmpt=0; cmpt < pTraits<T>::nComponents; ++cmpt)
45 {
46 component(ret, cmpt) = std::floor(component(x, cmpt));
47 }
48 return ret;
49 }
50 };
51 \endcode
52
53SourceFiles
54 Scalar.C
55
56\*---------------------------------------------------------------------------*/
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60namespace Foam
61{
62
63// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64
65// Template specialisation for pTraits<Scalar>
66template<>
68{
69 Scalar p_;
70
71public:
72
73 // Typedefs
74
75 //- Component type
77
78 //- Magnitude type
79 typedef Scalar magType;
80
81 //- Equivalent type of labels used for valid component indexing
82 typedef label labelType;
83
84
85 // Member Constants
86
87 //- Dimensionality of space
88 static constexpr direction dim = 3;
89
90 //- Rank of Scalar is 0
91 static constexpr direction rank = 0;
92
93 //- Number of components in Scalar is 1
94 static constexpr direction nComponents = 1;
95
96
97 // Static Data Members
98
99 static const char* const typeName;
100 static const char* const componentNames[];
101 static const Scalar zero;
102 static const Scalar one;
103 static const Scalar max;
104 static const Scalar min;
105 static const Scalar rootMax;
106 static const Scalar rootMin;
107 static const Scalar vsmall;
108
109
110 // Constructors
111
112 //- Copy construct from primitive
113 explicit pTraits(Scalar val) noexcept
114 :
115 p_(val)
116 {}
117
118 //- Read construct from Istream
119 explicit pTraits(Istream& is);
120
121
122 // Member Functions
123
124 //- Return the value
125 operator Scalar() const noexcept
126 {
127 return p_;
128 }
129
130 //- Access the value
131 operator Scalar&() noexcept
132 {
133 return p_;
134 }
135};
136
137
138// * * * * * * * * * * * * * * * IO/Conversion * * * * * * * * * * * * * * * //
139
140//- A word representation of a floating-point value.
141// Uses stringstream instead of std::to_string for more consistent formatting.
142word name(const Scalar val);
143
144
145//- A word representation of a floating-point value.
146template<>
148{
149 word operator()(const Scalar val) const
150 {
151 return Foam::name(val);
152 }
153};
154
155
156//- Parse entire buffer as a float/double, skipping leading/trailing whitespace.
157// \return Parsed value or FatalIOError on any problem
158Scalar ScalarRead(const char* buf);
159
160//- Parse entire buffer as a float/double, skipping leading/trailing whitespace.
161// \return True if successful.
162bool ScalarRead(const char* buf, Scalar& val);
163
164//- Parse entire string as a float/double, skipping leading/trailing whitespace.
165// \return Parsed value or FatalIOError on any problem
166inline Scalar ScalarRead(const std::string& str)
167{
168 return ScalarRead(str.c_str());
169}
170
171//- Parse entire string as a float/double, skipping leading/trailing whitespace.
172// \return True if successful.
173inline bool ScalarRead(const std::string& str, Scalar& val)
174{
175 return ScalarRead(str.c_str(), val);
176}
177
178
179Scalar ScalarRead(Istream& is);
180Istream& operator>>(Istream& is, Scalar& val);
181Ostream& operator<<(Ostream& os, const Scalar val);
182
183
184// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185
186// Standard C++ transcendental functions
188
205
206// Standard ANSI-C (but not in <cmath>) transcendental functions
207
211transFunc(tgamma)
212
219
220
221//- Non-const access to scalar-type (has no components)
223{
224 return val;
225}
226
227
228//- Return scalar value (has no components)
229inline constexpr Scalar component(const Scalar val, const direction) noexcept
230{
231 return val;
232}
233
234
235//- Return 1 if s is greater_equal zero, or otherwise -1
236inline Scalar sign(const Scalar s)
237{
238 return (s >= 0)? 1: -1;
239}
240
241
242//- Return 1 if s is greater than zero, otherwise 1
243inline Scalar pos(const Scalar s)
244{
245 return (s > 0)? 1: 0;
246}
247
248
249//- Return 1 if s is greater_equal zero, or otherwise 0
250inline Scalar pos0(const Scalar s)
251{
252 return (s >= 0)? 1: 0;
253}
254
255
256//- Return 1 if s is less than zero, or otherwise 0
257inline Scalar neg(const Scalar s)
258{
259 return (s < 0)? 1: 0;
260}
261
262
263//- Return 1 if s is less_equal zero, or otherwise 0
264inline Scalar neg0(const Scalar s)
265{
266 return (s <= 0)? 1: 0;
267}
268
269
270//- Return the positive part of s, otherwise zero. Same as max(0, s).
271inline Scalar posPart(const Scalar s)
272{
273 return (s > 0)? s: 0;
274}
275
276
277//- Return the negative part of s, otherwise zero. Same as min(0, s).
278// Does not change the sign
279inline Scalar negPart(const Scalar s)
280{
281 return (s < 0)? s: 0;
282}
283
284
285inline bool equal(const Scalar& s1, const Scalar& s2)
286{
287 return mag(s1 - s2) <= ScalarVSMALL;
288}
289
290
291inline bool notEqual(const Scalar s1, const Scalar s2)
292{
293 return mag(s1 - s2) > ScalarVSMALL;
294}
295
296
297inline Scalar limit(const Scalar s1, const Scalar s2)
298{
299 return (mag(s1) < mag(s2)) ? s1: 0.0;
300}
301
302
303inline Scalar minMod(const Scalar s1, const Scalar s2)
304{
305 return (mag(s1) < mag(s2)) ? s1: s2;
306}
307
308
309inline Scalar magSqr(const Scalar s)
310{
311 return s*s;
312}
313
314
315inline Scalar sqr(const Scalar s)
316{
317 return s*s;
318}
319
320
321inline Scalar pow3(const Scalar s)
322{
323 return s*sqr(s);
324}
325
326
327inline Scalar pow4(const Scalar s)
328{
329 return sqr(sqr(s));
330}
331
332
333inline Scalar pow5(const Scalar s)
334{
335 return s*pow4(s);
336}
337
338
339inline Scalar pow6(const Scalar s)
340{
341 return pow3(sqr(s));
342}
343
344
345inline Scalar pow025(const Scalar s)
346{
347 return sqrt(sqrt(s));
348}
349
350
351inline Scalar inv(const Scalar s)
352{
353 return 1.0/s;
354}
355
356
357inline Scalar dot(const Scalar s1, const Scalar s2)
358{
359 return s1*s2;
360}
361
362
363inline Scalar cmptMultiply(const Scalar s1, const Scalar s2)
364{
365 return s1*s2;
366}
367
368
369inline Scalar cmptPow(const Scalar s1, const Scalar s2)
370{
371 return pow(s1, s2);
372}
373
374
375inline Scalar cmptDivide(const Scalar s1, const Scalar s2)
376{
377 return s1/s2;
378}
379
380
381inline Scalar cmptMax(const Scalar s)
382{
383 return s;
384}
385
386
387inline Scalar cmptMin(const Scalar s)
388{
389 return s;
390}
391
392
393inline Scalar cmptAv(const Scalar s)
394{
395 return s;
396}
397
398
399inline Scalar cmptSqr(const Scalar s)
400{
401 return sqr(s);
402}
403
404
405inline Scalar cmptMag(const Scalar s)
406{
407 return mag(s);
408}
409
410
412{
413 return magSqr(s);
414}
415
416
417inline Scalar sqrtSumSqr(const Scalar a, const Scalar b)
418{
419 const Scalar maga = mag(a);
420 const Scalar magb = mag(b);
421
422 if (maga > magb)
423 {
424 return maga*sqrt(1.0 + sqr(magb/maga));
425 }
426 else
427 {
428 return magb < ScalarVSMALL ? 0.0 : magb*sqrt(1.0 + sqr(maga/magb));
429 }
430}
431
432
433//- Stabilisation around zero for division
434inline Scalar stabilise(const Scalar s, const Scalar tol)
435{
436 if (s >= 0)
437 {
438 return s + tol;
439 }
440 else
441 {
442 return s - tol;
443 }
444}
445
446
447// Specializations
448
449// Default definition in ops.H
450template<class T> struct compareOp;
451
452//- Compare scalar values
453template<>
455{
457
458 //- Construct with specified tolerance (non-negative value)
460 :
461 tolerance(tol)
462 {}
463
464 Scalar operator()(const Scalar& a, const Scalar& b) const
465 {
466 return (mag(a - b) <= tolerance) ? 0 : (a - b);
467 }
468};
469
470
471// Default definition in ops.H
472template<class T> struct equalOp;
473
474//- Compare scalar values for equality
475template<>
477{
479
480 //- Construct with specified tolerance (non-negative value)
482 :
483 tolerance(tol)
484 {}
485
486 bool operator()(const Scalar& a, const Scalar& b) const
487 {
488 return mag(a - b) <= tolerance;
489 }
490};
491
492
493// Default definition in ops.H
494template<class T> struct notEqualOp;
495
496//- Compare scalar values for inequality
497template<>
499{
501
502 //- Construct with specified tolerance (non-negative value)
504 :
505 tolerance(tol)
506 {}
507
508 bool operator()(const Scalar& a, const Scalar& b) const
509 {
510 return mag(a - b) > tolerance;
511 }
512};
513
514
515
516// Default definition in ops.H (future?)
517template<class T> struct floorOp;
518
519//- Round scalar downwards - functor version of std::floor
520template<>
522{
523 Scalar operator()(const Scalar& x) const
524 {
525 return std::floor(x);
526 }
527};
528
529
530// Default definition in ops.H (future?)
531template<class T> struct ceilOp;
532
533//- Round scalar upwards - functor version of std::ceil
534template<>
536{
537 Scalar operator()(const Scalar& x) const
538 {
539 return std::ceil(x);
540 }
541};
542
543
544// Default definition in ops.H (future?)
545template<class T> struct roundOp;
546
547//- Round scalar - functor version of std::round
548template<>
550{
551 Scalar operator()(const Scalar& x) const
552 {
553 return std::round(x);
554 }
555};
556
557
558// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
559
560} // End namespace Foam
561
562// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Scalar cmptType
Component type.
Definition: Scalar.H:76
pTraits(Scalar val) noexcept
Copy construct from primitive.
Definition: Scalar.H:113
static const Scalar one
Definition: Scalar.H:102
static const Scalar min
Definition: Scalar.H:104
static const Scalar rootMax
Definition: Scalar.H:105
static const Scalar max
Definition: Scalar.H:103
static const Scalar rootMin
Definition: Scalar.H:106
Scalar magType
Magnitude type.
Definition: Scalar.H:79
label labelType
Equivalent type of labels used for valid component indexing.
Definition: Scalar.H:82
static const Scalar vsmall
Definition: Scalar.H:107
static const Scalar zero
Definition: Scalar.H:101
static const char *const typeName
Definition: Scalar.H:99
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define transFunc(func)
#define ScalarVSMALL
Definition: doubleScalar.C:43
#define Scalar
Definition: doubleScalar.C:41
#define ScalarRead
Definition: doubleScalar.C:46
#define besselFunc2(func)
Definition: doubleScalar.H:132
#define besselFunc(func)
Definition: doubleScalar.H:127
OBJstream os(runTime.globalPath()/outputName)
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.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:291
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
Scalar minMod(const Scalar s1, const Scalar s2)
Definition: Scalar.H:303
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
dimensionedScalar cosh(const dimensionedScalar &ds)
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:369
dimensionedScalar sin(const dimensionedScalar &ds)
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components)
Definition: label.H:123
dimensionedScalar tanh(const dimensionedScalar &ds)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar erf(const dimensionedScalar &ds)
Scalar cmptSqr(const Scalar s)
Definition: Scalar.H:399
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar log10(const dimensionedScalar &ds)
complex limit(const complex &, const complex &)
Definition: complexI.H:263
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
void cmptMagSqr(Field< Type > &res, const UList< Type > &f)
void dot(FieldField< Field1, typename innerProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Istream & operator>>(Istream &, directionInfo &)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:56
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar atanh(const dimensionedScalar &ds)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const direction noexcept
Definition: Scalar.H:223
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &df)
dimensionedScalar atan(const dimensionedScalar &ds)
Scalar sqrtSumSqr(const Scalar a, const Scalar b)
Definition: Scalar.H:417
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)
volScalarField & b
Definition: createFields.H:27
Scalar operator()(const Scalar &x) const
Definition: Scalar.H:537
const Scalar tolerance
Definition: Scalar.H:456
Scalar operator()(const Scalar &a, const Scalar &b) const
Definition: Scalar.H:464
compareOp(Scalar tol=ScalarVSMALL)
Construct with specified tolerance (non-negative value)
Definition: Scalar.H:459
Three-way comparison operation of two parameters,.
Definition: ops.H:269
const Scalar tolerance
Definition: Scalar.H:478
equalOp(Scalar tol=ScalarVSMALL)
Construct with specified tolerance (non-negative value)
Definition: Scalar.H:481
bool operator()(const Scalar &a, const Scalar &b) const
Definition: Scalar.H:486
Scalar operator()(const Scalar &x) const
Definition: Scalar.H:523
word operator()(const Scalar val) const
Definition: Scalar.H:149
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:238
const Scalar tolerance
Definition: Scalar.H:500
notEqualOp(Scalar tol=ScalarVSMALL)
Construct with specified tolerance (non-negative value)
Definition: Scalar.H:503
bool operator()(const Scalar &a, const Scalar &b) const
Definition: Scalar.H:508
Scalar operator()(const Scalar &x) const
Definition: Scalar.H:551