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-2015 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  // Member constants
54 
55  enum
56  {
57  dim = 3,
58  rank = 0,
59  nComponents = 1
60  };
61 
62  // Static data members
63 
64  static const char* const typeName;
65  static const char* componentNames[];
66  static const Scalar zero;
67  static const Scalar one;
68  static const Scalar max;
69  static const Scalar min;
70  static const Scalar rootMax;
71  static const Scalar rootMin;
72 
73  // Constructors
74 
75  //- Construct from primitive
76  explicit pTraits(const Scalar&);
77 
78  //- Construct from Istream
79  pTraits(Istream&);
80 
81 
82  // Member Functions
83 
84  //- Access to the Scalar value
85  operator Scalar() const
86  {
87  return p_;
88  }
89 
90  //- Access to the Scalar value
91  operator Scalar&()
92  {
93  return p_;
94  }
95 };
96 
97 
98 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100 //- Return a string representation of a Scalar
101 word name(const Scalar);
102 
103 
105 {
106  return s;
107 }
108 
109 
110 inline Scalar component(const Scalar s, const direction)
111 {
112  return s;
113 }
114 
115 
116 inline Scalar sign(const Scalar s)
117 {
118  return (s >= 0)? 1: -1;
119 }
120 
121 
122 inline Scalar pos(const Scalar s)
123 {
124  return (s >= 0)? 1: 0;
125 }
126 
127 
128 inline Scalar neg(const Scalar s)
129 {
130  return (s < 0)? 1: 0;
131 }
132 
133 
134 inline Scalar posPart(const Scalar s)
135 {
136  return (s > 0)? s: 0;
137 }
138 
139 
140 inline Scalar negPart(const Scalar s)
141 {
142  return (s < 0)? s: 0;
143 }
144 
145 
146 inline bool equal(const Scalar& s1, const Scalar& s2)
147 {
148  return mag(s1 - s2) <= ScalarVSMALL;
149 }
150 
151 
152 inline bool notEqual(const Scalar s1, const Scalar s2)
153 {
154  return mag(s1 - s2) > ScalarVSMALL;
155 }
156 
157 
158 inline Scalar limit(const Scalar s1, const Scalar s2)
159 {
160  return (mag(s1) < mag(s2)) ? s1: 0.0;
161 }
162 
163 
164 inline Scalar minMod(const Scalar s1, const Scalar s2)
165 {
166  return (mag(s1) < mag(s2)) ? s1: s2;
167 }
168 
169 
170 inline Scalar magSqr(const Scalar s)
171 {
172  return s*s;
173 }
174 
175 
176 inline Scalar sqr(const Scalar s)
177 {
178  return s*s;
179 }
180 
181 
182 inline Scalar sqrtSumSqr(const Scalar a, const Scalar b)
183 {
184  Scalar maga = mag(a);
185  Scalar magb = mag(b);
186 
187  if (maga > magb)
188  {
189  return maga*sqrt(1.0 + sqr(magb/maga));
190  }
191  else
192  {
193  return magb < ScalarVSMALL ? 0.0 : magb*sqrt(1.0 + sqr(maga/magb));
194  }
195 }
196 
197 
198 inline Scalar pow3(const Scalar s)
199 {
200  return s*sqr(s);
201 }
202 
203 
204 inline Scalar pow4(const Scalar s)
205 {
206  return sqr(sqr(s));
207 }
208 
209 
210 inline Scalar pow5(const Scalar s)
211 {
212  return s*pow4(s);
213 }
214 
215 
216 inline Scalar pow6(const Scalar s)
217 {
218  return pow3(sqr(s));
219 }
220 
221 
222 inline Scalar pow025(const Scalar s)
223 {
224  return sqrt(sqrt(s));
225 }
226 
227 
228 inline Scalar inv(const Scalar s)
229 {
230  return 1.0/s;
231 }
232 
233 
234 inline Scalar dot(const Scalar s1, const Scalar s2)
235 {
236  return s1*s2;
237 }
238 
239 
240 inline Scalar cmptMultiply(const Scalar s1, const Scalar s2)
241 {
242  return s1*s2;
243 }
244 
245 
246 inline Scalar cmptPow(const Scalar s1, const Scalar s2)
247 {
248  return pow(s1, s2);
249 }
250 
251 
252 inline Scalar cmptDivide(const Scalar s1, const Scalar s2)
253 {
254  return s1/s2;
255 }
256 
257 
258 inline Scalar cmptMax(const Scalar s)
259 {
260  return s;
261 }
262 
263 
264 inline Scalar cmptMin(const Scalar s)
265 {
266  return s;
267 }
268 
269 
270 inline Scalar cmptAv(const Scalar s)
271 {
272  return s;
273 }
274 
275 
276 inline Scalar cmptMag(const Scalar s)
277 {
278  return mag(s);
279 }
280 
281 
282 // Standard C++ transcendental functions
285 transFunc(exp)
286 transFunc(log)
288 transFunc(sin)
289 transFunc(cos)
290 transFunc(tan)
300 
301 // Standard ANSI-C (but not in <cmath>) transcendental functions
302 
303 transFunc(erf)
306 
307 transFunc(j0)
308 transFunc(j1)
309 
310 transFunc(y0)
311 transFunc(y1)
312 
313 
314 // Stabilisation around zero for division
315 inline Scalar stabilise(const Scalar s, const Scalar small)
316 {
317  if (s >= 0)
318  {
319  return s + small;
320  }
321  else
322  {
323  return s - small;
324  }
325 }
326 
327 
328 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
329 
332 Ostream& operator<<(Ostream&, const Scalar);
333 
334 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
335 
336 } // End namespace Foam
337 
338 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
Scalar sqrtSumSqr(const Scalar a, const Scalar b)
Definition: Scalar.H:182
dimensionedScalar pow3(const dimensionedScalar &ds)
unsigned char direction
Definition: direction.H:43
dimensionedScalar erfc(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
Scalar cmptPow(const Scalar s1, const Scalar s2)
Definition: Scalar.H:246
dimensionedScalar j0(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar neg(const dimensionedScalar &ds)
pTraits(const PrimitiveType &p)
Construct from primitive.
Definition: pTraits.H:60
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 ))
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar j1(const dimensionedScalar &ds)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Scalar minMod(const Scalar s1, const Scalar s2)
Definition: Scalar.H:164
static const Scalar zero
Definition: Scalar.H:66
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const Scalar one
Definition: Scalar.H:67
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
dimensionedScalar posPart(const dimensionedScalar &ds)
static const Scalar rootMin
Definition: Scalar.H:71
label & setComponent(label &l, const direction)
Definition: label.H:79
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar log10(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
Scalar cmptType
Component type.
Definition: Scalar.H:51
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
dimensionedScalar y0(const dimensionedScalar &ds)
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:133
Namespace for OpenFOAM.
dimensionedScalar pow6(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)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
dimensionedScalar negPart(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
static const char *const typeName
Definition: Scalar.H:64
dimensionedScalar asinh(const dimensionedScalar &ds)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dimensionedScalar lgamma(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
#define Scalar
Definition: doubleScalar.C:33
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedScalar atan(const dimensionedScalar &ds)
tmp< DimensionedField< scalar, GeoMesh > > stabilise(const DimensionedField< scalar, GeoMesh > &dsf, const dimensioned< scalar > &ds)
bool notEqual(const Scalar s1, const Scalar s2)
Definition: Scalar.H:152
Traits class for primitives.
Definition: pTraits.H:50
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Istream & operator>>(Istream &, edgeMesh &)
Definition: edgeMeshIO.C:144
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
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
complex limit(const complex &, const complex &)
Definition: complexI.H:211
static const Scalar rootMax
Definition: Scalar.H:70
#define ScalarVSMALL
Definition: doubleScalar.C:35
dimensionedScalar pos(const dimensionedScalar &ds)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
const Scalar small
Definition: Scalar.H:316
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar acosh(const dimensionedScalar &ds)
dimensionedScalar cosh(const dimensionedScalar &ds)
dimensionedScalar asin(const dimensionedScalar &ds)
bool equal(const T &s1, const T &s2)
Definition: doubleFloat.H:62
dimensionedScalar atanh(const dimensionedScalar &ds)
static const Scalar max
Definition: Scalar.H:68
static const Scalar min
Definition: Scalar.H:69
dimensionedScalar sinh(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)