reducedUnits.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-2019 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 "reducedUnits.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 const Foam::scalar Foam::reducedUnits::kb = 1.3806504e-23;
31 
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::reducedUnits::calcRefValues()
36 {
37  if
38  (
39  refTime_ < vSmall
40  || refLength_ < vSmall
41  || refMass_ < vSmall
42  )
43  {
45  << "One of more reference values too small for floating point "
46  << "calculation: "
47  << "refTime_ = " << refTime_
48  << ", refLength = " << refTemp_
49  << ", refMass = " << refMass_
50  << nl << abort(FatalError);
51  }
52 
53  refEnergy_ = refLength_*refLength_*refMass_/(refTime_*refTime_);
54 
55  refTemp_ = refEnergy_ / kb;
56 
57  refForce_ = refEnergy_/refLength_;
58 
59  refVelocity_ = Foam::sqrt(refEnergy_/refMass_);
60 
61  refVolume_ = Foam::pow(refLength_,3.0);
62 
63  refPressure_ = refEnergy_/refVolume_;
64 
65  refMassDensity_ = refMass_/refVolume_;
66 
67  refNumberDensity_ = 1.0/refVolume_;
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
74 :
75  refLength_(1e-9),
76  refTime_(1e-12),
77  refMass_(1.660538782e-27)
78 {
79  calcRefValues();
80 }
81 
82 
84 (
85  scalar refLength,
86  scalar refTime,
87  scalar refMass
88 )
89 :
90  refLength_(refLength),
91  refTime_(refTime),
92  refMass_(refMass)
93 {
94  calcRefValues();
95 }
96 
97 
99 :
100  refLength_(),
101  refTime_(),
102  refMass_()
103 {
104  setRefValues(reducedUnitsDict);
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
109 
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 (
118  scalar refLength,
119  scalar refTime,
120  scalar refMass
121 )
122 {
123  refLength_ = refLength;
124 
125  refTime_ = refTime;
126 
127  refMass_ = refMass;
128 
129  calcRefValues();
130 }
131 
132 
134 (
135  const IOdictionary& reducedUnitsDict
136 )
137 {
138  refLength_ = reducedUnitsDict.template lookup<scalar>("refLength");
139 
140  refTime_ = reducedUnitsDict.template lookup<scalar>("refTime");
141 
142  refMass_ = reducedUnitsDict.template lookup<scalar>("refMass");
143 
144  calcRefValues();
145 }
146 
147 
148 // ************************************************************************* //
scalar refLength() const
Definition: reducedUnitsI.H:28
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
dimensionedScalar sqrt(const dimensionedScalar &ds)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
scalar refTime() const
Definition: reducedUnitsI.H:34
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const char nl
Definition: Ostream.H:260
scalar refMass() const
Definition: reducedUnitsI.H:40
reducedUnits()
Construct with no argument, uses default values:
Definition: reducedUnits.C:73
~reducedUnits()
Destructor.
Definition: reducedUnits.C:110
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static const scalar kb
Static data someStaticData.
Definition: reducedUnits.H:100
void setRefValues(scalar refLength, scalar refTime, scalar refMass)
Definition: reducedUnits.C:117
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105