reducedUnits.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "reducedUnits.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 const Foam::scalar Foam::reducedUnits::kb = 1.3806504e-23;
33 
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 void Foam::reducedUnits::calcRefValues()
38 {
39  if
40  (
41  refTime_ < VSMALL
42  || refLength_ < VSMALL
43  || refMass_ < VSMALL
44  )
45  {
47  << "One of more referencence values too small for floating point "
48  << "calculation: "
49  << "refTime_ = " << refTime_
50  << ", refLength = " << refTemp_
51  << ", refMass = " << refMass_
52  << nl << abort(FatalError);
53  }
54 
55  refEnergy_ = refLength_*refLength_*refMass_/(refTime_*refTime_);
56 
57  refTemp_ = refEnergy_ / kb;
58 
59  refForce_ = refEnergy_/refLength_;
60 
61  refVelocity_ = Foam::sqrt(refEnergy_/refMass_);
62 
63  refVolume_ = Foam::pow(refLength_,3.0);
64 
65  refPressure_ = refEnergy_/refVolume_;
66 
67  refMassDensity_ = refMass_/refVolume_;
68 
69  refNumberDensity_ = 1.0/refVolume_;
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74 
76 :
77  refLength_(1e-9),
78  refTime_(1e-12),
79  refMass_(1.660538782e-27)
80 {
81  calcRefValues();
82 }
83 
84 
86 (
87  scalar refLength,
88  scalar refTime,
89  scalar refMass
90 )
91 :
92  refLength_(refLength),
93  refTime_(refTime),
94  refMass_(refMass)
95 {
96  calcRefValues();
97 }
98 
99 
101 :
102  refLength_(),
103  refTime_(),
104  refMass_()
105 {
106  setRefValues(reducedUnitsDict);
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 (
114  scalar refLength,
115  scalar refTime,
116  scalar refMass
117 )
118 {
119  refLength_ = refLength;
120  refTime_ = refTime;
121  refMass_ = refMass;
122 
123  calcRefValues();
124 }
125 
126 
128 (
129  const IOdictionary& reducedUnitsDict
130 )
131 {
132  reducedUnitsDict.readEntry("refLength", refLength_);
133  reducedUnitsDict.readEntry("refTime", refTime_);
134  reducedUnitsDict.readEntry("refMass", refMass_);
135 
136  calcRefValues();
137 }
138 
139 
140 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
reducedUnits.H
Foam::reducedUnits::reducedUnits
reducedUnits()
Construct with no argument, uses default values:
Definition: reducedUnits.C:75
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::reducedUnits::setRefValues
void setRefValues(scalar refLength, scalar refTime, scalar refMass)
Definition: reducedUnits.C:113
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::reducedUnits::kb
static const scalar kb
Static data someStaticData.
Definition: reducedUnits.H:100