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-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 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 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef dimensionSet_H
40 #define dimensionSet_H
41 
42 #include "bool.H"
43 #include "className.H"
44 #include "NamedEnum.H"
45 #include "InfoProxy.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward declaration of friend functions and operators
53 
54 class dimensionSet;
55 class unitConversion;
56 template<class Type> class dimensioned;
57 typedef dimensioned<scalar> dimensionedScalar;
58 
59 // Friend Functions
60 
61 dimensionSet max(const dimensionSet&, const dimensionSet&);
62 dimensionSet min(const dimensionSet&, const dimensionSet&);
63 dimensionSet cmptMultiply(const dimensionSet&, const dimensionSet&);
64 dimensionSet cmptDivide(const dimensionSet&, const dimensionSet&);
65 dimensionSet cmptMag(const dimensionSet&);
66 
67 dimensionSet pow(const dimensionSet&, const scalar);
68 dimensionSet pow(const dimensionSet&, const dimensionedScalar&);
69 dimensionSet pow(const dimensionedScalar&, const dimensionSet&);
70 
71 dimensionSet sqr(const dimensionSet&);
72 dimensionSet pow3(const dimensionSet&);
73 dimensionSet pow4(const dimensionSet&);
74 dimensionSet pow5(const dimensionSet&);
75 dimensionSet pow6(const dimensionSet&);
76 dimensionSet pow025(const dimensionSet&);
77 
78 dimensionSet sqrt(const dimensionSet&);
79 dimensionSet cbrt(const dimensionSet&);
80 dimensionSet magSqr(const dimensionSet&);
81 dimensionSet mag(const dimensionSet&);
82 dimensionSet sign(const dimensionSet&);
83 dimensionSet pos(const dimensionSet&);
84 dimensionSet pos0(const dimensionSet&);
85 dimensionSet neg(const dimensionSet&);
86 dimensionSet neg0(const dimensionSet&);
87 dimensionSet posPart(const dimensionSet&);
88 dimensionSet negPart(const dimensionSet&);
89 dimensionSet inv(const dimensionSet&);
90 
91 // Function to check the argument is dimensionless
92 // for transcendental functions
93 dimensionSet trans(const dimensionSet&);
94 
95 dimensionSet atan2(const dimensionSet&, const dimensionSet&);
96 
97 // Return the argument; transformations do not change the dimensions
98 dimensionSet transform(const dimensionSet&);
99 
100 dimensionSet normalised(const dimensionSet&);
101 dimensionSet perpendicular(const dimensionSet&);
102 
103 // Friend operators
104 
105 dimensionSet operator-(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 dimensionSet operator^(const dimensionSet&, const dimensionSet&);
112 dimensionSet operator&&(const dimensionSet&, const dimensionSet&);
113 
114 // IOstream Operators
115 
116 Istream& operator>>(Istream&, dimensionSet&);
117 Ostream& operator<<(Ostream&, const dimensionSet&);
118 Ostream& operator<<(Ostream&, const InfoProxy<dimensionSet>&);
119 
120 
121 /*---------------------------------------------------------------------------*\
122  Class dimensionSet Declaration
123 \*---------------------------------------------------------------------------*/
124 
125 class dimensionSet
126 {
127 
128 public:
129 
130  // Member constants
131 
132  //- Define an enumeration for the number of dimensions
133  enum
134  {
135  nDimensions = 7 // Number of dimensions in SI is 7
136  };
137 
138  //- Define an enumeration for the names of the dimension exponents
139  enum dimensionType
140  {
141  MASS, // kilogram kg
142  LENGTH, // metre m
143  TIME, // second s
144  TEMPERATURE, // Kelvin K
145  MOLES, // kilomole kmol
146  CURRENT, // Ampere A
147  LUMINOUS_INTENSITY // Candela Cd
148  };
149 
150  //- Names of the dimensions
152 
153 
154  // Static Data Members
155 
156  //- A small exponent with which to perform inexact comparisons
157  static const scalar smallExponent;
158 
159 
160 private:
161 
162  // Private Data
163 
164  //- Array of dimension exponents
165  scalar exponents_[nDimensions];
166 
167 
168  // Private Member Functions
169 
170  //- Reset exponents to nearest integer if close to it. Used to
171  // handle reading with insufficient precision.
172  void round(const scalar tol);
173 
174 
175 public:
176 
177  // Declare name of the class and its debug switch
178  ClassName("dimensionSet");
179 
180 
181  // Constructors
182 
183  //- Construct given individual dimension exponents for all
184  // seven dimensions
186  (
187  const scalar mass,
188  const scalar length,
189  const scalar time,
190  const scalar temperature,
191  const scalar moles,
192  const scalar current,
193  const scalar luminousIntensity
194  );
195 
196  //- Construct given individual dimension exponents for first
197  // five dimensions
199  (
200  const scalar mass,
201  const scalar length,
202  const scalar time,
203  const scalar temperature,
204  const scalar moles
205  );
206 
207  //- Copy constructor
208  dimensionSet(const dimensionSet& ds);
209 
210  //- Construct and return a clone
212  {
213  return autoPtr<dimensionSet>(new dimensionSet(*this));
214  }
215 
216  //- Construct from Istream
218 
219 
220  // Member Functions
221 
222  //- Return true if it is dimensionless
223  bool dimensionless() const;
224 
225  void reset(const dimensionSet&);
226 
227 
228  // I/O
229 
230  //- Read without the square brackets
232 
233  //- Read
234  Istream& read(Istream& is);
235 
236  //- Write without the square brackets
237  Ostream& writeNoBeginOrEnd(Ostream& os) const;
238 
239  //- Write
240  Ostream& write(Ostream& os) const;
241 
242  //- Write info without the square brackets
244 
245  //- Return info proxy
246  inline InfoProxy<dimensionSet> info() const
247  {
248  return *this;
249  }
250 
251 
252  // Member Operators
253 
254  scalar operator[](const dimensionType) const;
255  scalar& operator[](const dimensionType);
256 
257  scalar operator[](const label) const;
258  scalar& operator[](const label);
259 
260  bool operator==(const dimensionSet&) const;
261  bool operator!=(const dimensionSet&) const;
262 
263  bool operator=(const dimensionSet&) const;
264 
265  bool operator+=(const dimensionSet&) const;
266  bool operator-=(const dimensionSet&) const;
267  bool operator*=(const dimensionSet&);
268  bool operator/=(const dimensionSet&);
269 
270 
271  // Friend Functions
272 
273  friend dimensionSet max(const dimensionSet&, const dimensionSet&);
274  friend dimensionSet min(const dimensionSet&, const dimensionSet&);
276  (
277  const dimensionSet&,
278  const dimensionSet&
279  );
280  friend dimensionSet cmptDivide
281  (
282  const dimensionSet&,
283  const dimensionSet&
284  );
285  friend dimensionSet cmptMag(const dimensionSet&);
286 
287  friend dimensionSet pow(const dimensionSet&, const scalar);
289  friend dimensionSet pow(const dimensionedScalar&, const dimensionSet&);
290 
291  friend dimensionSet sqr(const dimensionSet&);
296  friend dimensionSet pow025(const dimensionSet&);
297 
300  friend dimensionSet mag(const dimensionSet&);
303  friend dimensionSet neg(const dimensionSet&);
304  friend dimensionSet inv(const dimensionSet&);
305 
306  //- Function to check the argument is dimensionless
307  // for transcendental functions
308  friend dimensionSet trans(const dimensionSet&);
309 
310  friend dimensionSet atan2(const dimensionSet&, const dimensionSet&);
311 
312  friend dimensionSet transform(const dimensionSet&);
313 
315  friend dimensionSet perpendicular(const dimensionSet&);
316 
317 
318  // Friend Operators
319 
320  friend dimensionSet operator-(const dimensionSet&);
321 
322  friend dimensionSet operator+
323  (
324  const dimensionSet&,
325  const dimensionSet&
326  );
327 
328  friend dimensionSet operator-
329  (
330  const dimensionSet&,
331  const dimensionSet&
332  );
333 
334  friend dimensionSet operator*
335  (
336  const dimensionSet&,
337  const dimensionSet&
338  );
339 
340  friend dimensionSet operator/
341  (
342  const dimensionSet&,
343  const dimensionSet&
344  );
345 
346  friend dimensionSet operator&
347  (
348  const dimensionSet&,
349  const dimensionSet&
350  );
351 
352  friend dimensionSet operator^
353  (
354  const dimensionSet&,
355  const dimensionSet&
356  );
357 
358  friend dimensionSet operator&&
359  (
360  const dimensionSet&,
361  const dimensionSet&
362  );
363 
364 
365  // IOstream Operators
366 
370 };
371 
372 
373 void writeEntry(Ostream& os, const dimensionSet& value);
374 
375 
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377 
378 } // End namespace Foam
379 
380 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381 
382 #include "dimensionSets.H"
383 
384 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
385 
386 #endif
387 
388 // ************************************************************************* //
System bool.
A helper class for outputting values to Ostream.
Definition: InfoProxy.H:50
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Dimension set for the base types.
Definition: dimensionSet.H:125
friend dimensionSet pow5(const dimensionSet &)
friend dimensionSet operator-(const dimensionSet &)
bool operator+=(const dimensionSet &) const
Definition: dimensionSet.C:199
friend dimensionSet pow(const dimensionSet &, const dimensionedScalar &)
InfoProxy< dimensionSet > info() const
Return info proxy.
Definition: dimensionSet.H:245
bool operator==(const dimensionSet &) const
Definition: dimensionSet.C:161
Ostream & writeInfoNoBeginOrEnd(Ostream &os) const
Write info without the square brackets.
friend dimensionSet transform(const dimensionSet &)
friend dimensionSet min(const dimensionSet &, const dimensionSet &)
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
friend dimensionSet magSqr(const dimensionSet &)
friend dimensionSet pos0(const dimensionSet &)
friend dimensionSet neg(const dimensionSet &)
static const scalar smallExponent
A small exponent with which to perform inexact comparisons.
Definition: dimensionSet.H:156
friend dimensionSet sqr(const dimensionSet &)
friend dimensionSet pow025(const dimensionSet &)
friend dimensionSet atan2(const dimensionSet &, const dimensionSet &)
friend dimensionSet sqrt(const dimensionSet &)
Ostream & writeNoBeginOrEnd(Ostream &os) const
Write without the square brackets.
ClassName("dimensionSet")
bool operator/=(const dimensionSet &)
Definition: dimensionSet.C:235
friend dimensionSet cmptMag(const dimensionSet &)
Ostream & write(Ostream &os) const
Write.
friend Istream & operator>>(Istream &, dimensionSet &)
friend Ostream & operator<<(Ostream &, const dimensionSet &)
friend Ostream & operator<<(Ostream &, const InfoProxy< dimensionSet > &)
friend dimensionSet cmptMultiply(const dimensionSet &, const dimensionSet &)
bool operator-=(const dimensionSet &) const
Definition: dimensionSet.C:213
friend dimensionSet pow4(const dimensionSet &)
scalar operator[](const dimensionType) const
Definition: dimensionSet.C:137
friend dimensionSet mag(const dimensionSet &)
void reset(const dimensionSet &)
Definition: dimensionSet.C:126
friend dimensionSet pow3(const dimensionSet &)
friend dimensionSet sign(const dimensionSet &)
friend dimensionSet pow6(const dimensionSet &)
static const NamedEnum< dimensionType, 7 > dimensionTypeNames_
Names of the dimensions.
Definition: dimensionSet.H:150
Istream & readNoBeginOrEnd(Istream &is)
Read without the square brackets.
bool dimensionless() const
Return true if it is dimensionless.
Definition: dimensionSet.C:107
friend dimensionSet normalised(const dimensionSet &)
autoPtr< dimensionSet > clone() const
Construct and return a clone.
Definition: dimensionSet.H:210
friend dimensionSet pow(const dimensionSet &, const scalar)
dimensionType
Define an enumeration for the names of the dimension exponents.
Definition: dimensionSet.H:139
bool operator*=(const dimensionSet &)
Definition: dimensionSet.C:227
friend dimensionSet trans(const dimensionSet &)
Function to check the argument is dimensionless.
friend dimensionSet cmptDivide(const dimensionSet &, const dimensionSet &)
friend dimensionSet inv(const dimensionSet &)
friend dimensionSet perpendicular(const dimensionSet &)
friend dimensionSet pow(const dimensionedScalar &, const dimensionSet &)
friend dimensionSet max(const dimensionSet &, const dimensionSet &)
bool operator=(const dimensionSet &) const
Definition: dimensionSet.C:185
Istream & read(Istream &is)
Read.
Macro definitions for declaring ClassName(), NamespaceName(), etc.
Useful dimension sets.
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
void cmptDivide(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1, const LagrangianPatchField< Type > &f2)
HashSet< Key, Hash > operator^(const HashSet< Key, Hash > &hash1, const HashSet< Key, Hash > &hash2)
Create a HashSet that only contains unique entries (xor)
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar negPart(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
Istream & operator>>(Istream &, pistonPointEdgeData &)
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
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)
tmp< VolField< Type > > operator&(const fvMatrix< Type > &, const DimensionedField< Type, volMesh > &)
tmp< fvMatrix< Type > > operator*(const volScalarField::Internal &, const fvMatrix< Type > &)
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
tmp< fvMatrix< Type > > operator+(const fvMatrix< Type > &, const fvMatrix< Type > &)
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)
tmp< fvMatrix< Type > > operator/(const fvMatrix< Type > &, const volScalarField::Internal &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename scalarProduct< Type1, Type2 >::type > operator&&(const dimensioned< Type1 > &, const dimensioned< Type2 > &)