pairPotential.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-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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 "pairPotential.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
37}
38
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42void Foam::pairPotential::scaleEnergy(scalar& e, const scalar r) const
43{
44 if (!esfPtr_)
45 {
47 (
49 ).ptr();
50 }
51
53}
54
55
56// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57
59(
60 const word& name,
61 const dictionary& pairPotentialProperties
62)
63:
64 name_(name),
65 pairPotentialProperties_(pairPotentialProperties),
66 rCut_(pairPotentialProperties_.get<scalar>("rCut")),
67 rCutSqr_(rCut_*rCut_),
68 rMin_(pairPotentialProperties_.get<scalar>("rMin")),
69 dr_(pairPotentialProperties_.get<scalar>("dr")),
70 forceLookup_(0),
71 energyLookup_(0),
72 esfPtr_(nullptr),
73 writeTables_(pairPotentialProperties_.get<bool>("writeTables"))
74{}
75
76
77// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78
80{
81 label N = label((rCut_ - rMin_)/dr_) + 1;
82
83 forceLookup_.setSize(N);
84
85 energyLookup_.setSize(N);
86
87 forAll(forceLookup_, k)
88 {
89 energyLookup_[k] = scaledEnergy(k*dr_ + rMin_);
90
91 forceLookup_[k] = -energyDerivative((k*dr_ + rMin_), true);
92 }
93}
94
95
96Foam::scalar Foam::pairPotential::force(const scalar r) const
97{
98 scalar k_rIJ = (r - rMin_)/dr_;
99
100 label k = label(k_rIJ);
101
102 if (k < 0)
103 {
105 << "r less than rMin in pair potential " << name_ << nl
106 << abort(FatalError);
107 }
108
109 scalar f =
110 (k_rIJ - k)*forceLookup_[k+1]
111 + (k + 1 - k_rIJ)*forceLookup_[k];
112
113 return f;
114}
115
116
119{
120 List<Pair<scalar>> forceTab(forceLookup_.size());
121
122 forAll(forceLookup_,k)
123 {
124 forceTab[k].first() = rMin_ + k*dr_;
125
126 forceTab[k].second() = forceLookup_[k];
127 }
128
129 return forceTab;
130}
131
132
133Foam::scalar Foam::pairPotential::energy(const scalar r) const
134{
135 scalar k_rIJ = (r - rMin_)/dr_;
136
137 label k = label(k_rIJ);
138
139 if (k < 0)
140 {
142 << "r less than rMin in pair potential " << name_ << nl
143 << abort(FatalError);
144 }
145
146 scalar e =
147 (k_rIJ - k)*energyLookup_[k+1]
148 + (k + 1 - k_rIJ)*energyLookup_[k];
149
150 return e;
151}
152
153
156{
157 List<Pair<scalar>> energyTab(energyLookup_.size());
158
159 forAll(energyLookup_,k)
160 {
161 energyTab[k].first() = rMin_ + k*dr_;
162
163 energyTab[k].second() = energyLookup_[k];
164 }
165
166 return energyTab;
167}
168
169
170Foam::scalar Foam::pairPotential::scaledEnergy(const scalar r) const
171{
172 scalar e = unscaledEnergy(r);
173
174 scaleEnergy(e, r);
175
176 return e;
177}
178
179
181(
182 const scalar r,
183 const bool scaledEnergyDerivative
184) const
185{
186 // Local quadratic fit to energy: E = a0 + a1*r + a2*r^2
187 // Differentiate to give f = -dE/dr = -a1 - 2*a2*r
188
189 scalar ra = r - dr_;
190 scalar rf = r;
191 scalar rb = r + dr_;
192
193 scalar Ea, Ef, Eb;
194
195 if (scaledEnergyDerivative)
196 {
197 Ea = scaledEnergy(ra);
198 Ef = scaledEnergy(rf);
199 Eb = scaledEnergy(rb);
200 }
201 else
202 {
203 Ea = unscaledEnergy(ra);
204 Ef = unscaledEnergy(rf);
205 Eb = unscaledEnergy(rb);
206 }
207
208 scalar denominator = (ra - rf)*(ra - rb)*(rf - rb);
209
210 scalar a1 =
211 (
212 rb*rb*(Ea - Ef) + ra*ra*(Ef - Eb) + rf*rf*(Eb - Ea)
213 ) / denominator;
214
215 scalar a2 =
216 (
217 rb*(Ef - Ea) + rf*(Ea - Eb) + ra*(Eb - Ef)
218 ) / denominator;
219
220 return a1 + 2.0*a2*r;
221}
222
223
224bool Foam::pairPotential::read(const dictionary& pairPotentialProperties)
225{
226 pairPotentialProperties_ = pairPotentialProperties;
227
228 return true;
229}
230
231
232// ************************************************************************* //
scalar Ea(const scalar p, const scalar T) const
Definition: HtoEthermo.H:32
label k
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual void scaleEnergy(scalar &e, const scalar r) const =0
scalar energy(const scalar r) const
energyScalingFunction * esfPtr_
Definition: pairPotential.H:79
scalar energyDerivative(const scalar r, const bool scaledEnergyDerivative=true) const
void scaleEnergy(scalar &e, const scalar r) const
Definition: pairPotential.C:42
List< Pair< scalar > > energyTable() const
List< Pair< scalar > > forceTable() const
scalar scaledEnergy(const scalar r) const
dictionary pairPotentialProperties_
Definition: pairPotential.H:68
Base class for film (stress-based) force models.
Definition: force.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
bool
Definition: EEqn.H:20
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Namespace for OpenFOAM.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
const Vector< label > N(dict.get< Vector< label > >("N"))