FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

pairPotential.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 2008-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "pairPotential.H"
00027 #include <potential/energyScalingFunction.H>
00028 
00029 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033     defineTypeNameAndDebug(pairPotential, 0);
00034     defineRunTimeSelectionTable(pairPotential, dictionary);
00035 }
00036 
00037 
00038 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00039 
00040 
00041 void Foam::pairPotential::scaleEnergy(scalar& e, const scalar r) const
00042 {
00043     if (!esfPtr_)
00044     {
00045         esfPtr_ = energyScalingFunction::New
00046         (
00047             name_, pairPotentialProperties_, *this
00048         ).ptr();
00049     }
00050 
00051     esfPtr_->scaleEnergy(e, r);
00052 }
00053 
00054 
00055 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00056 
00057 Foam::pairPotential::pairPotential
00058 (
00059     const word& name,
00060     const dictionary& pairPotentialProperties
00061 )
00062 :
00063     name_(name),
00064     pairPotentialProperties_(pairPotentialProperties),
00065     rCut_(readScalar(pairPotentialProperties_.lookup("rCut"))),
00066     rCutSqr_(rCut_*rCut_),
00067     rMin_(readScalar(pairPotentialProperties_.lookup("rMin"))),
00068     dr_(readScalar(pairPotentialProperties_.lookup("dr"))),
00069     forceLookup_(0),
00070     energyLookup_(0),
00071     esfPtr_(NULL),
00072     writeTables_(Switch(pairPotentialProperties_.lookup("writeTables")))
00073 {}
00074 
00075 
00076 // * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * * //
00077 
00078 void Foam::pairPotential::setLookupTables()
00079 {
00080     label N = label((rCut_ - rMin_)/dr_) + 1;
00081 
00082     forceLookup_.setSize(N);
00083 
00084     energyLookup_.setSize(N);
00085 
00086     forAll(forceLookup_, k)
00087     {
00088         energyLookup_[k] = scaledEnergy(k*dr_ + rMin_);
00089 
00090         forceLookup_[k] = -energyDerivative((k*dr_ + rMin_), true);
00091     }
00092 }
00093 
00094 
00095 Foam::scalar Foam::pairPotential::force(const scalar r) const
00096 {
00097     scalar k_rIJ = (r - rMin_)/dr_;
00098 
00099     label k = label(k_rIJ);
00100 
00101     if (k < 0)
00102     {
00103         FatalErrorIn("pairPotential.C") << nl
00104             << "r less than rMin in pair potential " << name_ << nl
00105             << abort(FatalError);
00106     }
00107 
00108     scalar f =
00109         (k_rIJ - k)*forceLookup_[k+1]
00110       + (k + 1 - k_rIJ)*forceLookup_[k];
00111 
00112     return f;
00113 }
00114 
00115 
00116 Foam::List< Foam::Pair< Foam::scalar > >
00117 Foam::pairPotential::forceTable() const
00118 {
00119     List<Pair<scalar> > forceTab(forceLookup_.size());
00120 
00121     forAll(forceLookup_,k)
00122     {
00123         forceTab[k].first() = rMin_ + k*dr_;
00124 
00125         forceTab[k].second() = forceLookup_[k];
00126     }
00127 
00128     return forceTab;
00129 }
00130 
00131 
00132 Foam::scalar Foam::pairPotential::energy(const scalar r) const
00133 {
00134     scalar k_rIJ = (r - rMin_)/dr_;
00135 
00136     label k = label(k_rIJ);
00137 
00138     if (k < 0)
00139     {
00140         FatalErrorIn("pairPotential.C") << nl
00141             << "r less than rMin in pair potential " << name_ << nl
00142             << abort(FatalError);
00143     }
00144 
00145     scalar e =
00146         (k_rIJ - k)*energyLookup_[k+1]
00147       + (k + 1 - k_rIJ)*energyLookup_[k];
00148 
00149     return e;
00150 }
00151 
00152 
00153 Foam::List< Foam::Pair< Foam::scalar > >
00154     Foam::pairPotential::energyTable() const
00155 {
00156     List<Pair<scalar> > energyTab(energyLookup_.size());
00157 
00158     forAll(energyLookup_,k)
00159     {
00160         energyTab[k].first() = rMin_ + k*dr_;
00161 
00162         energyTab[k].second() = energyLookup_[k];
00163     }
00164 
00165     return energyTab;
00166 }
00167 
00168 
00169 Foam::scalar Foam::pairPotential::scaledEnergy(const scalar r) const
00170 {
00171     scalar e = unscaledEnergy(r);
00172 
00173     scaleEnergy(e, r);
00174 
00175     return e;
00176 }
00177 
00178 
00179 Foam::scalar Foam::pairPotential::energyDerivative
00180 (
00181     const scalar r,
00182     const bool scaledEnergyDerivative
00183 ) const
00184 {
00185     // Local quadratic fit to energy: E = a0 + a1*r + a2*r^2
00186     // Differentiate to give f = -dE/dr = -a1 - 2*a2*r
00187 
00188     scalar ra = r - dr_;
00189     scalar rf = r;
00190     scalar rb = r + dr_;
00191 
00192     scalar Ea, Ef, Eb;
00193 
00194     if (scaledEnergyDerivative)
00195     {
00196         Ea = scaledEnergy(ra);
00197         Ef = scaledEnergy(rf);
00198         Eb = scaledEnergy(rb);
00199     }
00200     else
00201     {
00202         Ea = unscaledEnergy(ra);
00203         Ef = unscaledEnergy(rf);
00204         Eb = unscaledEnergy(rb);
00205     }
00206 
00207     scalar denominator = (ra - rf)*(ra - rb)*(rf - rb);
00208 
00209     scalar a1 =
00210     (
00211         rb*rb*(Ea - Ef) + ra*ra*(Ef - Eb) + rf*rf*(Eb - Ea)
00212     ) / denominator;
00213 
00214     scalar a2 =
00215     (
00216         rb*(Ef - Ea) + rf*(Ea - Eb) + ra*(Eb - Ef)
00217     ) / denominator;
00218 
00219     return a1 + 2.0*a2*r;
00220 }
00221 
00222 
00223 bool Foam::pairPotential::read(const dictionary& pairPotentialProperties)
00224 {
00225     pairPotentialProperties_ = pairPotentialProperties;
00226 
00227     return true;
00228 }
00229 
00230 
00231 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines