dimensionSet.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 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 "dimensionSet.H"
27 #include "dimensionedScalar.H"
28 #include "OStringStream.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
35 }
36 
37 
40 {
41  "mass",
42  "length",
43  "time",
44  "temperature",
45  "moles",
46  "current",
47  "luminousIntensity"
48 };
49 
50 
51 namespace Foam
52 {
53  const scalar dimensionSet::smallExponent = rootSmall;
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const scalar mass,
62  const scalar length,
63  const scalar time,
64  const scalar temperature,
65  const scalar moles,
66  const scalar current,
67  const scalar luminousIntensity
68 )
69 {
70  exponents_[MASS] = mass;
71  exponents_[LENGTH] = length;
72  exponents_[TIME] = time;
73  exponents_[TEMPERATURE] = temperature;
74  exponents_[MOLES] = moles;
75  exponents_[CURRENT] = current;
76  exponents_[LUMINOUS_INTENSITY] = luminousIntensity;
77 }
78 
79 
81 (
82  const scalar mass,
83  const scalar length,
84  const scalar time,
85  const scalar temperature,
86  const scalar moles
87 )
88 {
89  exponents_[MASS] = mass;
90  exponents_[LENGTH] = length;
91  exponents_[TIME] = time;
92  exponents_[TEMPERATURE] = temperature;
93  exponents_[MOLES] = moles;
94  exponents_[CURRENT] = 0;
95  exponents_[LUMINOUS_INTENSITY] = 0;
96 }
97 
98 
100 {
101  reset(ds);
102 }
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  for (int Dimension=0; Dimension<nDimensions; ++Dimension)
110  {
111  // ie, mag(exponents_[Dimension]) > smallExponent
112  if
113  (
114  exponents_[Dimension] > smallExponent
115  || exponents_[Dimension] < -smallExponent
116  )
117  {
118  return false;
119  }
120  }
121 
122  return true;
123 }
124 
125 
127 {
128  for (int Dimension=0; Dimension<nDimensions; ++Dimension)
129  {
130  exponents_[Dimension] = ds.exponents_[Dimension];
131  }
132 }
133 
134 
135 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
136 
138 {
139  return exponents_[type];
140 }
141 
142 
144 {
145  return exponents_[type];
146 }
147 
148 
149 Foam::scalar Foam::dimensionSet::operator[](const label type) const
150 {
151  return exponents_[type];
152 }
153 
154 
156 {
157  return exponents_[type];
158 }
159 
160 
162 {
163  for (int Dimension=0; Dimension < nDimensions; ++Dimension)
164  {
165  if
166  (
167  mag(exponents_[Dimension] - ds.exponents_[Dimension])
168  > smallExponent
169  )
170  {
171  return false;
172  }
173  }
174 
175  return true;
176 }
177 
178 
180 {
181  return !(operator==(ds));
182 }
183 
184 
186 {
187  if (dimensionSet::debug && *this != ds)
188  {
190  << "Different dimensions for =" << endl
191  << " dimensions : " << info() << " = " << ds.info() << endl
192  << abort(FatalError);
193  }
194 
195  return true;
196 }
197 
198 
200 {
201  if (dimensionSet::debug && *this != ds)
202  {
204  << "Different dimensions for +=" << endl
205  << " dimensions : " << info() << " = " << ds.info() << endl
206  << abort(FatalError);
207  }
208 
209  return true;
210 }
211 
212 
214 {
215  if (dimensionSet::debug && *this != ds)
216  {
218  << "Different dimensions for -=" << endl
219  << " dimensions : " << info() << " = " << ds.info() << endl
220  << abort(FatalError);
221  }
222 
223  return true;
224 }
225 
226 
228 {
229  reset((*this)*ds);
230 
231  return true;
232 }
233 
234 
236 {
237  reset((*this)/ds);
238 
239  return true;
240 }
241 
242 
243 // * * * * * * * * * * * * * * * Friend functions * * * * * * * * * * * * * * //
244 
246 {
247  if (dimensionSet::debug && ds1 != ds2)
248  {
250  << "Arguments of max have different dimensions" << endl
251  << " dimensions : " << ds1.info() << " and " << ds2.info()
252  << endl << abort(FatalError);
253  }
254 
255  return ds1;
256 }
257 
258 
260 {
261  if (dimensionSet::debug && ds1 != ds2)
262  {
264  << "Arguments of min have different dimensions" << endl
265  << " dimensions : " << ds1.info() << " and " << ds2.info()
266  << endl << abort(FatalError);
267  }
268 
269  return ds1;
270 }
271 
272 
274 (
275  const dimensionSet& ds1,
276  const dimensionSet& ds2
277 )
278 {
279  return ds1*ds2;
280 }
281 
282 
284 (
285  const dimensionSet& ds1,
286  const dimensionSet& ds2
287 )
288 {
289  return ds1/ds2;
290 }
291 
292 
294 {
295  return ds;
296 }
297 
298 
299 Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
300 {
301  dimensionSet dimPow
302  (
303  ds[dimensionSet::MASS]*p,
305  ds[dimensionSet::TIME]*p,
310  );
311 
312  return dimPow;
313 }
314 
315 
317 (
318  const dimensionSet& ds,
319  const dimensionedScalar& dS
320 )
321 {
322  if (dimensionSet::debug && !dS.dimensions().dimensionless())
323  {
325  << "Exponent of pow is not dimensionless"
326  << abort(FatalError);
327  }
328 
329  dimensionSet dimPow
330  (
331  ds[dimensionSet::MASS]*dS.value(),
332  ds[dimensionSet::LENGTH]*dS.value(),
333  ds[dimensionSet::TIME]*dS.value(),
335  ds[dimensionSet::MOLES]*dS.value(),
336  ds[dimensionSet::CURRENT]*dS.value(),
338  );
339 
340  return dimPow;
341 }
342 
343 
345 (
346  const dimensionedScalar& dS,
347  const dimensionSet& ds
348 )
349 {
350  if
351  (
352  dimensionSet::debug
353  && !dS.dimensions().dimensionless()
354  && !ds.dimensionless()
355  )
356  {
358  << "Argument or exponent of pow not dimensionless" << endl
359  << abort(FatalError);
360  }
361 
362  return ds;
363 }
364 
365 
367 {
368  return pow(ds, 2);
369 }
370 
371 
373 {
374  return pow(ds, 3);
375 }
376 
377 
379 {
380  return pow(ds, 4);
381 }
382 
383 
385 {
386  return pow(ds, 5);
387 }
388 
389 
391 {
392  return pow(ds, 6);
393 }
394 
395 
397 {
398  return sqrt(sqrt(ds));
399 }
400 
401 
403 {
404  return pow(ds, 0.5);
405 }
406 
407 
409 {
410  return pow(ds, 1.0/3.0);
411 }
412 
413 
415 {
416  return pow(ds, 2);
417 }
418 
419 
421 {
422  return ds;
423 }
424 
425 
427 {
428  return dimless;
429 }
430 
431 
433 {
434  return dimless;
435 }
436 
437 
439 {
440  return dimless;
441 }
442 
443 
445 {
446  return dimless;
447 }
448 
449 
451 {
452  return dimless;
453 }
454 
455 
457 {
458  return ds;
459 }
460 
461 
463 {
464  return ds;
465 }
466 
467 
469 {
470  return dimless/ds;
471 }
472 
473 
475 {
476  if (dimensionSet::debug && !ds.dimensionless())
477  {
479  << "Argument of trancendental function not dimensionless"
480  << abort(FatalError);
481  }
482 
483  return ds;
484 }
485 
486 
488 {
489  if (dimensionSet::debug && ds1 != ds2)
490  {
492  << "Arguments of atan2 have different dimensions" << endl
493  << " dimensions : " << ds1.info() << " and " << ds2.info()
494  << endl << abort(FatalError);
495  }
496 
497  return dimless;
498 }
499 
500 
502 {
503  return ds;
504 }
505 
506 
508 {
509  return dimless;
510 }
511 
512 
514 {
515  return dimless;
516 }
517 
518 
519 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
520 
522 {
523  return ds;
524 }
525 
526 
527 Foam::dimensionSet Foam::operator+
528 (
529  const dimensionSet& ds1,
530  const dimensionSet& ds2
531 )
532 {
533  dimensionSet dimSum(ds1);
534 
535  if (dimensionSet::debug && ds1 != ds2)
536  {
538  << "LHS and RHS of + have different dimensions" << endl
539  << " dimensions : " << ds1.info() << " + " << ds2.info()
540  << endl << abort(FatalError);
541  }
542 
543  return dimSum;
544 }
545 
546 
547 Foam::dimensionSet Foam::operator-
548 (
549  const dimensionSet& ds1,
550  const dimensionSet& ds2
551 )
552 {
553  dimensionSet dimDifference(ds1);
554 
555  if (dimensionSet::debug && ds1 != ds2)
556  {
558  << "LHS and RHS of - have different dimensions" << endl
559  << " dimensions : " << ds1.info() << " - " << ds2.info()
560  << endl << abort(FatalError);
561  }
562 
563  return dimDifference;
564 }
565 
566 
567 Foam::dimensionSet Foam::operator*
568 (
569  const dimensionSet& ds1,
570  const dimensionSet& ds2
571 )
572 {
573  dimensionSet dimProduct(ds1);
574 
575  for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
576  {
577  dimProduct.exponents_[Dimension] += ds2.exponents_[Dimension];
578  }
579 
580  return dimProduct;
581 }
582 
583 
584 Foam::dimensionSet Foam::operator/
585 (
586  const dimensionSet& ds1,
587  const dimensionSet& ds2
588 )
589 {
590  dimensionSet dimQuotient(ds1);
591 
592  for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
593  {
594  dimQuotient.exponents_[Dimension] -= ds2.exponents_[Dimension];
595  }
596 
597  return dimQuotient;
598 }
599 
600 
601 Foam::dimensionSet Foam::operator&
602 (
603  const dimensionSet& ds1,
604  const dimensionSet& ds2
605 )
606 {
607  return ds1*ds2;
608 }
609 
610 
611 Foam::dimensionSet Foam::operator^
612 (
613  const dimensionSet& ds1,
614  const dimensionSet& ds2
615 )
616 {
617  return ds1*ds2;
618 }
619 
620 
621 Foam::dimensionSet Foam::operator&&
622 (
623  const dimensionSet& ds1,
624  const dimensionSet& ds2
625 )
626 {
627  return ds1*ds2;
628 }
629 
630 
631 // ************************************************************************* //
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Dimension set for the base types.
Definition: dimensionSet.H:125
bool operator+=(const dimensionSet &) const
Definition: dimensionSet.C:199
InfoProxy< dimensionSet > info() const
Return info proxy.
Definition: dimensionSet.H:245
bool operator==(const dimensionSet &) const
Definition: dimensionSet.C:161
bool operator!=(const dimensionSet &) const
Definition: dimensionSet.C:179
dimensionSet(const scalar mass, const scalar length, const scalar time, const scalar temperature, const scalar moles, const scalar current, const scalar luminousIntensity)
Construct given individual dimension exponents for all.
Definition: dimensionSet.C:60
static const scalar smallExponent
A small exponent with which to perform inexact comparisons.
Definition: dimensionSet.H:156
bool operator/=(const dimensionSet &)
Definition: dimensionSet.C:235
bool operator-=(const dimensionSet &) const
Definition: dimensionSet.C:213
scalar operator[](const dimensionType) const
Definition: dimensionSet.C:137
void reset(const dimensionSet &)
Definition: dimensionSet.C:126
static const NamedEnum< dimensionType, 7 > dimensionTypeNames_
Names of the dimensions.
Definition: dimensionSet.H:150
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:107
dimensionType
Define an enumeration for the names of the dimension exponents.
Definition: dimensionSet.H:139
bool operator*=(const dimensionSet &)
Definition: dimensionSet.C:227
bool operator=(const dimensionSet &) const
Definition: dimensionSet.C:185
const dimensionSet & dimensions() const
Return const reference to dimensions.
const Type & value() const
Return const reference to value.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
void cmptMultiply(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< Type > &f2)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void cmptDivide(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< Type > &f2)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void pow6(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar negPart(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:501
void pow5(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensionedScalar neg(const dimensionedScalar &ds)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
error FatalError
dimensionedScalar neg0(const dimensionedScalar &ds)
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
dimensionSet perpendicular(const dimensionSet &)
Definition: dimensionSet.C:513
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
tmp< fvMatrix< Type > > operator-(const fvMatrix< Type > &)
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:474
void cmptMag(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
dimensionedScalar posPart(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
volScalarField & p