initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / lagrangian / molecularDynamics / molecule / reducedUnits / reducedUnits.C
blob6836c2bc277e3f13ee476df76d7aca0b250a5102
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2008-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "reducedUnits.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 const  Foam::scalar Foam::reducedUnits::kb = 1.3806504e-23;
34 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
36 void Foam::reducedUnits::calcRefValues()
38     if
39     (
40         refTime_ < VSMALL
41      || refLength_ < VSMALL
42      || refMass_ < VSMALL
43     )
44     {
45         FatalErrorIn("Foam::reducedUnits::calcRefValues() ")
46             << "One of more referencence values too small for floating point "
47             << "calculation: "
48             << "refTime_ = " << refTime_
49             << ", refLength = " << refTemp_
50             << ", refMass = " << refMass_
51             << nl << abort(FatalError);
52     }
54     refEnergy_ = refLength_*refLength_*refMass_/(refTime_*refTime_);
56     refTemp_ = refEnergy_ / kb;
58     refForce_ = refEnergy_/refLength_;
60     refVelocity_ = Foam::sqrt(refEnergy_/refMass_);
62     refVolume_ = Foam::pow(refLength_,3.0);
64     refPressure_ = refEnergy_/refVolume_;
66     refMassDensity_ = refMass_/refVolume_;
68     refNumberDensity_ = 1.0/refVolume_;
72 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
74 Foam::reducedUnits::reducedUnits()
76     refLength_(1e-9),
77     refTime_(1e-12),
78     refMass_(1.660538782e-27)
80     calcRefValues();
84 Foam::reducedUnits::reducedUnits
86     scalar refLength,
87     scalar refTime,
88     scalar refMass
91     refLength_(refLength),
92     refTime_(refTime),
93     refMass_(refMass)
95     calcRefValues();
99 Foam::reducedUnits::reducedUnits(const IOdictionary& reducedUnitsDict)
101     refLength_(),
102     refTime_(),
103     refMass_()
105     setRefValues(reducedUnitsDict);
109 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
111 Foam::reducedUnits::~reducedUnits()
115 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
117 void Foam::reducedUnits::setRefValues
119     scalar refLength,
120     scalar refTime,
121     scalar refMass
124     refLength_ = refLength;
126     refTime_ = refTime;
128     refMass_ = refMass;
130     calcRefValues();
134 void Foam::reducedUnits::setRefValues
136     const IOdictionary& reducedUnitsDict
139     refLength_ = readScalar(reducedUnitsDict.lookup("refLength"));
141     refTime_ = readScalar(reducedUnitsDict.lookup("refTime"));
143     refMass_  = readScalar(reducedUnitsDict.lookup("refMass"));
145     calcRefValues();
149 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
151 void Foam::reducedUnits::operator=(const reducedUnits& rhs)
153     // Check for assignment to self
154     if (this == &rhs)
155     {
156         FatalErrorIn
157         (
158             "Foam::reducedUnits::operator=(const Foam::reducedUnits&)"
159         )   << "Attempted assignment to self"
160             << abort(FatalError);
161     }
165 // ************************************************************************* //