unitConversions.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) 2011-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 "unitConversions.H"
27 #include "demandDrivenData.H"
28 #include "dictionary.H"
29 
30 using namespace Foam::constant::mathematical;
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
38 
40 {
41  if (!unitsDictPtr_)
42  {
43  dictionary* cachedPtr = nullptr;
44 
46  (
48  (
49  debug::controlDict().found("UnitConversions")
50  ? "UnitConversions"
51  : debug::controlDict().found("DimensionSets")
52  ? "DimensionSets"
53  : "UnitConversions",
54  cachedPtr
55  )
56  );
57  }
58 
59  return *unitsDictPtr_;
60 }
61 
64 
65 // Delete the above data at the end of the run
67 {
69  {
73  }
74 };
75 
77 
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83 namespace Foam
84 {
85 
87 {
88  return dimensionSet(0, 0, 0, 0, 0);
89 }
90 
92 {
93  return unitConversion(makeDimless(), 0, 0, 1);
94 }
95 
97 {
98  return unitConversion(makeDimless(), 0, 0, 0);
99 }
101 {
102  return unitConversion(makeDimless(), 0, 0, -1);
103 }
104 
106 {
107  return unitConversion(makeDimless(), 1, 0, 1);
108 }
110 {
111  return unitConversion(makeDimless(), 1, 0, 0.01);
112 }
113 
115 {
116  return unitConversion(makeDimless(), 0, 1, 1);
117 }
119 {
120  return unitConversion(makeDimless(), 0, 1, 2*pi);
121 }
123 {
124  return unitConversion(makeDimless(), 0, 1, pi/180);
125 }
126 
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
133 
136 
139 
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 Foam::scalar Foam::degToRad(const scalar deg)
148 {
149  return unitDegrees.toStandard(deg);
150 }
151 
152 
153 Foam::scalar Foam::radToDeg(const scalar rad)
154 {
155  return unitDegrees.toUser(rad);
156 }
157 
158 
160 {
162 
163  if (!addedUnitsPtr_)
164  {
166  }
167 
168  addedUnitsPtr_->insert(name, units);
169 
171 }
172 
173 
175 {
176  if (!unitsPtr_)
177  {
179 
180  unitsPtr_->insert("%", makeUnitPercent());
181 
182  unitsPtr_->insert("rad", makeUnitRadians());
183  unitsPtr_->insert("rot", makeUnitRotations());
184  unitsPtr_->insert("deg", makeUnitDegrees());
185 
186  // Get the relevant part of the control dictionary
187  const dictionary& unitSetDict =
188  unitsDict().subDict(unitsDict().lookup<word>("unitSet") + "Coeffs");
189 
190  // Add units from the control dictionary
191  forAllConstIter(dictionary, unitSetDict, iter)
192  {
193  ITstream& is = iter().stream();
194 
195  const unitConversion units(is);
196  const scalar multiplier = pTraits<scalar>(is);
197 
198  const bool ok =
199  unitsPtr_->insert
200  (
201  iter().keyword(),
202  units*unitConversion(dimless, 0, 0, multiplier)
203  );
204 
205  if (!ok)
206  {
208  << "Duplicate unit " << iter().keyword()
209  << " read from dictionary"
210  << exit(FatalIOError);
211  }
212  }
213 
214  // Add programmatically defined units
215  if (addedUnitsPtr_)
216  {
218  {
219  const bool ok = unitsPtr_->insert(iter.key(), iter());
220 
221  if (!ok)
222  {
224  << "Duplicate unit " << iter.key()
225  << " added to dictionary"
226  << exit(FatalIOError);
227  }
228  }
229  }
230  }
231 
232  return *unitsPtr_;
233 }
234 
235 
236 // ************************************************************************* //
bool found
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
An STL-conforming hash table.
Definition: HashTable.H:127
Input token stream.
Definition: ITstream.H:53
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:849
Dimension set for the base types.
Definition: dimensionSet.H:125
Traits class for primitives.
Definition: pTraits.H:53
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
T toUser(const T &) const
Convert a value to user units.
T toStandard(const T &) const
Convert a value to standard units.
A class for handling words, derived from string.
Definition: word.H:62
Template functions to aid in the implementation of demand driven data.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
mathematical constants.
dictionary & switchSet(const char *subDictName, dictionary *&subDictPtr)
Internal function to lookup a sub-dictionary from controlDict.
Definition: debug.C:164
dictionary & controlDict()
The central control dictionary.
Definition: debug.C:120
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
unitConversion makeUnitNone()
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
dimensionSet makeDimless()
const unitConversion unitAny
unitConversion makeUnitDegrees()
unitConversion makeUnitless()
HashTable< unitConversion > * addedUnitsPtr_(nullptr)
unitConversion makeUnitRotations()
void deleteDemandDrivenData(DataType *&dataPtr)
const dimensionSet dimless
const unitConversion unitPercent
const unitConversion unitNone
unitConversion makeUnitAny()
const HashTable< unitConversion > & units()
Get the table of unit conversions.
const dictionary & unitsDict()
void addUnits(const word &name, const unitConversion &units)
Add a unit conversion.
HashTable< unitConversion > * unitsPtr_(nullptr)
IOerror FatalIOError
const unitConversion unitRadians
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary * unitsDictPtr_(nullptr)
const unitConversion unitless
deleteUnitsPtr deleteUnitsPtr_
const unitConversion unitRotations
unitConversion makeUnitRadians()
const unitConversion unitDegrees
unitConversion makeUnitPercent()
unitConversion makeUnitFraction()
const unitConversion unitFraction
Useful unit conversions.