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