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 | Copyright (C) 2011-2016 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 
113 {
114  return s;
115 }
116 
117 
118 inline Scalar component(const Scalar s, const direction)
119 {
120  return s;
121 }
122 
123 
124 inline Scalar sign(const Scalar s)
125 {
126  return (s >= 0)? 1: -1;
127 }
128 
129 
130 inline Scalar pos(const Scalar s)
131 {
132  return (s >= 0)? 1: 0;
133 }
134 
135 
136 inline Scalar neg(const Scalar s)
137 {
138  return (s < 0)? 1: 0;
139 }
140 
141 
142 //- Return the positive part of s
143 inline Scalar posPart(const Scalar s)
144 {
145  return (s > 0)? s: 0;
146 }
147 
148 
149 //- Return the negative part of s.
150 // Note: this function returns the actual negative part of s as a
151 // negative number and does not change the sign
152 inline Scalar negPart(const Scalar s)
153 {
154  return (s < 0)? s: 0;
155 }
156 
157 
158 inline bool equal(const Scalar& s1, const Scalar& s2)
159 {
160  return mag(s1 - s2) <= ScalarVSMALL;
161 }
162 
163 
164 inline bool notEqual(const Scalar s1, const Scalar s2)
165 {
166  return mag(s1 - s2) > ScalarVSMALL;
167 }
168 
169 
170 inline Scalar limit(const Scalar s1, const Scalar s2)
171 {
172  return (mag(s1) < mag(s2)) ? s1: 0.0;
173 }
174 
175 
176 inline Scalar minMod(const Scalar s1, const Scalar s2)
177 {
178  return (mag(s1) < mag(s2)) ? s1: s2;
179 }
180 
181 
182 inline Scalar magSqr(const Scalar s)
183 {
184  return s*s;
185 }
186 
187 
188 inline Scalar sqr(const Scalar s)
189 {
190  return s*s;
191 }
192 
193 
194 inline Scalar sqrtSumSqr(const Scalar a, const Scalar b)
195 {
196  Scalar maga = mag(a);
197  Scalar magb = mag(b);
198 
199  if (maga > magb)
200  {
201  return maga*sqrt(1.0 + sqr(magb/maga));
202  }
203  else
204  {
205  return magb < ScalarVSMALL ? 0.0 : magb*sqrt(1.0 + sqr(maga/magb));
206  }
207 }
208 
209 
210 inline Scalar pow3(const Scalar s)
211 {
212  return s*sqr(s);
213 }
214 
215 
216 inline Scalar pow4(const Scalar s)
217 {
218  return sqr(sqr(s));
219 }
220 
221 
222 inline Scalar pow5(const Scalar s)
223 {
224  return s*pow4(s);
225 }
226 
227 
228 inline Scalar pow6(const Scalar s)
229 {
230  return pow3(sqr(s));
231 }
232 
233 
234 inline Scalar pow025(const Scalar s)
235 {
236  return sqrt(sqrt(s));
237 }
238 
239 
240 inline Scalar inv(const Scalar s)
241 {
242  return 1.0/s;
243 }
244 
245 
246 inline Scalar dot(const Scalar s1, const Scalar s2)
247 {
248  return s1*s2;
249 }
250 
251 
252 inline Scalar cmptMultiply(const Scalar s1, const Scalar s2)
253 {
254  return s1*s2;
255 }
256 
257 
258 inline Scalar cmptPow(const Scalar s1, const Scalar s2)
259 {
260  return pow(s1, s2);
261 }
262 
263 
264 inline Scalar cmptDivide(const Scalar s1, const Scalar s2)
265 {
266  return s1/s2;
267 }
268 
269 
270 inline Scalar cmptMax(const Scalar s)
271 {
272  return s;
273 }
274 
275 
276 inline Scalar cmptMin(const Scalar s)
277 {
278  return s;
279 }
280 
281 
282 inline Scalar cmptAv(const Scalar s)
283 {
284  return s;
285 }
286 
287 
288 inline Scalar cmptMag(const Scalar s)
289 {
290  return mag(s);
291 }
292 
293 
294 // Standard C++ transcendental functions
297 transFunc(exp)
298 transFunc(log)
300 transFunc(sin)
301 transFunc(cos)
302 transFunc(tan)
312 
313 // Standard ANSI-C (but not in <cmath>) transcendental functions
314 
315 transFunc(erf)
318 
319 transFunc(j0)
320 transFunc(j1)
321 
322 transFunc(y0)
323 transFunc(y1)
324 
325 
326 // Stabilisation around zero for division
327 inline Scalar stabilise(const Scalar s, const Scalar small)
328 {
329  if (s >= 0)
330  {
331  return s + small;
332  }
333  else
334  {
335  return s - small;
336  }
337 }
338 
339 
340 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
341 
344 Ostream& operator<<(Ostream&, const Scalar);
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 } // End namespace Foam
349 
350 // ************************************************************************* //
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:164
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
uint8_t direction
Definition: direction.H:46
dimensionedScalar log(const dimensionedScalar &ds)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
#define ScalarVSMALL
Definition: doubleScalar.C:35
static const char *const typeName
Definition: Scalar.H:71
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
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(j0) transFunc(j1) transFunc(y0) transFunc(y1) inline Scalar stabilise(const Scalar s
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
Traits class for primitives.
Definition: pTraits.H:50
const Scalar small
Definition: Scalar.H:328
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)
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:209
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
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
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:194
static const Scalar min
Definition: Scalar.H:76
dimensionedScalar atanh(const dimensionedScalar &ds)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
dimensionedScalar y1(const dimensionedScalar &ds)
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:258
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:53
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
Scalar cmptType
Component type.
Definition: Scalar.H:51
#define Scalar
Definition: doubleScalar.C:33
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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:176
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:79
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)