Loading...
Searching...
No Matches
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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "dimensionedScalar.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37}
38
39const Foam::scalar Foam::dimensionSet::smallExponent = SMALL;
41
42// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
43
44namespace Foam
45{
46
47static 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;
65
66} // End namespace Foam
67
68
69// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72:
73 exponents_(Foam::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
101:
102 exponents_(dims)
103{}
104
105
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 }
124
125 return true;
126}
127
128
131{
132 return exponents_;
133}
134
135
140}
141
144{
145 exponents_ = Foam::zero{};
146}
147
148
151 exponents_ = ds.exponents_;
152}
153
154
155// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
157Foam::scalar Foam::dimensionSet::operator[](const dimensionType type) const
158{
159 return exponents_[type];
160}
161
164{
165 return exponents_[type];
166}
167
169Foam::scalar Foam::dimensionSet::operator[](const int type) const
170{
171 return exponents_[type];
172}
173
175Foam::scalar& Foam::dimensionSet::operator[](const int 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;
193 }
194
195 return true;
196}
197
200{
201 return !(operator==(ds));
202}
203
204
206{
208 {
209 checkDims("(a = b)", *this, ds);
210 }
211
212 return true;
213}
214
215
216bool Foam::dimensionSet::operator+=(const dimensionSet& ds) const
217{
219 {
220 checkDims("(a += b)", *this, ds);
221 }
222
223 return true;
224}
225
226
227bool Foam::dimensionSet::operator-=(const dimensionSet& ds) const
228{
230 {
231 checkDims("(a -= b)", *this, ds);
232 }
233
234 return true;
235}
236
237
238bool Foam::dimensionSet::operator*=(const dimensionSet& ds)
240 reset((*this)*ds);
241
242 return true;
243}
244
245
247{
248 reset((*this)/ds);
250 return true;
251}
252
253
254// * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * * //
255
256Foam::dimensionSet Foam::min(const dimensionSet& a, const dimensionSet& b)
257{
259 {
260 checkDims("min(a, b)", a, b);
261 }
262
263 return a;
264}
265
266
267Foam::dimensionSet Foam::max(const dimensionSet& a, const dimensionSet& b)
268{
270 {
271 checkDims("max(a, b)", a, b);
272 }
273
274 return a;
275}
276
277
278Foam::dimensionSet Foam::clamp(const dimensionSet& a, const dimensionSet& range)
279{
280 // In may cases the min/max range will be created from raw values
281 // (no dimension) so accept those without error
282 if (dimensionSet::checking() && !range.dimensionless())
284 checkDims("clamp(a, b)", a, range);
285 }
286
287 return a;
288}
289
290Foam::dimensionSet Foam::lerp(const dimensionSet& a, const dimensionSet& b)
291{
293 {
294 checkDims("lerp(a, b)", a, b);
295 }
296
297 return a;
298}
299
300
301Foam::dimensionSet Foam::cmptMultiply
302(
303 const dimensionSet& ds1,
304 const dimensionSet& ds2
305)
306{
307 return ds1*ds2;
308}
309
310
312(
313 const dimensionSet& ds1,
314 const dimensionSet& ds2
315)
316{
317 return ds1/ds2;
318}
319
320
321Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
322{
323 return dimensionSet
324 (
332 );
333}
334
335
337(
338 const dimensionSet& ds,
339 const dimensionedScalar& p
340)
341{
342 if (dimensionSet::checking() && !p.dimensions().dimensionless())
343 {
345 << "Exponent of pow is not dimensionless" << endl
347 }
348
349 return pow(ds, p.value());
350}
351
354{
355 return pow(ds, 2);
356}
357
360{
361 return pow(ds, 2);
362}
363
366{
367 return pow(ds, 3);
368}
369
372{
373 return pow(ds, 4);
374}
375
378{
379 return pow(ds, 5);
380}
381
384{
385 return pow(ds, 6);
386}
387
390{
391 return pow(ds, 0.25);
392}
393
396{
397 return pow(ds, 0.5);
398}
399
402{
403 return pow(ds, 1.0/3.0);
404}
405
408{
409 return pow(ds, 2);
410}
411
414{
415 return ds;
416}
417
422}
423
428}
429
434}
435
440}
441
446}
447
450{
451 return ds;
452}
453
456{
457 return ds;
458}
459
460
461Foam::dimensionSet Foam::inv(const dimensionSet& ds)
462{
463 return dimensionSet
464 (
465 // Avoid signed zero values by subtracting from zero (with rounding)
466 // instead of using a simple '-value'
467 0.0-ds[dimensionSet::MASS],
468 0.0-ds[dimensionSet::LENGTH],
469 0.0-ds[dimensionSet::TIME],
472 0.0-ds[dimensionSet::CURRENT],
474 );
475}
476
477
478Foam::dimensionSet Foam::trans(const dimensionSet& ds)
479{
480 if (dimensionSet::checking() && !ds.dimensionless())
481 {
483 << "Argument of trancendental function not dimensionless" << nl
485 }
486
487 return ds;
488}
489
490
491Foam::dimensionSet Foam::atan2(const dimensionSet& ds1, const dimensionSet& ds2)
492{
494 {
495 checkDims("atan2(a, b)", ds1, ds2);
496 }
497
498 return dimless;
499}
500
501
502Foam::dimensionSet Foam::hypot(const dimensionSet& ds1, const dimensionSet& ds2)
503{
505 {
506 checkDims("hypot(a, b)", ds1, ds2);
507 }
508
509 return ds1;
510}
511
512
513Foam::dimensionSet Foam::stabilise
514(
515 const dimensionSet& ds1,
516 const dimensionSet& ds2
517)
518{
520 {
521 checkDims("stabilise(a, b)", ds1, ds2);
522 }
523
524 return ds1;
525}
526
529{
530 return ds;
531}
532
533
534Foam::dimensionSet Foam::invTransform(const dimensionSet& ds)
536 return ds;
537}
538
539
540// * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
543{
544 return inv(ds);
545}
546
549{
550 return ds;
551}
552
553
554Foam::dimensionSet Foam::operator+
555(
556 const dimensionSet& ds1,
557 const dimensionSet& ds2
558)
559{
561 {
562 checkDims("(a + b)", ds1, ds2);
563 }
564
565 return ds1;
566}
567
568
569Foam::dimensionSet Foam::operator-
570(
571 const dimensionSet& ds1,
572 const dimensionSet& ds2
573)
574{
576 {
577 checkDims("(a - b)", ds1, ds2);
578 }
579
580 return ds1;
581}
582
583
584Foam::dimensionSet Foam::operator*
585(
586 const dimensionSet& ds1,
587 const dimensionSet& ds2
588)
589{
590 dimensionSet result(ds1);
591
592 auto rhs = ds2.values().begin();
593
594 for (scalar& val : result.values())
595 {
596 val += *rhs;
598 }
599
600 return result;
601}
602
603
604Foam::dimensionSet Foam::operator/
605(
606 const dimensionSet& ds1,
607 const dimensionSet& ds2
608)
609{
610 dimensionSet result(ds1);
611
612 auto rhs = ds2.values().begin();
613
614 for (scalar& val : result.values())
615 {
616 val -= *rhs;
618 }
619
620 return result;
621}
622
623
624Foam::dimensionSet Foam::operator&
625(
626 const dimensionSet& ds1,
627 const dimensionSet& ds2
628)
629{
630 return ds1*ds2;
631}
632
633
634Foam::dimensionSet Foam::operator^
635(
636 const dimensionSet& ds1,
637 const dimensionSet& ds2
638)
639{
640 return ds1*ds2;
641}
642
643
644Foam::dimensionSet Foam::operator&&
645(
646 const dimensionSet& ds1,
647 const dimensionSet& ds2
648)
649{
650 return ds1*ds2;
651}
652
653
654// ************************************************************************* //
scalar range
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
bool operator+=(const dimensionSet &) const
Addition operation, checks for identical dimensions.
bool operator==(const dimensionSet &) const
static constexpr int nDimensions
There are 7 base dimensions.
bool operator!=(const dimensionSet &) const
dimensionSet()
Default construct (dimensionless).
static bool checking(bool on) noexcept
Turn dimension checking on/off.
static bool checking() noexcept
True if dimension checking is enabled (the usual default).
bool operator/=(const dimensionSet &)
Division, modifies the exponents.
bool operator-=(const dimensionSet &) const
Subtraction operation, checks for identical dimensions.
scalar operator[](const dimensionType) const
const FixedList< scalar, 7 > & values() const noexcept
Const access to the exponents as a list.
static const scalar smallExponent
Tolerance for 'small' exponents, for near-zero rounding.
bool dimensionless() const
Return true if it is dimensionless.
void clear()
Clear exponents - resets to be dimensionless.
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
dimensionType
Enumeration for the dimension exponents.
@ LUMINOUS_INTENSITY
Candela cd.
@ MASS
kilogram kg
@ TEMPERATURE
Kelvin K.
bool operator*=(const dimensionSet &)
Multiplication, modifies the exponents.
bool operator=(const dimensionSet &) const
Assignment operation, checks for identical dimensions. Use reset() to change the dimensions.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for OpenFOAM.
tmp< faMatrix< Type > > operator-(const faMatrix< Type > &)
Unary negation.
bitSet operator~(const bitSet &bitset)
Bitset complement, returns a copy of the bitset with all its bits flipped.
Definition bitSetI.H:674
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionSet invTransform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar negPart(const dimensionedScalar &ds)
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static bool checkDims(const char *what, const dimensionSet &a, const dimensionSet &b)
dimensionedScalar hypot(const dimensionedScalar &x, const dimensionedScalar &y)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensionedScalar pow4(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition errorManip.H:139
dimensionSet trans(const dimensionSet &ds)
Check the argument is dimensionless (for transcendental functions).
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar stabilise(const dimensionedScalar &x, const dimensionedScalar &y)
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const direction noexcept
Definition scalarImpl.H:265
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionSet pow2(const dimensionSet &ds)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & b