dimensionedScalar.C
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) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "dimensionedScalar.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 dimensionedScalar operator+(const dimensionedScalar& ds1, const scalar s2)
39 {
40  return ds1 + dimensionedScalar(s2);
41 }
42 
43 
44 dimensionedScalar operator+(const scalar s1, const dimensionedScalar& ds2)
45 {
46  return dimensionedScalar(s1) + ds2;
47 }
48 
49 
50 dimensionedScalar operator-(const dimensionedScalar& ds1, const scalar s2)
51 {
52  return ds1 - dimensionedScalar(s2);
53 }
54 
55 
56 dimensionedScalar operator-(const scalar s1, const dimensionedScalar& ds2)
57 {
58  return dimensionedScalar(s1) - ds2;
59 }
60 
61 
62 dimensionedScalar operator*(const dimensionedScalar& ds1, const scalar s2)
63 {
64  return ds1 * dimensionedScalar(s2);
65 }
66 
67 
68 dimensionedScalar operator/(const scalar s1, const dimensionedScalar& ds2)
69 {
70  return dimensionedScalar(s1)/ds2;
71 }
72 
73 
75 (
76  const dimensionedScalar& ds,
77  const dimensionedScalar& expt
78 )
79 {
80  return dimensionedScalar
81  (
82  "pow(" + ds.name() + ',' + expt.name() + ')',
83  pow(ds.dimensions(), expt),
84  ::pow(ds.value(), expt.value())
85  );
86 }
87 
88 
90 {
91  return dimensionedScalar
92  (
93  "pow3(" + ds.name() + ')',
94  pow3(ds.dimensions()),
95  pow3(ds.value())
96  );
97 }
98 
99 
101 {
102  return dimensionedScalar
103  (
104  "pow4(" + ds.name() + ')',
105  pow4(ds.dimensions()),
106  pow4(ds.value())
107  );
108 }
109 
110 
112 {
113  return dimensionedScalar
114  (
115  "pow5(" + ds.name() + ')',
116  pow5(ds.dimensions()),
117  pow5(ds.value())
118  );
119 }
120 
121 
123 {
124  return dimensionedScalar
125  (
126  "pow6(" + ds.name() + ')',
127  pow6(ds.dimensions()),
128  pow6(ds.value())
129  );
130 }
131 
132 
134 {
135  return dimensionedScalar
136  (
137  "pow025(" + ds.name() + ')',
138  pow025(ds.dimensions()),
139  pow025(ds.value())
140  );
141 }
142 
143 
145 {
146  return dimensionedScalar
147  (
148  "sqrt(" + ds.name() + ')',
149  pow(ds.dimensions(), dimensionedScalar("0.5", dimless, 0.5)),
150  ::sqrt(ds.value())
151  );
152 }
153 
154 
156 {
157  return dimensionedScalar
158  (
159  "cbrt(" + ds.name() + ')',
160  pow(ds.dimensions(), dimensionedScalar("(1|3)", dimless, 1.0/3.0)),
161  ::cbrt(ds.value())
162  );
163 }
164 
165 
167 {
168  return dimensionedScalar
169  (
170  "sign(" + ds.name() + ')',
171  sign(ds.dimensions()),
172  ::Foam::sign(ds.value())
173  );
174 }
175 
176 
178 {
179  return dimensionedScalar
180  (
181  "pos(" + ds.name() + ')',
182  pos(ds.dimensions()),
183  ::Foam::pos(ds.value())
184  );
185 }
186 
187 
189 {
190  return dimensionedScalar
191  (
192  "pos0(" + ds.name() + ')',
193  pos0(ds.dimensions()),
194  ::Foam::pos0(ds.value())
195  );
196 }
197 
198 
200 {
201  return dimensionedScalar
202  (
203  "neg(" + ds.name() + ')',
204  neg(ds.dimensions()),
205  ::Foam::neg(ds.value())
206  );
207 }
208 
209 
211 {
212  return dimensionedScalar
213  (
214  "neg0(" + ds.name() + ')',
215  neg0(ds.dimensions()),
216  ::Foam::neg0(ds.value())
217  );
218 }
219 
220 
222 {
223  return dimensionedScalar
224  (
225  "posPart(" + ds.name() + ')',
226  posPart(ds.dimensions()),
227  ::Foam::pos0(ds.value())
228  );
229 }
230 
231 
233 {
234  return dimensionedScalar
235  (
236  "negPart(" + ds.name() + ')',
237  negPart(ds.dimensions()),
238  ::Foam::neg(ds.value())
239  );
240 }
241 
242 
243 #define transFunc(func) \
244 dimensionedScalar func(const dimensionedScalar& ds) \
245 { \
246  if (!ds.dimensions().dimensionless()) \
247  { \
248  FatalErrorInFunction \
249  << "ds not dimensionless" \
250  << abort(FatalError); \
251  } \
252  \
253  return dimensionedScalar \
254  ( \
255  #func "(" + ds.name() + ')', \
256  dimless, \
257  ::func(ds.value()) \
258  ); \
259 }
260 
283 
284 #undef transFunc
285 
286 
287 #define transFunc(func) \
288 dimensionedScalar func(const int n, const dimensionedScalar& ds) \
289 { \
290  if (!ds.dimensions().dimensionless()) \
291  { \
292  FatalErrorInFunction \
293  << "ds not dimensionless" \
294  << abort(FatalError); \
295  } \
296  \
297  return dimensionedScalar \
298  ( \
299  #func "(" + name(n) + ',' + ds.name() + ')', \
300  dimless, \
301  ::func(n, ds.value()) \
302  ); \
303 }
304 
307 
308 #undef transFunc
309 
310 
312 (
313  const dimensionedScalar& x,
314  const dimensionedScalar& y
315 )
316 {
317  return dimensionedScalar
318  (
319  "atan2(" + x.name() + ',' + y.name() + ')',
320  atan2(x.dimensions(), y.dimensions()),
321  ::atan2(x.value(), y.value())
322  );
323 }
324 
325 
327 (
328  const dimensionedScalar& x,
329  const dimensionedScalar& y
330 )
331 {
332  return dimensionedScalar
333  (
334  "hypot(" + x.name() + ',' + y.name() + ')',
335  hypot(x.dimensions(), y.dimensions()),
336  ::hypot(x.value(), y.value())
337  );
338 }
339 
340 
341 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
342 
343 } // End namespace Foam
344 
345 // ************************************************************************* //
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::cosh
dimensionedScalar cosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:271
Foam::y1
dimensionedScalar y1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:282
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::jn
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:305
Foam::posPart
dimensionedScalar posPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:221
Foam::atan2
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:312
Foam::operator-
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Foam::neg0
dimensionedScalar neg0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:210
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::erf
dimensionedScalar erf(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:276
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::hypot
dimensionedScalar hypot(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:327
Foam::atanh
dimensionedScalar atanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:275
Foam::lgamma
dimensionedScalar lgamma(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:278
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:122
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::log10
dimensionedScalar log10(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::y0
dimensionedScalar y0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:281
Foam::erfc
dimensionedScalar erfc(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:277
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::asinh
dimensionedScalar asinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:273
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:111
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::operator/
dimensionedScalar operator/(const scalar s1, const dimensionedScalar &ds2)
Definition: dimensionedScalar.C:68
Foam::yn
dimensionedScalar yn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:306
dimensionedScalar.H
Foam::acosh
dimensionedScalar acosh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:274
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:232
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::j0
dimensionedScalar j0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:279
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::operator+
tmp< faMatrix< Type > > operator+(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::operator*
tmp< faMatrix< Type > > operator*(const areaScalarField &, const faMatrix< Type > &)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::transFunc
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 &s
Foam::atan
dimensionedScalar atan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:269
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:420
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::j1
dimensionedScalar j1(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:280
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::asin
dimensionedScalar asin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:267
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177
Foam::sinh
dimensionedScalar sinh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:270