dimensionSet.H
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-2021 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 Class
25  Foam::dimensionSet
26 
27 Description
28  Dimension set for the base types.
29 
30  This type may be used to implement rigorous dimension checking
31  for algebraic manipulation.
32 
33 SourceFiles
34  dimensionSet.C
35  dimensionSetIO.C
36  dimensionSets.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef dimensionSet_H
41 #define dimensionSet_H
42 
43 #include "bool.H"
44 #include "dimensionedScalarFwd.H"
45 #include "className.H"
46 #include "scalarField.H"
47 #include "PtrList.H"
48 #include "HashTable.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 // Forward declaration of friend functions and operators
56 
57 class dimensionSet;
58 class dimensionSets;
59 
60 // Friend Functions
61 
62 dimensionSet max(const dimensionSet&, const dimensionSet&);
63 dimensionSet min(const dimensionSet&, const dimensionSet&);
64 dimensionSet cmptMultiply(const dimensionSet&, const dimensionSet&);
65 dimensionSet cmptDivide(const dimensionSet&, const dimensionSet&);
66 dimensionSet cmptMag(const dimensionSet&);
67 
68 dimensionSet pow(const dimensionSet&, const scalar);
69 dimensionSet pow(const dimensionSet&, const dimensionedScalar&);
70 dimensionSet pow(const dimensionedScalar&, const dimensionSet&);
71 
72 dimensionSet sqr(const dimensionSet&);
73 dimensionSet pow3(const dimensionSet&);
74 dimensionSet pow4(const dimensionSet&);
75 dimensionSet pow5(const dimensionSet&);
76 dimensionSet pow6(const dimensionSet&);
77 dimensionSet pow025(const dimensionSet&);
78 
79 dimensionSet sqrt(const dimensionSet&);
80 dimensionSet cbrt(const dimensionSet&);
81 dimensionSet magSqr(const dimensionSet&);
82 dimensionSet mag(const dimensionSet&);
83 dimensionSet sign(const dimensionSet&);
84 dimensionSet pos(const dimensionSet&);
85 dimensionSet pos0(const dimensionSet&);
86 dimensionSet neg(const dimensionSet&);
87 dimensionSet neg0(const dimensionSet&);
88 dimensionSet posPart(const dimensionSet&);
89 dimensionSet negPart(const dimensionSet&);
90 dimensionSet inv(const dimensionSet&);
91 
92 // Function to check the argument is dimensionless
93 // for transcendental functions
94 dimensionSet trans(const dimensionSet&);
95 
96 dimensionSet atan2(const dimensionSet&, const dimensionSet&);
97 
98 // Return the argument; transformations do not change the dimensions
99 dimensionSet transform(const dimensionSet&);
100 
101 // Friend operators
102 
103 dimensionSet operator-(const dimensionSet&);
104 dimensionSet operator+(const dimensionSet&, const dimensionSet&);
105 dimensionSet operator-(const dimensionSet&, const dimensionSet&);
106 dimensionSet operator*(const dimensionSet&, const dimensionSet&);
107 dimensionSet operator/(const dimensionSet&, const dimensionSet&);
108 dimensionSet operator&(const dimensionSet&, const dimensionSet&);
109 dimensionSet operator^(const dimensionSet&, const dimensionSet&);
110 dimensionSet operator&&(const dimensionSet&, const dimensionSet&);
111 
112 // IOstream Operators
113 
114 Istream& operator>>(Istream&, dimensionSet&);
115 Ostream& operator<<(Ostream&, const dimensionSet&);
116 
117 
118 /*---------------------------------------------------------------------------*\
119  Class dimensionSet Declaration
120 \*---------------------------------------------------------------------------*/
122 class dimensionSet
123 {
124 
125 public:
126 
127  // Member constants
128 
129  enum
130  {
131  nDimensions = 7 // Number of dimensions in SI is 7
132  };
133 
134  //- Define an enumeration for the names of the dimension exponents
135  enum dimensionType
136  {
137  MASS, // kilogram kg
138  LENGTH, // metre m
139  TIME, // second s
140  TEMPERATURE, // Kelvin K
141  MOLES, // mole mol
142  CURRENT, // Ampere A
143  LUMINOUS_INTENSITY // Candela Cd
144  };
145 
146 
147  // Static Data Members
149  static const scalar smallExponent;
150 
151 
152 private:
153 
154  // Private classes
155 
156  class tokeniser
157  {
158  // Private Data
159 
160  Istream& is_;
161 
162  List<token> tokens_;
163 
164  label start_;
165 
166  label size_;
167 
168 
169  // Private Member Functions
170 
171  void push(const token&);
172 
173  token pop();
174 
175  void unpop(const token&);
176 
177  public:
178 
179  // Constructors
180 
181  tokeniser(Istream&);
182 
183 
184  // Member Functions
185 
186  Istream& stream()
187  {
188  return is_;
189  }
190 
191  bool hasToken() const;
192 
193  token nextToken();
194 
195  void putBack(const token&);
196 
197  void splitWord(const word&);
198 
199  static bool valid(char c);
200 
201  static label priority(const token& t);
202  };
203 
204 
205  //- Reset exponents to nearest integer if close to it. Used to
206  // handle reading with insufficient precision.
207  void round(const scalar tol);
208 
209  dimensionedScalar parse
210  (
211  const label lastPrior,
212  tokeniser& tis,
214  ) const;
215 
216 
217  // private data
218 
219  // dimensionSet stored as an array of dimension exponents
220  scalar exponents_[nDimensions];
221 
222 
223 public:
224 
225  // Declare name of the class and its debug switch
226  ClassName("dimensionSet");
227 
228 
229  // Constructors
230 
231  //- Construct given individual dimension exponents for all
232  // seven dimensions
234  (
235  const scalar mass,
236  const scalar length,
237  const scalar time,
238  const scalar temperature,
239  const scalar moles,
240  const scalar current,
241  const scalar luminousIntensity
242  );
243 
244  //- Construct given individual dimension exponents for first
245  // five dimensions
247  (
248  const scalar mass,
249  const scalar length,
250  const scalar time,
251  const scalar temperature,
252  const scalar moles
253  );
254 
255  //- Copy constructor
256  dimensionSet(const dimensionSet& ds);
257 
258  //- Construct and return a clone
260  {
261  return autoPtr<dimensionSet>(new dimensionSet(*this));
262  }
263 
264  //- Construct from Istream
266 
267 
268  // Member Functions
269 
270  //- Return true if it is dimensionless
271  bool dimensionless() const;
272 
273  void reset(const dimensionSet&);
274 
275 
276  // I/O
277 
278  //- Read using provided units. Used only in initial parsing
279  Istream& read
280  (
281  Istream& is,
282  scalar& multiplier,
283  const dictionary&
284  );
285 
286  //- Read using provided units
287  Istream& read
288  (
289  Istream& is,
290  scalar& multiplier,
292  );
293 
294  //- Read using system units
295  Istream& read
296  (
297  Istream& is,
298  scalar& multiplier
299  );
300 
301  //- Write using provided units
302  Ostream& write
303  (
304  Ostream& os,
305  scalar& multiplier,
306  const dimensionSets&
307  ) const;
308 
309  //- Write using system units
310  Ostream& write
311  (
312  Ostream& os,
313  scalar& multiplier
314  ) const;
315 
316 
317  // Member Operators
318 
319  scalar operator[](const dimensionType) const;
320  scalar& operator[](const dimensionType);
321 
322  scalar operator[](const label) const;
323  scalar& operator[](const label);
324 
325  bool operator==(const dimensionSet&) const;
326  bool operator!=(const dimensionSet&) const;
327 
328  bool operator=(const dimensionSet&) const;
329 
330  bool operator+=(const dimensionSet&) const;
331  bool operator-=(const dimensionSet&) const;
332  bool operator*=(const dimensionSet&);
333  bool operator/=(const dimensionSet&);
334 
335 
336  // Friend Functions
337 
338  friend dimensionSet max(const dimensionSet&, const dimensionSet&);
339  friend dimensionSet min(const dimensionSet&, const dimensionSet&);
340  friend dimensionSet cmptMultiply
341  (
342  const dimensionSet&,
343  const dimensionSet&
344  );
345  friend dimensionSet cmptDivide
346  (
347  const dimensionSet&,
348  const dimensionSet&
349  );
350  friend dimensionSet cmptMag(const dimensionSet&);
351 
352  friend dimensionSet pow(const dimensionSet&, const scalar);
353  friend dimensionSet pow(const dimensionSet&, const dimensionedScalar&);
354  friend dimensionSet pow(const dimensionedScalar&, const dimensionSet&);
355 
356  friend dimensionSet sqr(const dimensionSet&);
357  friend dimensionSet pow3(const dimensionSet&);
358  friend dimensionSet pow4(const dimensionSet&);
359  friend dimensionSet pow5(const dimensionSet&);
360  friend dimensionSet pow6(const dimensionSet&);
361  friend dimensionSet pow025(const dimensionSet&);
362 
363  friend dimensionSet sqrt(const dimensionSet&);
364  friend dimensionSet magSqr(const dimensionSet&);
365  friend dimensionSet mag(const dimensionSet&);
366  friend dimensionSet sign(const dimensionSet&);
367  friend dimensionSet pos0(const dimensionSet&);
368  friend dimensionSet neg(const dimensionSet&);
369  friend dimensionSet inv(const dimensionSet&);
370 
371  //- Function to check the argument is dimensionless
372  // for transcendental functions
373  friend dimensionSet trans(const dimensionSet&);
374 
375  friend dimensionSet atan2(const dimensionSet&, const dimensionSet&);
376 
377  //- Return the argument; transformations do not change the dimensions
378  friend dimensionSet transform(const dimensionSet&);
379 
380 
381  // Friend operators
382 
383  friend dimensionSet operator-(const dimensionSet&);
384 
385  friend dimensionSet operator+
386  (
387  const dimensionSet&,
388  const dimensionSet&
389  );
390 
391  friend dimensionSet operator-
392  (
393  const dimensionSet&,
394  const dimensionSet&
395  );
396 
397  friend dimensionSet operator*
398  (
399  const dimensionSet&,
400  const dimensionSet&
401  );
402 
403  friend dimensionSet operator/
404  (
405  const dimensionSet&,
406  const dimensionSet&
407  );
408 
409  friend dimensionSet operator&
410  (
411  const dimensionSet&,
412  const dimensionSet&
413  );
414 
415  friend dimensionSet operator^
416  (
417  const dimensionSet&,
418  const dimensionSet&
419  );
420 
421  friend dimensionSet operator&&
422  (
423  const dimensionSet&,
424  const dimensionSet&
425  );
426 
427 
428  // IOstream Operators
429 
430  friend Istream& operator>>(Istream&, dimensionSet&);
431  friend Ostream& operator<<(Ostream&, const dimensionSet&);
432 };
433 
434 
435 void writeEntry(Ostream& os, const dimensionSet& value);
436 
437 
438 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439 
440 } // End namespace Foam
441 
442 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
443 
444 #include "dimensionSets.H"
445 
446 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
447 
448 #endif
449 
450 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionSet trans(const dimensionSet &)
Definition: dimensionSet.C:456
tmp< fvMatrix< Type > > operator*(const volScalarField::Internal &, const fvMatrix< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool operator/=(const dimensionSet &)
Definition: dimensionSet.C:217
autoPtr< dimensionSet > clone() const
Construct and return a clone.
Definition: dimensionSet.H:258
friend dimensionSet pow(const dimensionSet &, const scalar)
HashSet< Key, Hash > operator^(const HashSet< Key, Hash > &hash1, const HashSet< Key, Hash > &hash2)
Create a HashSet that only contains unique entries (xor)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
friend dimensionSet mag(const dimensionSet &)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
friend dimensionSet neg(const dimensionSet &)
friend dimensionSet pow025(const dimensionSet &)
A token holds items read from Istream.
Definition: token.H:72
bool operator-=(const dimensionSet &) const
Definition: dimensionSet.C:195
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
friend dimensionSet pow4(const dimensionSet &)
bool operator=(const dimensionSet &) const
Definition: dimensionSet.C:167
ClassName("dimensionSet")
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:42
friend Istream & operator>>(Istream &, dimensionSet &)
tmp< fvMatrix< Type > > operator/(const fvMatrix< Type > &, const volScalarField::Internal &)
dimensionedScalar posPart(const dimensionedScalar &ds)
Ostream & write(Ostream &os, scalar &multiplier, const dimensionSets &) const
Write using provided units.
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedScalar neg(const dimensionedScalar &ds)
dimensionedScalar pow5(const dimensionedScalar &ds)
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
bool operator*=(const dimensionSet &)
Definition: dimensionSet.C:209
tmp< GeometricField< Type, fvPatchField, volMesh > > operator &(const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &)
friend dimensionSet sign(const dimensionSet &)
dimensionedScalar pos(const dimensionedScalar &ds)
friend dimensionSet pow6(const dimensionSet &)
friend dimensionSet magSqr(const dimensionSet &)
Dimension set for the base types.
Definition: dimensionSet.H:121
friend Ostream & operator<<(Ostream &, const dimensionSet &)
dimensionSet operator &&(const dimensionSet &, const dimensionSet &)
A class for handling words, derived from string.
Definition: word.H:59
Istream & operator>>(Istream &, directionInfo &)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensionedScalar neg0(const dimensionedScalar &ds)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
tmp< fvMatrix< Type > > operator-(const fvMatrix< Type > &)
bool operator+=(const dimensionSet &) const
Definition: dimensionSet.C:181
tmp< fvMatrix< Type > > operator+(const fvMatrix< Type > &, const fvMatrix< Type > &)
An STL-conforming hash table.
Definition: HashTable.H:61
friend dimensionSet max(const dimensionSet &, const dimensionSet &)
friend dimensionSet cmptDivide(const dimensionSet &, const dimensionSet &)
bool operator!=(const dimensionSet &) const
Definition: dimensionSet.C:161
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
friend dimensionSet min(const dimensionSet &, const dimensionSet &)
friend dimensionSet operator-(const dimensionSet &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
friend dimensionSet pow3(const dimensionSet &)
friend dimensionSet cmptMultiply(const dimensionSet &, const dimensionSet &)
friend dimensionSet trans(const dimensionSet &)
Function to check the argument is dimensionless.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void reset(const dimensionSet &)
Definition: dimensionSet.C:108
dimensionedScalar pow3(const dimensionedScalar &ds)
scalar operator[](const dimensionType) const
Definition: dimensionSet.C:119
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Ostream & operator<<(Ostream &, const ensightPart &)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar pow6(const dimensionedScalar &ds)
dimensionType
Define an enumeration for the names of the dimension exponents.
Definition: dimensionSet.H:134
dimensioned< scalar > mag(const dimensioned< Type > &)
friend dimensionSet transform(const dimensionSet &)
Return the argument; transformations do not change the dimensions.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Useful dimension sets.
friend dimensionSet sqr(const dimensionSet &)
friend dimensionSet pow5(const dimensionSet &)
friend dimensionSet pos0(const dimensionSet &)
bool operator==(const dimensionSet &) const
Definition: dimensionSet.C:143
static const scalar smallExponent
Definition: dimensionSet.H:148
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:89
friend dimensionSet inv(const dimensionSet &)
friend dimensionSet atan2(const dimensionSet &, const dimensionSet &)
friend dimensionSet cmptMag(const dimensionSet &)
Istream & read(Istream &is, scalar &multiplier, const dictionary &)
Read using provided units. Used only in initial parsing.
System bool.
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:483
friend dimensionSet sqrt(const dimensionSet &)