Scalar.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-2023 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 Typedef
25  Foam::Scalar
26 
27 Description
28  Single floating point number (float or double)
29 
30 SourceFiles
31  Scalar.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 // template specialisation for pTraits<Scalar>
43 template<>
45 {
46  Scalar p_;
47 
48 public:
49 
50  //- Component type
51  typedef Scalar cmptType;
52 
53  //- Equivalent type of labels used for valid component indexing
54  typedef label labelType;
55 
56 
57  // Member constants
58 
59  //- Dimensionality of space
60  static const direction dim = 3;
61 
62  //- Rank of Scalar is 0
63  static const direction rank = 0;
64 
65  //- Number of components in Scalar is 1
66  static const direction nComponents = 1;
67 
68 
69  // Static Data Members
70 
71  static const char* const typeName;
72  static const char* const componentNames[];
73  static const Scalar zero;
74  static const Scalar one;
75  static const Scalar max;
76  static const Scalar min;
77  static const Scalar rootMax;
78  static const Scalar rootMin;
79  static const Scalar nan;
80 
81 
82  // Constructors
83 
84  //- Construct from primitive
85  explicit pTraits(const Scalar&);
86 
87  //- Construct from Istream
88  pTraits(Istream&);
89 
90 
91  // Member Functions
92 
93  //- Access to the Scalar value
94  operator Scalar() const
95  {
96  return p_;
97  }
98 
99  //- Access to the Scalar value
100  operator Scalar&()
101  {
102  return p_;
103  }
104 };
105 
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 //- Return a string representation of a Scalar
110 word name(const Scalar);
111 
112 
113 // Standard C++ transcendental functions
115 
117 transFunc(exp)
118 transFunc(log)
120 transFunc(sin)
121 transFunc(cos)
122 transFunc(tan)
132 
133 // Standard ANSI-C (but not in <cmath>) transcendental functions
134 
135 transFunc(erf)
138 transFunc(tgamma)
139 
140 transFunc(j0)
141 transFunc(j1)
142 
143 transFunc(y0)
144 transFunc(y1)
145 
146 
147 inline Scalar& setComponent(Scalar& s, const direction)
148 {
149  return s;
150 }
151 
152 
153 inline Scalar component(const Scalar s, const direction)
154 {
155  return s;
156 }
157 
158 
159 //- Return 1 if s is positive or 0 otherwise -1
160 inline int sign(const Scalar s)
161 {
162  return (s >= 0) ? 1: -1;
163 }
164 
165 
166 //- Return 1 if s is positive but not 0
167 inline int pos(const Scalar s)
168 {
169  return (s > 0) ? 1: 0;
170 }
171 
172 
173 //- Return 1 if s is positive or 0
174 inline int pos0(const Scalar s)
175 {
176  return (s >= 0) ? 1: 0;
177 }
178 
179 
180 //- Return 1 if s is negative but not 0
181 inline int neg(const Scalar s)
182 {
183  return (s < 0) ? 1: 0;
184 }
185 
186 
187 //- Return 1 if s is negative or 0
188 inline int neg0(const Scalar s)
189 {
190  return (s <= 0) ? 1: 0;
191 }
192 
193 
194 //- Return the positive part of s
195 inline Scalar posPart(const Scalar s)
196 {
197  return (s > 0)? s: 0;
198 }
199 
200 
201 //- Return the negative part of s.
202 // Note: this function returns the actual negative part of s as a
203 // negative number and does not change the sign
204 inline Scalar negPart(const Scalar s)
205 {
206  return (s < 0)? s: 0;
207 }
208 
209 
210 inline bool equal(const Scalar& s1, const Scalar& s2)
211 {
212  return mag(s1 - s2) <= ScalarVSmall;
213 }
214 
215 
216 inline bool notEqual(const Scalar s1, const Scalar s2)
217 {
218  return mag(s1 - s2) > ScalarVSmall;
219 }
220 
221 
222 inline Scalar limit(const Scalar s1, const Scalar s2)
223 {
224  return (mag(s1) < mag(s2)) ? s1: 0.0;
225 }
226 
227 
228 inline Scalar minMod(const Scalar s1, const Scalar s2)
229 {
230  return (mag(s1) < mag(s2)) ? s1: s2;
231 }
232 
233 
234 inline Scalar magSqr(const Scalar s)
235 {
236  return s*s;
237 }
238 
239 
240 inline Scalar sqr(const Scalar s)
241 {
242  return s*s;
243 }
244 
245 
246 inline Scalar pow3(const Scalar s)
247 {
248  return s*sqr(s);
249 }
250 
251 
252 inline Scalar pow4(const Scalar s)
253 {
254  return sqr(sqr(s));
255 }
256 
257 
258 inline Scalar pow5(const Scalar s)
259 {
260  return s*pow4(s);
261 }
262 
263 
264 inline Scalar pow6(const Scalar s)
265 {
266  return pow3(sqr(s));
267 }
268 
269 
270 inline Scalar pow025(const Scalar s)
271 {
272  return sqrt(sqrt(s));
273 }
274 
275 
276 inline Scalar inv(const Scalar s)
277 {
278  return 1.0/s;
279 }
280 
281 
282 template<class Type>
283 inline Type dot(const Scalar s, const Type& t)
284 {
285  return s*t;
286 }
287 
288 
289 template<class Type>
290 inline Type dot(const Type& t, const Scalar s)
291 {
292  return t*s;
293 }
294 
295 
296 inline Scalar dot(const Scalar s1, const Scalar s2)
297 {
298  return s1*s2;
299 }
300 
301 
302 inline Scalar cmptMultiply(const Scalar s1, const Scalar s2)
303 {
304  return s1*s2;
305 }
306 
307 
308 inline Scalar cmptPow(const Scalar s1, const Scalar s2)
309 {
310  return pow(s1, s2);
311 }
312 
313 
314 inline Scalar cmptDivide(const Scalar s1, const Scalar s2)
315 {
316  return s1/s2;
317 }
318 
319 
320 inline Scalar cmptMax(const Scalar s)
321 {
322  return s;
323 }
324 
325 
326 inline Scalar cmptMin(const Scalar s)
327 {
328  return s;
329 }
330 
331 
332 inline Scalar cmptAv(const Scalar s)
333 {
334  return s;
335 }
336 
337 
338 inline Scalar cmptSqr(const Scalar s)
339 {
340  return sqr(s);
341 }
342 
343 
344 inline Scalar cmptMag(const Scalar s)
345 {
346  return mag(s);
347 }
348 
349 
350 inline Scalar sqrtSumSqr(const Scalar a, const Scalar b)
351 {
352  Scalar maga = mag(a);
353  Scalar magb = mag(b);
354 
355  if (maga > magb)
356  {
357  return maga*sqrt(1 + sqr(magb/maga));
358  }
359  else
360  {
361  return magb < ScalarVSmall ? 0 : magb*sqrt(1 + sqr(maga/magb));
362  }
363 }
364 
365 
366 //- Stabilisation around zero for division
367 inline Scalar stabilise(const Scalar s, const Scalar small)
368 {
369  if (s >= 0)
370  {
371  return s + small;
372  }
373  else
374  {
375  return s - small;
376  }
377 }
378 
379 
380 // * * * * * * * * * * * * * * * IOstream Functions * * * * * * * * * * * * //
381 
382 Scalar readScalar(Istream&);
383 
384 void writeEntry(Ostream& os, const Scalar value);
385 
386 
387 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
388 
389 Istream& operator>>(Istream&, Scalar&);
390 Ostream& operator<<(Ostream&, const Scalar);
391 
392 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393 
394 } // End namespace Foam
395 
396 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
Scalar cmptType
Component type.
Definition: Scalar.H:51
static const Scalar one
Definition: Scalar.H:74
static const Scalar min
Definition: Scalar.H:76
static const Scalar rootMax
Definition: Scalar.H:77
static const Scalar max
Definition: Scalar.H:75
static const Scalar rootMin
Definition: Scalar.H:78
static const Scalar nan
Definition: Scalar.H:79
label labelType
Equivalent type of labels used for valid component indexing.
Definition: Scalar.H:54
static const Scalar zero
Definition: Scalar.H:73
static const char *const typeName
Definition: Scalar.H:71
Traits class for primitives.
Definition: pTraits.H:53
pTraits(const PrimitiveType &p)
Construct from primitive.
Definition: pTraits.H:60
#define ScalarVSmall
Definition: doubleScalar.C:35
#define Scalar
Definition: doubleScalar.C:33
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:25
Namespace for OpenFOAM.
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensionedScalar erfc(const dimensionedScalar &ds)
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
dimensionedScalar asin(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tan(const dimensionedScalar &ds)
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:216
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar j1(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedScalar y0(const dimensionedScalar &ds)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
label & setComponent(label &l, const direction)
Definition: label.H:86
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar log10(const dimensionedScalar &ds)
complex limit(const complex &, const complex &)
Definition: complexI.H:202
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)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Scalar minMod(const Scalar s1, const Scalar s2)
Definition: Scalar.H:228
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar acosh(const dimensionedScalar &ds)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
dimensionedScalar sqrt(const dimensionedScalar &ds)
Istream & operator>>(Istream &, pistonPointEdgeData &)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar pow4(const dimensionedScalar &ds)
Scalar sqrtSumSqr(const Scalar a, const Scalar b)
Definition: Scalar.H:350
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:296
dimensionedScalar atanh(const dimensionedScalar &ds)
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:308
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) transFunc(j0) transFunc(j1) transFunc(y0) transFunc(y1) inline Scalar &setComponent(Scalar &s
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar posPart(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:45
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar j0(const dimensionedScalar &ds)
Scalar cmptSqr(const Scalar s)
Definition: Scalar.H:338
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar asinh(const dimensionedScalar &ds)