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-2021 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 
131 {
132  return exponents_;
133 }
134 
135 
138 {
139  return exponents_;
140 }
141 
142 
144 {
145  exponents_ = Zero;
146 }
147 
148 
150 {
151  exponents_ = ds.exponents_;
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
156 
158 {
159  return exponents_[type];
160 }
161 
162 
164 {
165  return exponents_[type];
166 }
167 
168 
169 Foam::scalar Foam::dimensionSet::operator[](const label type) const
170 {
171  return exponents_[type];
172 }
173 
174 
175 Foam::scalar& Foam::dimensionSet::operator[](const label type)
176 {
177  return exponents_[type];
178 }
179 
180 
182 {
183  for (int d=0; d<nDimensions; ++d)
184  {
185  if
186  (
187  mag(exponents_[d] - ds.exponents_[d])
188  > smallExponent
189  )
190  {
191  return false;
192  }
193  }
194 
195  return true;
196 }
197 
198 
200 {
201  return !(operator==(ds));
202 }
203 
204 
206 {
208  {
209  checkDims("(a = b)", *this, ds);
210  }
211 
212  return true;
213 }
214 
215 
217 {
219  {
220  checkDims("(a += b)", *this, ds);
221  }
222 
223  return true;
224 }
225 
226 
228 {
230  {
231  checkDims("(a -= b)", *this, ds);
232  }
233 
234  return true;
235 }
236 
237 
239 {
240  reset((*this)*ds);
241 
242  return true;
243 }
244 
245 
247 {
248  reset((*this)/ds);
249 
250  return true;
251 }
252 
253 
254 // * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * * //
255 
257 {
259  {
260  checkDims("min(a, b)", ds1, ds2);
261  }
262 
263  return ds1;
264 }
265 
266 
268 {
270  {
271  checkDims("max(a, b)", ds1, ds2);
272  }
273 
274  return ds1;
275 }
276 
277 
279 {
281  {
282  checkDims("clip(a, b)", ds1, ds2);
283  }
284 
285  return ds1;
286 }
287 
288 
290 (
291  const dimensionSet& ds1,
292  const dimensionSet& ds2
293 )
294 {
295  return ds1*ds2;
296 }
297 
298 
300 (
301  const dimensionSet& ds1,
302  const dimensionSet& ds2
303 )
304 {
305  return ds1/ds2;
306 }
307 
308 
309 Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
310 {
311  return dimensionSet
312  (
313  ds[dimensionSet::MASS]*p,
315  ds[dimensionSet::TIME]*p,
320  );
321 }
322 
323 
325 (
326  const dimensionSet& ds,
327  const dimensionedScalar& dS
328 )
329 {
331  {
333  << "Exponent of pow is not dimensionless" << endl
334  << abort(FatalError);
335  }
336 
337  return pow(ds, dS.value());
338 }
339 
340 
342 (
343  const dimensionedScalar& dS,
344  const dimensionSet& ds
345 )
346 {
347  if
348  (
350  && !dS.dimensions().dimensionless()
351  && !ds.dimensionless()
352  )
353  {
355  << "Argument or exponent of pow not dimensionless" << endl
356  << abort(FatalError);
357  }
358 
359  return ds;
360 }
361 
362 
364 {
365  return pow(ds, 2);
366 }
367 
368 
370 {
371  return pow(ds, 2);
372 }
373 
374 
376 {
377  return pow(ds, 3);
378 }
379 
380 
382 {
383  return pow(ds, 4);
384 }
385 
386 
388 {
389  return pow(ds, 5);
390 }
391 
392 
394 {
395  return pow(ds, 6);
396 }
397 
398 
400 {
401  return pow(ds, 0.25);
402 }
403 
404 
406 {
407  return pow(ds, 0.5);
408 }
409 
410 
412 {
413  return pow(ds, 1.0/3.0);
414 }
415 
416 
418 {
419  return pow(ds, 2);
420 }
421 
422 
424 {
425  return ds;
426 }
427 
428 
430 {
431  return dimless;
432 }
433 
434 
436 {
437  return dimless;
438 }
439 
440 
442 {
443  return dimless;
444 }
445 
446 
448 {
449  return dimless;
450 }
451 
452 
454 {
455  return dimless;
456 }
457 
458 
460 {
461  return ds;
462 }
463 
464 
466 {
467  return ds;
468 }
469 
470 
472 {
473  return dimensionSet
474  (
475  0.0-ds[dimensionSet::MASS],
476  0.0-ds[dimensionSet::LENGTH],
477  0.0-ds[dimensionSet::TIME],
479  0.0-ds[dimensionSet::MOLES],
480  0.0-ds[dimensionSet::CURRENT],
482  );
483 }
484 
485 
487 {
488  if (dimensionSet::checking() && !ds.dimensionless())
489  {
491  << "Argument of trancendental function not dimensionless" << nl
492  << abort(FatalError);
493  }
494 
495  return ds;
496 }
497 
498 
500 {
502  {
503  checkDims("atan2(a, b)", ds1, ds2);
504  }
505 
506  return dimless;
507 }
508 
509 
511 {
513  {
514  checkDims("hypot(a, b)", ds1, ds2);
515  }
516 
517  return ds1;
518 }
519 
520 
522 {
523  return ds;
524 }
525 
526 
528 {
529  return ds;
530 }
531 
532 
533 // * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
534 
536 {
537  return inv(ds);
538 }
539 
540 
542 {
543  return ds;
544 }
545 
546 
547 Foam::dimensionSet Foam::operator+
548 (
549  const dimensionSet& ds1,
550  const dimensionSet& ds2
551 )
552 {
554  {
555  checkDims("(a + b)", ds1, ds2);
556  }
557 
558  return ds1;
559 }
560 
561 
562 Foam::dimensionSet Foam::operator-
563 (
564  const dimensionSet& ds1,
565  const dimensionSet& ds2
566 )
567 {
569  {
570  checkDims("(a - b)", ds1, ds2);
571  }
572 
573  return ds1;
574 }
575 
576 
577 Foam::dimensionSet Foam::operator*
578 (
579  const dimensionSet& ds1,
580  const dimensionSet& ds2
581 )
582 {
583  dimensionSet result(ds1);
584 
585  auto rhs = ds2.values().begin();
586 
587  for (scalar& val : result.values())
588  {
589  val += *rhs;
590  ++rhs;
591  }
592 
593  return result;
594 }
595 
596 
597 Foam::dimensionSet Foam::operator/
598 (
599  const dimensionSet& ds1,
600  const dimensionSet& ds2
601 )
602 {
603  dimensionSet result(ds1);
604 
605  auto rhs = ds2.values().begin();
606 
607  for (scalar& val : result.values())
608  {
609  val -= *rhs;
610  ++rhs;
611  }
612 
613  return result;
614 }
615 
616 
617 Foam::dimensionSet Foam::operator&
618 (
619  const dimensionSet& ds1,
620  const dimensionSet& ds2
621 )
622 {
623  return ds1*ds2;
624 }
625 
626 
627 Foam::dimensionSet Foam::operator^
628 (
629  const dimensionSet& ds1,
630  const dimensionSet& ds2
631 )
632 {
633  return ds1*ds2;
634 }
635 
636 
637 Foam::dimensionSet Foam::operator&&
638 (
639  const dimensionSet& ds1,
640  const dimensionSet& ds2
641 )
642 {
643  return ds1*ds2;
644 }
645 
646 
647 // ************************************************************************* //
Foam::trans
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions)
Definition: dimensionSet.C:486
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dimensionSet::operator+=
bool operator+=(const dimensionSet &) const
Definition: dimensionSet.C:216
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::dimensionSet::operator/=
bool operator/=(const dimensionSet &)
Definition: dimensionSet.C:246
Foam::pow2
dimensionSet pow2(const dimensionSet &ds)
Definition: dimensionSet.C:369
Foam::dimensionSet::LUMINOUS_INTENSITY
Candela cd.
Definition: dimensionSet.H:130
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:137
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:238
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
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, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
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:521
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::FixedList::begin
iterator begin() noexcept
Return an iterator to begin traversing the FixedList.
Definition: FixedListI.H:524
Foam::dimensionSet::dimensionType
dimensionType
Enumeration for the dimension exponents.
Definition: dimensionSet.H:122
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:127
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:181
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::dimensionSet::clear
void clear()
Clear exponents - resets to be dimensionless.
Definition: dimensionSet.C:143
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:199
dimensionSet.H
Foam::dimensionSet::operator-=
bool operator-=(const dimensionSet &) const
Definition: dimensionSet.C:227
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:124
Foam::dimensionSet::values
const FixedList< scalar, 7 > & values() const noexcept
Const access to the exponents as a list.
Definition: dimensionSet.C:130
reset
meshPtr reset(new Foam::fvMesh(Foam::IOobject(regionName, runTime.timeName(), runTime, Foam::IOobject::MUST_READ), false))
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:453
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:278
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::dimensionSet::CURRENT
Ampere A.
Definition: dimensionSet.H:129
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:157
Foam::invTransform
dimensionSet invTransform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:527
Foam::dimensionSet::TIME
second s
Definition: dimensionSet.H:126
Foam::dimensionSet::operator=
bool operator=(const dimensionSet &) const
Definition: dimensionSet.C:205
Foam::dimensionSet::LENGTH
metre m
Definition: dimensionSet.H:125
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:128
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::checking
static bool checking() noexcept
True if dimension checking is enabled (the usual default)
Definition: dimensionSet.H:229
Foam::dimensionSet::reset
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:149
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177