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-2024 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  template<>
38  {
39  "mass",
40  "length",
41  "time",
42  "temperature",
43  "moles",
44  "current",
45  "luminousIntensity"
46  };
47 }
48 
49 
52 
53 
54 namespace Foam
55 {
56  const scalar dimensionSet::smallExponent = rootSmall;
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const scalar mass,
65  const scalar length,
66  const scalar time,
67  const scalar temperature,
68  const scalar moles,
69  const scalar current,
70  const scalar luminousIntensity
71 )
72 {
73  exponents_[MASS] = mass;
74  exponents_[LENGTH] = length;
75  exponents_[TIME] = time;
76  exponents_[TEMPERATURE] = temperature;
77  exponents_[MOLES] = moles;
78  exponents_[CURRENT] = current;
79  exponents_[LUMINOUS_INTENSITY] = luminousIntensity;
80 }
81 
82 
84 (
85  const scalar mass,
86  const scalar length,
87  const scalar time,
88  const scalar temperature,
89  const scalar moles
90 )
91 {
92  exponents_[MASS] = mass;
93  exponents_[LENGTH] = length;
94  exponents_[TIME] = time;
95  exponents_[TEMPERATURE] = temperature;
96  exponents_[MOLES] = moles;
97  exponents_[CURRENT] = 0;
98  exponents_[LUMINOUS_INTENSITY] = 0;
99 }
100 
101 
103 {
104  reset(ds);
105 }
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
112  for (int Dimension=0; Dimension<nDimensions; ++Dimension)
113  {
114  // ie, mag(exponents_[Dimension]) > smallExponent
115  if
116  (
117  exponents_[Dimension] > smallExponent
118  || exponents_[Dimension] < -smallExponent
119  )
120  {
121  return false;
122  }
123  }
124 
125  return true;
126 }
127 
128 
130 {
131  for (int Dimension=0; Dimension<nDimensions; ++Dimension)
132  {
133  exponents_[Dimension] = ds.exponents_[Dimension];
134  }
135 }
136 
137 
138 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
139 
141 {
142  return exponents_[type];
143 }
144 
145 
147 {
148  return exponents_[type];
149 }
150 
151 
152 Foam::scalar Foam::dimensionSet::operator[](const label type) const
153 {
154  return exponents_[type];
155 }
156 
157 
159 {
160  return exponents_[type];
161 }
162 
163 
165 {
166  for (int Dimension=0; Dimension < nDimensions; ++Dimension)
167  {
168  if
169  (
170  mag(exponents_[Dimension] - ds.exponents_[Dimension])
171  > smallExponent
172  )
173  {
174  return false;
175  }
176  }
177 
178  return true;
179 }
180 
181 
183 {
184  return !(operator==(ds));
185 }
186 
187 
189 {
190  if (dimensionSet::debug && *this != ds)
191  {
193  << "Different dimensions for =" << endl
194  << " dimensions : " << info() << " = " << ds.info() << endl
195  << abort(FatalError);
196  }
197 
198  return true;
199 }
200 
201 
203 {
204  if (dimensionSet::debug && *this != ds)
205  {
207  << "Different dimensions for +=" << endl
208  << " dimensions : " << info() << " = " << ds.info() << endl
209  << abort(FatalError);
210  }
211 
212  return true;
213 }
214 
215 
217 {
218  if (dimensionSet::debug && *this != ds)
219  {
221  << "Different dimensions for -=" << endl
222  << " dimensions : " << info() << " = " << ds.info() << endl
223  << abort(FatalError);
224  }
225 
226  return true;
227 }
228 
229 
231 {
232  reset((*this)*ds);
233 
234  return true;
235 }
236 
237 
239 {
240  reset((*this)/ds);
241 
242  return true;
243 }
244 
245 
246 // * * * * * * * * * * * * * * * Friend functions * * * * * * * * * * * * * * //
247 
249 {
250  if (dimensionSet::debug && ds1 != ds2)
251  {
253  << "Arguments of max have different dimensions" << endl
254  << " dimensions : " << ds1.info() << " and " << ds2.info()
255  << endl << abort(FatalError);
256  }
257 
258  return ds1;
259 }
260 
261 
263 {
264  if (dimensionSet::debug && ds1 != ds2)
265  {
267  << "Arguments of min have different dimensions" << endl
268  << " dimensions : " << ds1.info() << " and " << ds2.info()
269  << endl << abort(FatalError);
270  }
271 
272  return ds1;
273 }
274 
275 
277 (
278  const dimensionSet& ds1,
279  const dimensionSet& ds2
280 )
281 {
282  return ds1*ds2;
283 }
284 
285 
287 (
288  const dimensionSet& ds1,
289  const dimensionSet& ds2
290 )
291 {
292  return ds1/ds2;
293 }
294 
295 
297 {
298  return ds;
299 }
300 
301 
302 Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
303 {
304  dimensionSet dimPow
305  (
306  ds[dimensionSet::MASS]*p,
308  ds[dimensionSet::TIME]*p,
313  );
314 
315  return dimPow;
316 }
317 
318 
320 (
321  const dimensionSet& ds,
322  const dimensionedScalar& dS
323 )
324 {
325  if (dimensionSet::debug && !dS.dimensions().dimensionless())
326  {
328  << "Exponent of pow is not dimensionless"
329  << abort(FatalError);
330  }
331 
332  dimensionSet dimPow
333  (
334  ds[dimensionSet::MASS]*dS.value(),
335  ds[dimensionSet::LENGTH]*dS.value(),
336  ds[dimensionSet::TIME]*dS.value(),
338  ds[dimensionSet::MOLES]*dS.value(),
339  ds[dimensionSet::CURRENT]*dS.value(),
341  );
342 
343  return dimPow;
344 }
345 
346 
348 (
349  const dimensionedScalar& dS,
350  const dimensionSet& ds
351 )
352 {
353  if
354  (
355  dimensionSet::debug
356  && !dS.dimensions().dimensionless()
357  && !ds.dimensionless()
358  )
359  {
361  << "Argument or exponent of pow not dimensionless" << endl
362  << abort(FatalError);
363  }
364 
365  return ds;
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 sqrt(sqrt(ds));
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 dimless/ds;
474 }
475 
476 
478 {
479  if (dimensionSet::debug && !ds.dimensionless())
480  {
482  << "Argument of trancendental function not dimensionless"
483  << abort(FatalError);
484  }
485 
486  return ds;
487 }
488 
489 
491 {
492  if (dimensionSet::debug && ds1 != ds2)
493  {
495  << "Arguments of atan2 have different dimensions" << endl
496  << " dimensions : " << ds1.info() << " and " << ds2.info()
497  << endl << abort(FatalError);
498  }
499 
500  return dimless;
501 }
502 
503 
505 {
506  return ds;
507 }
508 
509 
511 {
512  return dimless;
513 }
514 
515 
517 {
518  return dimless;
519 }
520 
521 
522 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
523 
525 {
526  return ds;
527 }
528 
529 
530 Foam::dimensionSet Foam::operator+
531 (
532  const dimensionSet& ds1,
533  const dimensionSet& ds2
534 )
535 {
536  dimensionSet dimSum(ds1);
537 
538  if (dimensionSet::debug && ds1 != ds2)
539  {
541  << "LHS and RHS of + have different dimensions" << endl
542  << " dimensions : " << ds1.info() << " + " << ds2.info()
543  << endl << abort(FatalError);
544  }
545 
546  return dimSum;
547 }
548 
549 
550 Foam::dimensionSet Foam::operator-
551 (
552  const dimensionSet& ds1,
553  const dimensionSet& ds2
554 )
555 {
556  dimensionSet dimDifference(ds1);
557 
558  if (dimensionSet::debug && ds1 != ds2)
559  {
561  << "LHS and RHS of - have different dimensions" << endl
562  << " dimensions : " << ds1.info() << " - " << ds2.info()
563  << endl << abort(FatalError);
564  }
565 
566  return dimDifference;
567 }
568 
569 
570 Foam::dimensionSet Foam::operator*
571 (
572  const dimensionSet& ds1,
573  const dimensionSet& ds2
574 )
575 {
576  dimensionSet dimProduct(ds1);
577 
578  for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
579  {
580  dimProduct.exponents_[Dimension] += ds2.exponents_[Dimension];
581  }
582 
583  return dimProduct;
584 }
585 
586 
587 Foam::dimensionSet Foam::operator/
588 (
589  const dimensionSet& ds1,
590  const dimensionSet& ds2
591 )
592 {
593  dimensionSet dimQuotient(ds1);
594 
595  for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
596  {
597  dimQuotient.exponents_[Dimension] -= ds2.exponents_[Dimension];
598  }
599 
600  return dimQuotient;
601 }
602 
603 
604 Foam::dimensionSet Foam::operator&
605 (
606  const dimensionSet& ds1,
607  const dimensionSet& ds2
608 )
609 {
610  return ds1*ds2;
611 }
612 
613 
614 Foam::dimensionSet Foam::operator^
615 (
616  const dimensionSet& ds1,
617  const dimensionSet& ds2
618 )
619 {
620  return ds1*ds2;
621 }
622 
623 
624 Foam::dimensionSet Foam::operator&&
625 (
626  const dimensionSet& ds1,
627  const dimensionSet& ds2
628 )
629 {
630  return ds1*ds2;
631 }
632 
633 
634 // ************************************************************************* //
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
Dimension set for the base types.
Definition: dimensionSet.H:125
bool operator+=(const dimensionSet &) const
Definition: dimensionSet.C:202
InfoProxy< dimensionSet > info() const
Return info proxy.
Definition: dimensionSet.H:239
bool operator==(const dimensionSet &) const
Definition: dimensionSet.C:164
bool operator!=(const dimensionSet &) const
Definition: dimensionSet.C:182
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:63
static const scalar smallExponent
A small exponent with which to perform inexact comparisons.
Definition: dimensionSet.H:156
bool operator/=(const dimensionSet &)
Definition: dimensionSet.C:238
bool operator-=(const dimensionSet &) const
Definition: dimensionSet.C:216
scalar operator[](const dimensionType) const
Definition: dimensionSet.C:140
void reset(const dimensionSet &)
Definition: dimensionSet.C:129
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:110
dimensionType
Define an enumeration for the names of the dimension exponents.
Definition: dimensionSet.H:139
bool operator*=(const dimensionSet &)
Definition: dimensionSet.C:230
bool operator=(const dimensionSet &) const
Definition: dimensionSet.C:188
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 pow6(const dimensionedScalar &ds)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:510
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:296
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionSet perpendicular(const dimensionSet &)
Definition: dimensionSet.C:516
tmp< fvMatrix< Type > > operator-(const fvMatrix< Type > &)
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:477
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & p