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-2025 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 
35 
38 {
39  "fraction",
40  "angle"
41 };
42 
43 
44 namespace Foam
45 {
46  const scalar unitConversion::smallExponent = rootSmall;
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 bool Foam::unitConversion::compare
53 (
54  const unitConversion& a,
55  const unitConversion& b,
56  const bool compareMultiplier
57 )
58 {
59  if (a.any()) return true;
60  if (b.any()) return true;
61  if (a.none()) return false;
62  if (b.none()) return false;
63 
64  // Check the dimensions are the same
65  if (a.dimensions_ != b.dimensions_) return false;
66 
67  // Check the dimensionless units are the same
68  for (int i = 0; i < unitConversion::nDimlessUnits; ++ i)
69  {
70  if
71  (
72  mag(a.exponents_[i] - b.exponents_[i])
74  )
75  {
76  return false;
77  }
78  }
79 
80  // If specified, check the unit conversion factors are the same
81  return !compareMultiplier || a.multiplier_ == b.multiplier_;
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
88 (
89  const dimensionSet& dimensions,
90  const scalar fraction,
91  const scalar angle,
92  const scalar multiplier
93 )
94 :
95  dimensions_(dimensions),
96  multiplier_(multiplier)
97 {
98  exponents_[FRACTION] = fraction;
99  exponents_[ANGLE] = angle;
100 }
101 
102 
104 :
105  dimensions_(dimensions),
106  multiplier_(1)
107 {
108  exponents_[FRACTION] = 0;
109  exponents_[ANGLE] = 0;
110 }
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 {
117  dimensions_.reset(units.dimensions_);
118  for (int i = 0; i < unitConversion::nDimlessUnits; ++ i)
119  {
120  exponents_[i] = units.exponents_[i];
121  }
122  multiplier_ = units.multiplier_;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
127 
129 {
130  if (units.any()) return units;
131  if (units.none()) return unitNone;
132 
133  return
135  (
136  pow(units.dimensions_, exp),
137  units.exponents_[0]*exp,
138  units.exponents_[1]*exp,
139  pow(units.multiplier_, exp)
140  );
141 }
142 
143 
144 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
145 
146 const Foam::unitConversion& Foam::operator+
147 (
148  const unitConversion& a,
149  const unitConversion& b
150 )
151 {
152  if (a.any()) return b;
153  if (b.any()) return a;
154  if (a.none()) return unitNone;
155  if (b.none()) return unitNone;
156 
157  if (!unitConversion::compare(a, b, true))
158  {
160  << "Different units for +" << endl
161  << " units : " << a << " + " << b << endl
162  << abort(FatalError);
163  }
164 
165  return a;
166 }
167 
168 
169 Foam::unitConversion Foam::operator*
170 (
171  const unitConversion& a,
172  const unitConversion& b
173 )
174 {
175  if (a.any()) return a;
176  if (b.any()) return b;
177  if (a.none()) return unitNone;
178  if (b.none()) return unitNone;
179 
180  return
182  (
183  a.dimensions_*b.dimensions_,
184  a.exponents_[0] + b.exponents_[0],
185  a.exponents_[1] + b.exponents_[1],
186  a.multiplier_*b.multiplier_
187  );
188 }
189 
190 
191 Foam::unitConversion Foam::operator/
192 (
193  const unitConversion& a,
194  const unitConversion& b
195 )
196 {
197  if (a.any()) return a;
198  if (b.any()) return b;
199  if (a.none()) return unitNone;
200  if (b.none()) return unitNone;
201 
202  return
204  (
205  a.dimensions_/b.dimensions_,
206  a.exponents_[0] - b.exponents_[0],
207  a.exponents_[1] - b.exponents_[1],
208  a.multiplier_/b.multiplier_
209  );
210 }
211 
212 
213 // ************************************************************************* //
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
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:27
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:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const unitConversion unitNone
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
defineTypeNameAndDebug(combustionModel, 0)
error FatalError