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-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 "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 
84 
87 
90 
94 
95 
96 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
97 
98 Foam::scalar Foam::degToRad(const scalar deg)
99 {
100  return unitDegrees.toStandard(deg);
101 }
102 
103 
104 Foam::scalar Foam::radToDeg(const scalar rad)
105 {
106  return unitDegrees.toUser(rad);
107 }
108 
109 
111 {
113 
114  if (!addedUnitsPtr_)
115  {
117  }
118 
119  addedUnitsPtr_->insert(name, units);
120 
122 }
123 
124 
126 {
127  if (!unitsPtr_)
128  {
130 
131  unitsPtr_->insert("%", unitPercent);
132 
133  unitsPtr_->insert("rad", unitRadians);
134  unitsPtr_->insert("rot", unitRotations);
135  unitsPtr_->insert("deg", unitDegrees);
136 
137  // Get the relevant part of the control dictionary
138  const dictionary& unitSetDict =
139  unitsDict().subDict(unitsDict().lookup<word>("unitSet") + "Coeffs");
140 
141  // Add units from the control dictionary
142  forAllConstIter(dictionary, unitSetDict, iter)
143  {
144  ITstream& is = iter().stream();
145 
146  const unitConversion units(is);
147  const scalar multiplier = pTraits<scalar>(is);
148 
149  const bool ok =
150  unitsPtr_->insert
151  (
152  iter().keyword(),
153  units*unitConversion(dimless, 0, 0, multiplier)
154  );
155 
156  if (!ok)
157  {
159  << "Duplicate unit " << iter().keyword()
160  << " read from dictionary"
161  << exit(FatalIOError);
162  }
163  }
164 
165  // Add programmatically defined units
166  if (addedUnitsPtr_)
167  {
169  {
170  const bool ok = unitsPtr_->insert(iter.key(), iter());
171 
172  if (!ok)
173  {
175  << "Duplicate unit " << iter.key()
176  << " added to dictionary"
177  << exit(FatalIOError);
178  }
179  }
180  }
181  }
182 
183  return *unitsPtr_;
184 }
185 
186 
187 // ************************************************************************* //
bool found
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
An STL-conforming hash table.
Definition: HashTable.H:127
Input token stream.
Definition: ITstream.H:53
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
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
scalar radToDeg(const scalar rad)
Convert radians to degrees.
scalar degToRad(const scalar deg)
Convert degrees to radians.
const unitConversion unitAny
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
HashTable< unitConversion > * addedUnitsPtr_(nullptr)
void deleteDemandDrivenData(DataType *&dataPtr)
const dimensionSet dimless
const unitConversion unitPercent
const unitConversion unitNone
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
dictionary * unitsDictPtr_(nullptr)
const unitConversion unitless
deleteUnitsPtr deleteUnitsPtr_
const unitConversion unitRotations
const unitConversion unitDegrees
const unitConversion unitFraction
Useful unit conversions.