1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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 "pairPotential.H"
28 #include "energyScalingFunction.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 defineTypeNameAndDebug(pairPotential, 0);
35 defineRunTimeSelectionTable(pairPotential, dictionary);
39 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 void Foam::pairPotential::scaleEnergy(scalar& e, const scalar r) const
46 esfPtr_ = energyScalingFunction::New
48 name_, pairPotentialProperties_, *this
52 esfPtr_->scaleEnergy(e, r);
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 Foam::pairPotential::pairPotential
61 const dictionary& pairPotentialProperties
65 pairPotentialProperties_(pairPotentialProperties),
66 rCut_(readScalar(pairPotentialProperties_.lookup("rCut"))),
67 rCutSqr_(rCut_*rCut_),
68 rMin_(readScalar(pairPotentialProperties_.lookup("rMin"))),
69 dr_(readScalar(pairPotentialProperties_.lookup("dr"))),
73 writeTables_(Switch(pairPotentialProperties_.lookup("writeTables")))
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
79 void Foam::pairPotential::setLookupTables()
81 label N = label((rCut_ - rMin_)/dr_) + 1;
83 forceLookup_.setSize(N);
85 energyLookup_.setSize(N);
87 forAll(forceLookup_, k)
89 energyLookup_[k] = scaledEnergy(k*dr_ + rMin_);
91 forceLookup_[k] = -energyDerivative((k*dr_ + rMin_), true);
96 Foam::scalar Foam::pairPotential::force(const scalar r) const
98 scalar k_rIJ = (r - rMin_)/dr_;
100 label k = label(k_rIJ);
104 FatalErrorIn("pairPotential.C") << nl
105 << "r less than rMin in pair potential " << name_ << nl
106 << abort(FatalError);
110 (k_rIJ - k)*forceLookup_[k+1]
111 + (k + 1 - k_rIJ)*forceLookup_[k];
117 Foam::List< Foam::Pair< Foam::scalar > >
118 Foam::pairPotential::forceTable() const
120 List<Pair<scalar> > forceTab(forceLookup_.size());
122 forAll(forceLookup_,k)
124 forceTab[k].first() = rMin_ + k*dr_;
126 forceTab[k].second() = forceLookup_[k];
133 Foam::scalar Foam::pairPotential::energy(const scalar r) const
135 scalar k_rIJ = (r - rMin_)/dr_;
137 label k = label(k_rIJ);
141 FatalErrorIn("pairPotential.C") << nl
142 << "r less than rMin in pair potential " << name_ << nl
143 << abort(FatalError);
147 (k_rIJ - k)*energyLookup_[k+1]
148 + (k + 1 - k_rIJ)*energyLookup_[k];
154 Foam::List< Foam::Pair< Foam::scalar > >
155 Foam::pairPotential::energyTable() const
157 List<Pair<scalar> > energyTab(energyLookup_.size());
159 forAll(energyLookup_,k)
161 energyTab[k].first() = rMin_ + k*dr_;
163 energyTab[k].second() = energyLookup_[k];
170 Foam::scalar Foam::pairPotential::scaledEnergy(const scalar r) const
172 scalar e = unscaledEnergy(r);
180 Foam::scalar Foam::pairPotential::energyDerivative
183 const bool scaledEnergyDerivative
186 // Local quadratic fit to energy: E = a0 + a1*r + a2*r^2
187 // Differentiate to give f = -dE/dr = -a1 - 2*a2*r
195 if (scaledEnergyDerivative)
197 Ea = scaledEnergy(ra);
198 Ef = scaledEnergy(rf);
199 Eb = scaledEnergy(rb);
203 Ea = unscaledEnergy(ra);
204 Ef = unscaledEnergy(rf);
205 Eb = unscaledEnergy(rb);
208 scalar denominator = (ra - rf)*(ra - rb)*(rf - rb);
212 rb*rb*(Ea - Ef) + ra*ra*(Ef - Eb) + rf*rf*(Eb - Ea)
217 rb*(Ef - Ea) + rf*(Ea - Eb) + ra*(Eb - Ef)
220 return a1 + 2.0*a2*r;
224 bool Foam::pairPotential::read(const dictionary& pairPotentialProperties)
226 pairPotentialProperties_ = pairPotentialProperties;
232 // ************************************************************************* //