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