unitConversion.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) 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::unitConversion
26 
27 Description
28  Unit conversion structure. Contains the associated dimensions and the
29  multiplier with which to convert values.
30 
31 SourceFiles
32  unitConversion.C
33  unitConversionIO.C
34  unitConversionI.H
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef unitConversion_H
39 #define unitConversion_H
40 
41 #include "dimensionSet.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 // Forward declaration of classes
49 class Istream;
50 class Ostream;
51 template<class Type> class List;
52 template<class Type> class Field;
53 
54 // Forward declaration of friend functions and operators
55 class unitConversion;
56 unitConversion pow(const unitConversion&, const scalar);
57 const unitConversion& operator+(const unitConversion&, const unitConversion&);
58 unitConversion operator*(const unitConversion&, const unitConversion&);
59 unitConversion operator/(const unitConversion&, const unitConversion&);
60 Istream& operator>>(Istream&, unitConversion&);
61 Ostream& operator<<(Ostream&, const unitConversion&);
62 Ostream& operator<<(Ostream&, const InfoProxy<unitConversion>&);
63 
64 /*---------------------------------------------------------------------------*\
65  Class unitConversion Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class unitConversion
69 {
70 public:
71 
72  // Member constants
73 
74  //- Define an enumeration for the number of dimensionless units
75  enum
76  {
78  };
79 
80  //- Define an enumeration for the names of the dimensionless unit
81  // exponents
82  enum dimlessUnitType
83  {
84  FRACTION, // fraction %
85  ANGLE // angle rad, rot, deg
86  };
87 
88  //- Names of the dimensionless units
90 
91 
92  // Static Data Members
93 
94  //- A small exponent with which to perform inexact comparisons
95  static const scalar smallExponent;
96 
97 
98 private:
99 
100  // Private Data
101 
102  //- The dimensions
103  dimensionSet dimensions_;
104 
105  //- Array of dimensionless unit exponents
106  scalar exponents_[nDimlessUnits];
107 
108  //- The conversion multiplier with which to multiply quantities
109  // in this unit in order to recover the quantity in standard units
110  scalar multiplier_;
111 
112 
113  // Private Member Functions
114 
115  //- Compare two unit conversions and return if they are the same.
116  // An additional control determines whether the multiplier is compared
117  // as well as the dimensions and dimensionless units.
118  static bool compare
119  (
120  const unitConversion&,
121  const unitConversion&,
122  const bool compareMultiplier
123  );
124 
125 
126 public:
127 
128  // Static Data Members
129 
130  //- Run time type information
131  ClassName("unitConversion");
132 
133 
134  // Constructors
135 
136  //- Construct from components
138  (
139  const dimensionSet&,
140  const scalar fraction,
141  const scalar angle,
142  const scalar multiplier
143  );
144 
145  //- Construct from a dimension set. No dimensionless units. Unity
146  // multiplier.
148 
149  //- Copy constructor
150  unitConversion(const unitConversion&) = default;
151 
152  //- Move constructor
153  unitConversion(unitConversion&&) = default;
154 
155  //- Construct from stream
156  unitConversion(Istream& is);
157 
158 
159  // Member Functions
160 
161  //- Access the dimensions
162  inline const dimensionSet& dimensions() const;
163 
164  //- Convert a value to standard units
165  template<class T>
166  T toStandard(const T&) const;
167 
168  //- Convert a pair of values to standard units
169  template<class T>
170  Pair<T> toStandard(const Pair<T>&) const;
171 
172  //- Convert a list of values to standard units
173  template<class T>
174  List<T> toStandard(const List<T>&) const;
175 
176  //- Convert a field of values to standard units
177  template<class T>
178  tmp<Field<T>> toStandard(const Field<T>&) const;
179 
180  //- Convert a tmp field of values to standard units
181  template<class T>
182  tmp<Field<T>> toStandard(const tmp<Field<T>>&) const;
183 
184  //- Convert a value to standard units
185  template<class T>
186  void makeStandard(T&) const;
187 
188  //- Convert a pair of values to standard units
189  template<class T>
190  void makeStandard(Pair<T>&) const;
191 
192  //- Convert a list of values to standard units
193  template<class T>
194  void makeStandard(List<T>&) const;
195 
196  //- Convert a value to user units
197  template<class T>
198  T toUser(const T&) const;
199 
200  //- Convert a pair of values to user units
201  template<class T>
202  Pair<T> toUser(const Pair<T>&) const;
203 
204  //- Convert a list of values to user units
205  template<class T>
206  List<T> toUser(const List<T>&) const;
207 
208  //- Convert a field of values to user units
209  template<class T>
210  tmp<Field<T>> toUser(const Field<T>&) const;
211 
212  //- Convert a tmp field of values to user units
213  template<class T>
214  tmp<Field<T>> toUser(const tmp<Field<T>>&) const;
215 
216  //- Return whether this is the "any" unit. I.e., the case where
217  // dimensions and dimensionless units are not checked, and any
218  // conversion is permitted.
219  inline bool any() const;
220 
221  //- Return whether this is the "none" unit. I.e., the case where unit
222  // conversions are prohibited.
223  inline bool none() const;
224 
225  //- Return whether this unit is standard. I.e., is its multiplier one?
226  inline bool standard() const;
227 
228  //- Reset the unit conversion
229  void reset(const unitConversion&);
230 
231  //- Update
232  void read(const word& keyword, const dictionary&);
233 
234  //- Update
235  void read(Istream& is);
236 
237  //- Update
238  void read(const word& keyword, const dictionary&, Istream& is);
239 
240  //- Update if found in the dictionary
241  bool readIfPresent(const word& keyword, const dictionary&);
242 
243  //- Update if found on the stream
244  bool readIfPresent(Istream& is);
245 
246  //- Update if found on the dictionary stream
247  bool readIfPresent(const word& keyword, const dictionary&, Istream& is);
248 
249  //- Return info proxy
250  inline InfoProxy<unitConversion> info() const
251  {
252  return *this;
253  }
254 
255 
256  // Member Operators
257 
258  //- Disallow default bitwise assignment
259  void operator=(const unitConversion&) = delete;
260 
261  //- Disallow default bitwise move assignment
262  void operator=(const unitConversion&&) = delete;
263 
264 
265  // Friend Functions
266 
267  //- Raise to a power
268  friend unitConversion pow(const unitConversion&, const scalar);
269 
270 
271  // Friend Operators
272 
273  //- Combine
274  friend const unitConversion& operator+
275  (
276  const unitConversion&,
277  const unitConversion&
278  );
279 
280  //- Multiply
281  friend unitConversion operator*
282  (
283  const unitConversion&,
284  const unitConversion&
285  );
286 
287  //- Divide
288  friend unitConversion operator/
289  (
290  const unitConversion&,
291  const unitConversion&
292  );
293 
294 
295  // IOstream Operators
296 
297  //- Read from stream
299 
300  //- Write to stream
301  friend Ostream& operator<<(Ostream&, const unitConversion&);
302 
303  //- Write info to stream
305 };
306 
307 
308 // Global Functions
309 
310 //- Write a type with a given unit conversion
311 template<class Type>
312 void writeEntry(Ostream& os, const unitConversion&, const Type& t);
313 
314 
315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
316 
317 } // End namespace Foam
318 
319 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
320 
321 #include "unitConversionI.H"
322 
323 #ifdef NoRepository
324  #include "unitConversionTemplates.C"
325 #endif
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 #endif
330 
331 // ************************************************************************* //
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
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An ordered pair of two objects of type <Type> with first() and second() elements.
Definition: Pair.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
A class for managing temporary objects.
Definition: tmp.H:55
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
friend Ostream & operator<<(Ostream &, const unitConversion &)
Write to stream.
friend unitConversion pow(const unitConversion &, const scalar)
Raise to a power.
ClassName("unitConversion")
Run time type information.
friend Istream & operator>>(Istream &, unitConversion &)
Read from stream.
bool readIfPresent(const word &keyword, const dictionary &)
Update if found in the dictionary.
void operator=(const unitConversion &)=delete
Disallow default bitwise assignment.
T toUser(const T &) const
Convert a value to user units.
bool none() const
Return whether this is the "none" unit. I.e., the case where unit.
static const NamedEnum< dimlessUnitType, 2 > dimlessUnitTypeNames_
Names of the dimensionless units.
const dimensionSet & dimensions() const
Access the dimensions.
static const scalar smallExponent
A small exponent with which to perform inexact comparisons.
void makeStandard(T &) const
Convert a value to standard units.
unitConversion(const dimensionSet &, const scalar fraction, const scalar angle, const scalar multiplier)
Construct from components.
InfoProxy< unitConversion > info() const
Return info proxy.
bool standard() const
Return whether this unit is standard. I.e., is its multiplier one?
void read(const word &keyword, const dictionary &)
Update.
void reset(const unitConversion &)
Reset the unit conversion.
T toStandard(const T &) const
Convert a value to standard units.
dimlessUnitType
Define an enumeration for the names of the dimensionless unit.
bool any() const
Return whether this is the "any" unit. I.e., the case where.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
Istream & operator>>(Istream &, pistonPointEdgeData &)
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 > &)
tmp< fvMatrix< Type > > operator/(const fvMatrix< Type > &, const volScalarField::Internal &)