dimensionSet.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) 2017-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 "dimensionSet.H"
30 #include "dimensionedScalar.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(dimensionSet, 1);
37 }
38 
39 const Foam::scalar Foam::dimensionSet::smallExponent = SMALL;
40 
41 
42 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 static inline bool checkDims
48 (
49  const char* what,
50  const dimensionSet& a,
51  const dimensionSet& b
52 )
53 {
54  if (a != b)
55  {
57  << "Different dimensions for '" << what
58  << "'\n dimensions : " << a << " != " << b << nl
59  << abort(FatalError);
60  return false;
61  }
62 
63  return true;
64 }
65 
66 } // End namespace Foam
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 :
73  exponents_(Zero)
74 {}
75 
76 
78 (
79  const scalar mass,
80  const scalar length,
81  const scalar time,
82  const scalar temperature,
83  const scalar moles,
84  const scalar current,
85  const scalar luminousIntensity
86 )
87 :
88  exponents_()
89 {
90  exponents_[MASS] = mass;
91  exponents_[LENGTH] = length;
92  exponents_[TIME] = time;
93  exponents_[TEMPERATURE] = temperature;
94  exponents_[MOLES] = moles;
95  exponents_[CURRENT] = current;
96  exponents_[LUMINOUS_INTENSITY] = luminousIntensity;
97 }
98 
99 
101 :
102  exponents_(dims)
103 {}
104 
105 
107 :
108  exponents_(ds.exponents_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
115 {
116  for (const scalar& val : exponents_)
117  {
118  // ie, mag(val) > smallExponent
119  if ((val > smallExponent) || (val < -smallExponent))
120  {
121  return false;
122  }
123  }
124 
125  return true;
126 }
127 
128 
130 {
131  return exponents_;
132 }
133 
134 
136 {
137  return exponents_;
138 }
139 
140 
142 {
143  exponents_ = Zero;
144 }
145 
146 
148 {
149  exponents_ = ds.exponents_;
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
154 
156 {
157  return exponents_[type];
158 }
159 
160 
162 {
163  return exponents_[type];
164 }
165 
166 
167 Foam::scalar Foam::dimensionSet::operator[](const label type) const
168 {
169  return exponents_[type];
170 }
171 
172 
173 Foam::scalar& Foam::dimensionSet::operator[](const label type)
174 {
175  return exponents_[type];
176 }
177 
178 
180 {
181  for (int d=0; d<nDimensions; ++d)
182  {
183  if
184  (
185  mag(exponents_[d] - ds.exponents_[d])
186  > smallExponent
187  )
188  {
189  return false;
190  }
191  }
192 
193  return true;
194 }
195 
196 
198 {
199  return !(operator==(ds));
200 }
201 
202 
204 {
206  {
207  checkDims("(a = b)", *this, ds);
208  }
209 
210  return true;
211 }
212 
213 
215 {
217  {
218  checkDims("(a += b)", *this, ds);
219  }
220 
221  return true;
222 }
223 
224 
226 {
228  {
229  checkDims("(a -= b)", *this, ds);
230  }
231 
232  return true;
233 }
234 
235 
237 {
238  reset((*this)*ds);
239 
240  return true;
241 }
242 
243 
245 {
246  reset((*this)/ds);
247 
248  return true;
249 }
250 
251 
252 // * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * * //
253 
255 {
257  {
258  checkDims("min(a, b)", ds1, ds2);
259  }
260 
261  return ds1;
262 }
263 
264 
266 {
268  {
269  checkDims("max(a, b)", ds1, ds2);
270  }
271 
272  return ds1;
273 }
274 
275 
277 {
279  {
280  checkDims("clip(a, b)", ds1, ds2);
281  }
282 
283  return ds1;
284 }
285 
286 
288 (
289  const dimensionSet& ds1,
290  const dimensionSet& ds2
291 )
292 {
293  return ds1*ds2;
294 }
295 
296 
298 (
299  const dimensionSet& ds1,
300  const dimensionSet& ds2
301 )
302 {
303  return ds1/ds2;
304 }
305 
306 
307 Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
308 {
309  return dimensionSet
310  (
311  ds[dimensionSet::MASS]*p,
313  ds[dimensionSet::TIME]*p,
318  );
319 }
320 
321 
323 (
324  const dimensionSet& ds,
325  const dimensionedScalar& dS
326 )
327 {
329  {
331  << "Exponent of pow is not dimensionless" << endl
332  << abort(FatalError);
333  }
334 
335  return pow(ds, dS.value());
336 }
337 
338 
340 (
341  const dimensionedScalar& dS,
342  const dimensionSet& ds
343 )
344 {
345  if
346  (
348  && !dS.dimensions().dimensionless()
349  && !ds.dimensionless()
350  )
351  {
353  << "Argument or exponent of pow not dimensionless" << endl
354  << abort(FatalError);
355  }
356 
357  return ds;
358 }
359 
360 
362 {
363  return pow(ds, 2);
364 }
365 
366 
368 {
369  return pow(ds, 2);
370 }
371 
372 
374 {
375  return pow(ds, 3);
376 }
377 
378 
380 {
381  return pow(ds, 4);
382 }
383 
384 
386 {
387  return pow(ds, 5);
388 }
389 
390 
392 {
393  return pow(ds, 6);
394 }
395 
396 
398 {
399  return pow(ds, 0.25);
400 }
401 
402 
404 {
405  return pow(ds, 0.5);
406 }
407 
408 
410 {
411  return pow(ds, 1.0/3.0);
412 }
413 
414 
416 {
417  return pow(ds, 2);
418 }
419 
420 
422 {
423  return ds;
424 }
425 
426 
428 {
429  return dimless;
430 }
431 
432 
434 {
435  return dimless;
436 }
437 
438 
440 {
441  return dimless;
442 }
443 
444 
446 {
447  return dimless;
448 }
449 
450 
452 {
453  return dimless;
454 }
455 
456 
458 {
459  return ds;
460 }
461 
462 
464 {
465  return ds;
466 }
467 
468 
470 {
471  return dimensionSet
472  (
473  0.0-ds[dimensionSet::MASS],
474  0.0-ds[dimensionSet::LENGTH],
475  0.0-ds[dimensionSet::TIME],
477  0.0-ds[dimensionSet::MOLES],
478  0.0-ds[dimensionSet::CURRENT],
480  );
481 }
482 
483 
485 {
486  if (dimensionSet::debug && !ds.dimensionless())
487  {
489  << "Argument of trancendental function not dimensionless" << nl
490  << abort(FatalError);
491  }
492 
493  return ds;
494 }
495 
496 
498 {
500  {
501  checkDims("atan2(a, b)", ds1, ds2);
502  }
503 
504  return dimless;
505 }
506 
507 
509 {
511  {
512  checkDims("hypot(a, b)", ds1, ds2);
513  }
514 
515  return ds1;
516 }
517 
518 
520 {
521  return ds;
522 }
523 
524 
526 {
527  return ds;
528 }
529 
530 
531 // * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
532 
534 {
535  return inv(ds);
536 }
537 
538 
540 {
541  return ds;
542 }
543 
544 
545 Foam::dimensionSet Foam::operator+
546 (
547  const dimensionSet& ds1,
548  const dimensionSet& ds2
549 )
550 {
552  {
553  checkDims("(a + b)", ds1, ds2);
554  }
555 
556  return ds1;
557 }
558 
559 
560 Foam::dimensionSet Foam::operator-
561 (
562  const dimensionSet& ds1,
563  const dimensionSet& ds2
564 )
565 {
567  {
568  checkDims("(a - b)", ds1, ds2);
569  }
570 
571  return ds1;
572 }
573 
574 
575 Foam::dimensionSet Foam::operator*
576 (
577  const dimensionSet& ds1,
578  const dimensionSet& ds2
579 )
580 {
581  dimensionSet result(ds1);
582 
583  auto rhs = ds2.values().begin();
584 
585  for (scalar& val : result.values())
586  {
587  val += *rhs;
588  ++rhs;
589  }
590 
591  return result;
592 }
593 
594 
595 Foam::dimensionSet Foam::operator/
596 (
597  const dimensionSet& ds1,
598  const dimensionSet& ds2
599 )
600 {
601  dimensionSet result(ds1);
602 
603  auto rhs = ds2.values().begin();
604 
605  for (scalar& val : result.values())
606  {
607  val -= *rhs;
608  ++rhs;
609  }
610 
611  return result;
612 }
613 
614 
615 Foam::dimensionSet Foam::operator&
616 (
617  const dimensionSet& ds1,
618  const dimensionSet& ds2
619 )
620 {
621  return ds1*ds2;
622 }
623 
624 
625 Foam::dimensionSet Foam::operator^
626 (
627  const dimensionSet& ds1,
628  const dimensionSet& ds2
629 )
630 {
631  return ds1*ds2;
632 }
633 
634 
635 Foam::dimensionSet Foam::operator&&
636 (
637  const dimensionSet& ds1,
638  const dimensionSet& ds2
639 )
640 {
641  return ds1*ds2;
642 }
643 
644 
645 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::trans
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions)
Definition: dimensionSet.C:484
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::dimensionSet::operator+=
bool operator+=(const dimensionSet &) const
Definition: dimensionSet.C:214
Foam::FixedList::begin
iterator begin()
Return an iterator to begin traversing the FixedList.
Definition: FixedListI.H:506
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::dimensionSet::operator/=
bool operator/=(const dimensionSet &)
Definition: dimensionSet.C:244
Foam::pow2
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:367
Foam::dimensionSet::LUMINOUS_INTENSITY
Candela Cd.
Definition: dimensionSet.H:87
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimensionSet::smallExponent
static const scalar smallExponent
Tolerance for 'small' exponents, for near-zero rounding.
Definition: dimensionSet.H:94
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::dimensionSet::operator*=
bool operator*=(const dimensionSet &)
Definition: dimensionSet.C:236
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
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::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::dimensionSet::dimensionSet
dimensionSet()
Default construct (dimensionless).
Definition: dimensionSet.C:71
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
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::dimensionSet::dimensionType
dimensionType
Enumeration for the dimension exponents.
Definition: dimensionSet.H:79
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::dimensionSet::dimensionless
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:114
Foam::pow6
dimensionedScalar pow6(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:122
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::dimensionSet::TEMPERATURE
Kelvin K.
Definition: dimensionSet.H:84
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::dimensionSet::operator==
bool operator==(const dimensionSet &) const
Definition: dimensionSet.C:179
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::dimensionSet::clear
void clear()
Reset exponents to be dimensionless.
Definition: dimensionSet.C:141
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::dimensionSet::operator!=
bool operator!=(const dimensionSet &) const
Definition: dimensionSet.C:197
dimensionSet.H
Foam::dimensionSet::operator-=
bool operator-=(const dimensionSet &) const
Definition: dimensionSet.C:225
Foam::FatalError
error FatalError
Foam::pow5
dimensionedScalar pow5(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:111
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::dimensionSet::MASS
kilogram kg
Definition: dimensionSet.H:81
Foam::cmptDivide
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar.H
Foam::operator~
bitSet operator~(const bitSet &bitset)
Bitset complement, returns a copy of the bitset with all its bits flipped.
Definition: bitSetI.H:746
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::dimensionSet::values
const FixedList< scalar, 7 > & values() const
Return const access to the exponents as a list.
Definition: dimensionSet.C:129
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::negPart
dimensionedScalar negPart(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:232
Foam::clip
dimensionSet clip(const dimensionSet &ds1, const dimensionSet &ds2)
Definition: dimensionSet.C:276
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::dimensionSet::CURRENT
Ampere A.
Definition: dimensionSet.H:86
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::FixedList< scalar, 7 >
Foam::dimensionSet::operator[]
scalar operator[](const dimensionType) const
Definition: dimensionSet.C:155
Foam::invTransform
dimensionSet invTransform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:525
Foam::dimensionSet::TIME
second s
Definition: dimensionSet.H:83
Foam::dimensionSet::operator=
bool operator=(const dimensionSet &) const
Definition: dimensionSet.C:203
Foam::dimensionSet::LENGTH
metre m
Definition: dimensionSet.H:82
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::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::dimensionSet::MOLES
mole mol
Definition: dimensionSet.H:85
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::checkDims
static bool checkDims(const char *what, const dimensionSet &a, const dimensionSet &b)
Definition: dimensionSet.C:48
Foam::dimensionSet::reset
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:147
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177