unitConversion.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) 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 "unitConversion.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
33 
34  template<>
36  {
37  "fraction",
38  "angle"
39  };
40 }
41 
42 
45 
46 
47 namespace Foam
48 {
49  const scalar unitConversion::smallExponent = rootSmall;
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 bool Foam::unitConversion::compare
56 (
57  const unitConversion& a,
58  const unitConversion& b,
59  const bool compareMultiplier
60 )
61 {
62  if (a.any()) return true;
63  if (b.any()) return true;
64  if (a.none()) return false;
65  if (b.none()) return false;
66 
67  // Check the dimensions are the same
68  if (a.dimensions_ != b.dimensions_) return false;
69 
70  // Check the dimensionless units are the same
71  for (int i = unitConversion::nDimlessUnits; i < 0; ++ i)
72  {
73  if
74  (
75  mag(a.exponents_[i] - b.exponents_[i])
77  )
78  {
79  return false;
80  }
81  }
82 
83  // If specified, check the unit conversion factors are the same
84  return !compareMultiplier || a.multiplier_ == b.multiplier_;
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
91 (
92  const dimensionSet& dimensions,
93  const scalar fraction,
94  const scalar angle,
95  const scalar multiplier
96 )
97 :
98  dimensions_(dimensions),
99  multiplier_(multiplier)
100 {
101  exponents_[FRACTION] = fraction;
102  exponents_[ANGLE] = angle;
103 }
104 
105 
107 :
108  dimensions_(dimensions),
109  multiplier_(1)
110 {
111  exponents_[FRACTION] = 0;
112  exponents_[ANGLE] = 0;
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  dimensions_.reset(units.dimensions_);
121  for (int i = 0; i < unitConversion::nDimlessUnits; ++ i)
122  {
123  exponents_[i] = units.exponents_[i];
124  }
125  multiplier_ = units.multiplier_;
126 }
127 
128 
129 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
130 
132 {
133  if (units.any()) return units;
134  if (units.none()) return unitNone;
135 
136  return
138  (
139  pow(units.dimensions_, exp),
140  units.exponents_[0]*exp,
141  units.exponents_[1]*exp,
142  pow(units.multiplier_, exp)
143  );
144 }
145 
146 
147 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
148 
149 const Foam::unitConversion& Foam::operator+
150 (
151  const unitConversion& a,
152  const unitConversion& b
153 )
154 {
155  if (a.any()) return b;
156  if (b.any()) return a;
157  if (a.none()) return unitNone;
158  if (b.none()) return unitNone;
159 
160  if (!unitConversion::compare(a, b, true))
161  {
163  << "Different units for +" << endl
164  << " units : "
165  << a << " {" << a.multiplier_ << "} + "
166  << b << " {" << b.multiplier_ << "} + " << endl
167  << abort(FatalError);
168  }
169 
170  return a;
171 }
172 
173 
174 Foam::unitConversion Foam::operator*
175 (
176  const unitConversion& a,
177  const unitConversion& b
178 )
179 {
180  if (a.any()) return a;
181  if (b.any()) return b;
182  if (a.none()) return unitNone;
183  if (b.none()) return unitNone;
184 
185  return
187  (
188  a.dimensions_*b.dimensions_,
189  a.exponents_[0] + b.exponents_[0],
190  a.exponents_[1] + b.exponents_[1],
191  a.multiplier_*b.multiplier_
192  );
193 }
194 
195 
196 Foam::unitConversion Foam::operator/
197 (
198  const unitConversion& a,
199  const unitConversion& b
200 )
201 {
202  if (a.any()) return a;
203  if (b.any()) return b;
204  if (a.none()) return unitNone;
205  if (b.none()) return unitNone;
206 
207  return
209  (
210  a.dimensions_/b.dimensions_,
211  a.exponents_[0] - b.exponents_[0],
212  a.exponents_[1] - b.exponents_[1],
213  a.multiplier_/b.multiplier_
214  );
215 }
216 
217 
218 // ************************************************************************* //
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
Dimension set for the base types.
Definition: dimensionSet.H:125
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
static const NamedEnum< dimlessUnitType, 2 > dimlessUnitTypeNames_
Names of the dimensionless units.
static const scalar smallExponent
A small exponent with which to perform inexact comparisons.
unitConversion(const dimensionSet &, const scalar fraction, const scalar angle, const scalar multiplier)
Construct from components.
void reset(const unitConversion &)
Reset the unit conversion.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField & b
Definition: createFields.H:25
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const unitConversion unitNone
const HashTable< unitConversion > & units()
Get the table of unit conversions.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
error FatalError