Loading...
Searching...
No Matches
scalarImpl.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-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
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
42 {
43 T ret;
44 for (direction cmpt=0; cmpt < pTraits<T>::nComponents; ++cmpt)
45 {
46 setComponent(ret, cmpt) = std::floor(component(x, cmpt));
47 }
48 return ret;
49 }
50 };
51 \endcode
52
53SourceFiles
54 scalarImpl.txx
55
56\*---------------------------------------------------------------------------*/
57
58// Could also check if 'Scalar' is defined...
59#ifndef Foam_use_scalarImpl_header
60#error "scalarImpl.H" is only to be included internally
61#endif
62
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66namespace Foam
67{
68
69// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70
71// Template specialisation for pTraits<Scalar>
72template<>
73class pTraits<Scalar>
74{
75 Scalar p_;
76
77public:
78
79 // Typedefs
81 //- Component type
82 typedef Scalar cmptType;
83
84 //- Magnitude type
85 typedef Scalar magType;
86
87 //- Equivalent type of labels used for valid component indexing
88 typedef label labelType;
89
91 // Member Constants
92
93 //- Dimensionality of space
94 static constexpr direction dim = 3;
95
96 //- Rank of Scalar is 0
97 static constexpr direction rank = 0;
99 //- Number of components in Scalar is 1
100 static constexpr direction nComponents = 1;
101
102
103 // Limits
104
105 static constexpr Scalar min_ = -ScalarVGREAT;
106 static constexpr Scalar max_ = ScalarVGREAT;
107 static constexpr Scalar vsmall_ = ScalarVSMALL;
109
110 // Static Data Members
111
112 static const char* const typeName;
113 static const char* const componentNames[];
114 static const Scalar zero;
115 static const Scalar one;
116 static const Scalar max;
117 static const Scalar min;
118 static const Scalar rootMax;
119 static const Scalar rootMin;
120 static const Scalar vsmall;
123 // Constructors
125 //- Copy construct from primitive
126 explicit pTraits(Scalar val) noexcept : p_(val) {}
128 //- Read construct from Istream (uses tokenizer)
129 explicit pTraits(Istream& is);
130
131
132 // Member Functions
133
134 //- Return the value
135 operator Scalar() const noexcept { return p_; }
137 //- Access the value
138 operator Scalar&() noexcept { return p_; }
139};
140
142// * * * * * * * * * * * * * * * IO/Conversion * * * * * * * * * * * * * * * //
143
144//- A word representation of a floating-point value.
145// Uses stringstream instead of std::to_string for more consistent formatting.
146word name(const Scalar val);
147
148
149//- A word representation of a floating-point value.
150template<>
151struct nameOp<Scalar>
152{
153 word operator()(const Scalar val) const
155 return Foam::name(val);
156 }
157};
158
159
160//- Parse entire buffer as a float/double, skipping leading/trailing whitespace.
161// \return Parsed value or FatalIOError on any problem
162Scalar ScalarRead(const char* buf);
163
164//- Parse entire buffer as a float/double, skipping leading/trailing whitespace.
165// \return True if successful.
166bool ScalarRead(const char* buf, Scalar& val);
167
168//- Parse entire string as a float/double, skipping leading/trailing whitespace.
169// \return Parsed value or FatalIOError on any problem
170inline Scalar ScalarRead(const std::string& str)
171{
172 return ScalarRead(str.c_str());
173}
175//- Parse entire string as a float/double, skipping leading/trailing whitespace.
176// \return True if successful.
177inline bool ScalarRead(const std::string& str, Scalar& val)
178{
179 return ScalarRead(str.c_str(), val);
180}
181
182
183//- Read from Istream (uses tokenizer)
185Istream& operator>>(Istream& is, Scalar& val);
187
188
189// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190
191// Standard C++ transcendental functions
211// Standard ANSI-C (but not in <cmath>) transcendental functions
212
216transFunc(tgamma)
217
224
225
226//- Non-const access to scalar-type (has no components)
228{
229 return val;
230}
231
232
233//- Return scalar value (has no components)
234inline constexpr Scalar component(const Scalar val, const direction) noexcept
235{
236 return val;
237}
238
239
240//- Return 1 if s is greater_equal zero, or otherwise -1
241inline Scalar sign(const Scalar s) noexcept
242{
243 return (s >= 0)? 1: -1;
244}
245
246
247//- Return 1 if s is greater than zero, otherwise 0
248inline Scalar pos(const Scalar s) noexcept
249{
250 return (s > 0)? 1: 0;
251}
252
253
254//- Return 1 if s is greater_equal zero, or otherwise 0
255inline Scalar pos0(const Scalar s) noexcept
256{
257 return (s >= 0)? 1: 0;
258}
259
260
261//- Return 1 if s is less than zero, or otherwise 0
262inline Scalar neg(const Scalar s) noexcept
263{
264 return (s < 0)? 1: 0;
265}
266
267
268//- Return 1 if s is less_equal zero, or otherwise 0
269inline Scalar neg0(const Scalar s) noexcept
270{
271 return (s <= 0)? 1: 0;
272}
274
275//- Return the positive part of s, otherwise zero. Same as max(0, s).
276inline Scalar posPart(const Scalar s) noexcept
277{
278 return (s > 0)? s: 0;
279}
280
281
282//- Return the negative part of s, otherwise zero. Same as min(0, s).
283// Does not change the sign
284inline Scalar negPart(const Scalar s) noexcept
285{
286 return (s < 0)? s: 0;
287}
288
289
290inline bool equal(const Scalar& a, const Scalar& b)
292 return mag(a - b) <= ScalarVSMALL;
293}
294
295
296inline bool notEqual(const Scalar a, const Scalar b)
297{
298 return mag(a - b) > ScalarVSMALL;
299}
301
302//- Clamp scalar value to a 0-1 range
303inline Scalar clamp(const Scalar val, Foam::zero_one)
304{
305 return (val < 0) ? 0 : (1 < val) ? 1 : val;
306}
307
308
309// Fast implementation, and with scalar promotion of upper/lower limits
310inline Scalar clamp(const Scalar val, const Scalar lower, const Scalar upper)
311{
312 return (val < lower) ? lower : (upper < val) ? upper : val;
313}
314
315
316inline Scalar limit(const Scalar s1, const Scalar s2)
317{
318 return (mag(s1) < mag(s2)) ? s1: 0.0;
319}
320
321
322inline Scalar minMod(const Scalar s1, const Scalar s2)
323{
324 return (mag(s1) < mag(s2)) ? s1: s2;
325}
326
328//- Linear interpolation of scalars a and b by factor t
329// Explicit form (not fused multiply add etc) for exact values at end points.
330// Unlike std::lerp (C++20) there are no isfinite() checks etc.
331inline constexpr Scalar lerp(const Scalar a, const Scalar b, const Scalar t)
332{
333 return (Scalar{1}-t)*a + t*b;
334}
335
336
337inline Scalar magSqr(const Scalar s)
339 return s*s;
340}
341
342
343inline Scalar sqr(const Scalar s)
345 return s*s;
346}
347
348
349inline Scalar pow3(const Scalar s)
351 return s*sqr(s);
352}
353
354
355inline Scalar pow4(const Scalar s)
356{
357 return sqr(sqr(s));
358}
360
361inline Scalar pow5(const Scalar s)
362{
363 return s*pow4(s);
364}
365
367inline Scalar pow6(const Scalar s)
368{
369 return pow3(sqr(s));
370}
371
373inline Scalar pow025(const Scalar s)
374{
375 return sqrt(sqrt(s));
376}
377
379inline Scalar inv(const Scalar s)
380{
381 return 1.0/s;
382}
383
384
385inline Scalar dot(const Scalar s1, const Scalar s2)
386{
387 return s1*s2;
388}
389
391inline Scalar cmptMultiply(const Scalar s1, const Scalar s2)
392{
393 return s1*s2;
394}
395
397inline Scalar cmptPow(const Scalar s1, const Scalar s2)
398{
399 return pow(s1, s2);
400}
401
403inline Scalar cmptDivide(const Scalar s1, const Scalar s2)
404{
405 return s1/s2;
406}
407
409inline Scalar cmptMax(const Scalar s)
410{
411 return s;
412}
413
415inline Scalar cmptMin(const Scalar s)
416{
417 return s;
418}
419
421inline Scalar cmptAv(const Scalar s)
422{
423 return s;
424}
425
427inline Scalar cmptSqr(const Scalar s)
428{
429 return sqr(s);
430}
431
433inline Scalar cmptMag(const Scalar s)
434{
435 return mag(s);
436}
437
439inline Scalar cmptMagSqr(const Scalar s)
440{
441 return magSqr(s);
442}
443
445inline Scalar sqrtSumSqr(const Scalar a, const Scalar b)
446{
447 const Scalar maga = mag(a);
448 const Scalar magb = mag(b);
449
450 if (maga > magb)
451 {
452 return maga*sqrt(1.0 + sqr(magb/maga));
453 }
454 else
455 {
456 return magb < ScalarVSMALL ? 0.0 : magb*sqrt(1.0 + sqr(maga/magb));
457 }
458}
459
460
461//- Stabilisation around zero for division
462inline Scalar stabilise(const Scalar s, const Scalar tol)
463{
464 if (s >= 0)
465 {
466 return s + tol;
467 }
468 else
469 {
470 return s - tol;
471 }
472}
473
475// Specializations
476
477// Default definition in ops.H
478template<class T> struct compareOp;
479
480//- Compare scalar values
481template<>
482struct compareOp<Scalar>
483{
484 const Scalar tolerance;
485
486 //- Construct with specified tolerance (non-negative value)
488 :
489 tolerance(tol)
490 {}
491
492 Scalar operator()(const Scalar& a, const Scalar& b) const
493 {
494 return (mag(a - b) <= tolerance) ? 0 : (a - b);
495 }
496};
497
499// Default definition in ops.H
500template<class T> struct equalOp;
501
502//- Compare scalar values for equality
503template<>
505{
506 const Scalar tolerance;
507
508 //- Construct with specified tolerance (non-negative value)
510 :
511 tolerance(tol)
512 {}
513
514 bool operator()(const Scalar& a, const Scalar& b) const
515 {
516 return mag(a - b) <= tolerance;
517 }
518};
519
520
521// Default definition in ops.H
522template<class T> struct notEqualOp;
524//- Compare scalar values for inequality
525template<>
526struct notEqualOp<Scalar>
527{
528 const Scalar tolerance;
529
530 //- Construct with specified tolerance (non-negative value)
532 :
533 tolerance(tol)
534 {}
535
536 bool operator()(const Scalar& a, const Scalar& b) const
537 {
538 return mag(a - b) > tolerance;
539 }
540};
541
542
543
544// Default definition in ops.H (future?)
545template<class T> struct floorOp;
546
547//- Round scalar downwards - functor version of std::floor
548template<>
549struct floorOp<Scalar>
550{
551 Scalar operator()(const Scalar& x) const
553 return std::floor(x);
554 }
555};
556
558// Default definition in ops.H (future?)
559template<class T> struct ceilOp;
560
561//- Round scalar upwards - functor version of std::ceil
562template<>
563struct ceilOp<Scalar>
564{
565 Scalar operator()(const Scalar& x) const
566 {
567 return std::ceil(x);
568 }
569};
570
572// Default definition in ops.H (future?)
573template<class T> struct roundOp;
574
575//- Round scalar - functor version of std::round
576template<>
577struct roundOp<Scalar>
579 Scalar operator()(const Scalar& x) const
580 {
581 return std::round(x);
582 }
584
585
586// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
587
588} // End namespace Foam
589
590// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static const char *const typeName
Definition scalarImpl.H:120
Scalar cmptType
Component type.
Definition scalarImpl.H:80
pTraits(Scalar val) noexcept
Copy construct from primitive.
Definition scalarImpl.H:136
static const Scalar one
Definition scalarImpl.H:123
static const Scalar min
Definition scalarImpl.H:125
static const char *const componentNames[]
Definition scalarImpl.H:121
static const Scalar rootMax
Definition scalarImpl.H:126
static const Scalar max
Definition scalarImpl.H:124
static const Scalar rootMin
Definition scalarImpl.H:127
static constexpr Scalar vsmall_
Definition scalarImpl.H:115
static constexpr direction nComponents
Number of components in Scalar is 1.
Definition scalarImpl.H:108
static constexpr direction rank
Rank of Scalar is 0.
Definition scalarImpl.H:103
static constexpr direction dim
Dimensionality of space.
Definition scalarImpl.H:98
Scalar magType
Magnitude type.
Definition scalarImpl.H:85
pTraits(Istream &is)
Read construct from Istream (uses tokenizer).
label labelType
Equivalent type of labels used for valid component indexing.
Definition scalarImpl.H:90
static const Scalar vsmall
Definition scalarImpl.H:128
static constexpr Scalar max_
Definition scalarImpl.H:114
static const Scalar zero
Definition scalarImpl.H:122
static constexpr Scalar min_
Definition scalarImpl.H:113
pTraits(const Base &obj)
Copy construct from base class.
Definition pTraits.H:72
A class for handling words, derived from Foam::string.
Definition word.H:66
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
#define transFunc(func)
#define besselFunc2(func)
#define ScalarVSMALL
#define ScalarVGREAT
#define Scalar
#define besselFunc(func)
#define ScalarRead
OBJstream os(runTime.globalPath()/outputName)
auto & name
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)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
complex limit(const complex &c1, const complex &c2)
Definition complexI.H:235
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
Scalar minMod(const Scalar s1, const Scalar s2)
Definition scalarImpl.H:378
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 scalarImpl.H:456
dimensionedScalar sin(const dimensionedScalar &ds)
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components).
Definition label.H:160
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 scalarImpl.H:486
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
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).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Istream & operator>>(Istream &, directionInfo &)
void cmptMagSqr(Field< Type > &result, const UList< Type > &f1)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
uint8_t direction
Definition direction.H:49
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 scalarImpl.H:265
transFunc(sqrt) transFunc(cbrt) transFunc(exp) transFunc(log) transFunc(log10) transFunc(sin) transFunc(cos) transFunc(tan) transFunc(asin) transFunc(acos) transFunc(atan) transFunc(sinh) transFunc(cosh) transFunc(tanh) transFunc(asinh) transFunc(acosh) transFunc(atanh) transFunc(erf) transFunc(erfc) transFunc(lgamma) transFunc(tgamma) besselFunc(j0) besselFunc(j1) besselFunc(y0) besselFunc(y1) besselFunc2(jn) besselFunc2(yn) inline Scalar &setComponent(Scalar &val
Non-const access to scalar-type (has no components).
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 expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
bool notEqual(const Scalar a, const Scalar b)
Definition scalarImpl.H:350
dimensionedScalar atan(const dimensionedScalar &ds)
Scalar sqrtSumSqr(const Scalar a, const Scalar b)
Definition scalarImpl.H:504
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)
volScalarField & b
Scalar operator()(const Scalar &x) const
Definition scalarImpl.H:642
Scalar operator()(const Scalar &a, const Scalar &b) const
Definition scalarImpl.H:557
compareOp(Scalar tol=ScalarVSMALL)
Construct with specified tolerance (non-negative value).
Definition scalarImpl.H:552
Three-way comparison operation of two parameters,.
Definition ops.H:259
int operator()(const T &a, const T &b) const
Definition ops.H:260
equalOp(Scalar tol=ScalarVSMALL)
Construct with specified tolerance (non-negative value).
Definition scalarImpl.H:578
bool operator()(const Scalar &a, const Scalar &b) const
Definition scalarImpl.H:583
bool operator()(const T &x, const T &y) const
Definition ops.H:223
Scalar operator()(const Scalar &x) const
Definition scalarImpl.H:626
word operator()(const Scalar val) const
Definition scalarImpl.H:174
Extract name (as a word) from an object, typically using its name() method.
Definition word.H:341
word operator()(const T &obj) const
Definition word.H:342
notEqualOp(Scalar tol=ScalarVSMALL)
Construct with specified tolerance (non-negative value).
Definition scalarImpl.H:604
bool operator()(const Scalar &a, const Scalar &b) const
Definition scalarImpl.H:609
bool operator()(const T &x, const T &y) const
Definition ops.H:224
Scalar operator()(const Scalar &x) const
Definition scalarImpl.H:658