ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/SlaterIntegrals.hpp
Revision: 1718
Committed: Thu May 24 01:29:59 2012 UTC (13 years, 2 months ago) by gezelter
File size: 18540 byte(s)
Log Message:
Adding fluctuating charges, still a work in progress

File Contents

# User Rev Content
1 gezelter 1718 /***************************************************************************
2     * This program is free software; you can redistribute it and/or modify *
3     * it under the terms of the GNU General Public License as published by *
4     * the Free Software Foundation; either version 3 of the License, or *
5     * (at your option) any later version. *
6     * *
7     * This program is distributed in the hope that it will be useful, *
8     * but WITHOUT ANY WARRANTY; without even the implied warranty of *
9     * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
10     * GNU General Public License for more details. *
11     * *
12     * You should have received a copy of the GNU General Public License *
13     * along with this program; if not, see <http://www.gnu.org/licenses/>. *
14     ***************************************************************************/
15    
16     /***
17     * This file was imported from qtpie, found at http://code.google.com/p/qtpie
18     * and was modified minimally for use in OpenMD.
19     *
20     * No author attribution was found in the code, but it presumably is
21     * the work of J. Chen and Todd J. Martinez.
22     *
23     * QTPIE (charge transfer with polarization current equilibration) is
24     * a new charge model, similar to other charge models like QEq,
25     * fluc-q, EEM or ABEEM. Unlike other existing charge models, however,
26     * it is capable of describing both charge transfer and polarization
27     * phenomena. It is also unique for its ability to describe
28     * intermolecular charge transfer at reasonable computational cost.
29     *
30     * Good references to cite when using this code are:
31     *
32     * J. Chen and T. J. Martinez, "QTPIE: Charge transfer with
33     * polarization current equalization. A fluctuating charge model with
34     * correct asymptotics", Chemical Physics Letters 438 (2007),
35     * 315-320. DOI: 10.1016/j.cplett.2007.02.065. Erratum: ibid, 463
36     * (2008), 288. DOI: 10.1016/j.cplett.2008.08.060
37     *
38     * J. Chen, D. Hundertmark and T. J. Martinez, "A unified theoretical
39     * framework for fluctuating-charge models in atom-space and in
40     * bond-space", Journal of Chemical Physics 129 (2008), 214113. DOI:
41     * 10.1063/1.3021400.
42     *
43     * J. Chen and T. J. Martinez, "Charge conservation in
44     * electronegativity equalization and its implications for the
45     * electrostatic properties of fluctuating-charge models", Journal of
46     * Chemical Physics 131 (2009), 044114. DOI: 10.1063/1.3183167
47     *
48     * J. Chen and T. J. Martinez, "The dissociation catastrophe in
49     * fluctuating-charge models and its implications for the concept of
50     * atomic electronegativity", Progress in Theoretical Chemistry and
51     * Physics, to appear. arXiv:0812.1543
52     *
53     * J. Chen, "Theory and applications of fluctuating-charge models",
54     * PhD (chemical physics) thesis, University of Illinois at
55     * Urbana-Champaign, Department of Chemistry, 2009.
56     *
57     * J. Chen and T. J. Martinez, "Size-extensive polarizabilities with
58     * intermolecular charge transfer in a fluctuating-charge model", in
59     * preparation. arXiv:0812.1544
60     */
61    
62     #include "config.h"
63     #include <cmath>
64     #include "math/Factorials.hpp"
65    
66     #ifndef NONBONDED_SLATERINTEGRALS_HPP
67     #define NONBONDED_SLATERINTEGRALS_HPP
68    
69     template <typename T> inline T sqr(T t) { return t*t; }
70     template <typename T> inline T mod(T x, T m)
71     { return x<0 ? m - 1 - ((-x) - 1)%m : x%m; }
72    
73     // #include "parameters.h"
74    
75     /**
76     * @brief Computes Rosen's Guillimer-Zener function A
77     * Computes Rosen's A integral, an auxiliary quantity needed to
78     * compute integrals involving Slater-type orbitals of s symmetry.
79     * \f[
80     * A_n(\alpha) = \int_1^\infty x^n e^{-\alpha x}dx
81     * = \frac{n! e^{-\alpha}}{\alpha^{n+1}}\sum_{\nu=0}^n
82     * \frac{\alpha^\nu}{\nu!}
83     * \f]
84     * @param n - principal quantum number
85     * @param alpha - Slater exponent
86     * @return the value of Rosen's A integral
87     * @note N. Rosen, Phys. Rev., 38 (1931), 255
88     */
89     inline RealType RosenA(int n, RealType a)
90     {
91     RealType RosenA_ = 0.;
92     if (a != 0.)
93     {
94     RealType Term = 1.;
95     RosenA_ = Term;
96     for (unsigned nu=1; nu<=n; nu++)
97     {
98     Term *= a / nu;
99     RosenA_ += Term;
100     }
101     RosenA_ = (RosenA_/Term) * (exp(-a)/a);
102     }
103     return RosenA_;
104     }
105    
106     /**
107     * @brief Computes Rosen's Guillimer-Zener function B
108     * Computes Rosen's B integral, an auxiliary quantity needed to
109     * compute integrals involving Slater-type orbitals of s symmetry.
110     * \f[
111     * B_n(\alpha) = \int_{-1}^1 x^n e^{-\alpha x} dx
112     * = \frac{n!}{\alpha^{n+1}}
113     * \sum_{\nu=0}^n \frac{e^\alpha(-\alpha)^\nu
114     * - e^{-\alpha} \alpha^\nu}{\nu!}
115     * \f]
116     * @param n - principal quantum number
117     * @param alpha - Slater exponent
118     * @return the value of Rosen's B integral
119     * @note N. Rosen, Phys. Rev., 38 (1931), 255
120     */
121     inline RealType RosenB(int n, RealType alpha)
122     {
123     RealType TheSum, Term;
124     RealType RosenB_, PSinhRosenA, PCoshRosenA, PHyperRosenA;
125     int nu;
126     bool IsPositive;
127     if (alpha != 0.)
128     {
129     Term = 1.;
130     TheSum = 1.;
131     IsPositive = true;
132    
133     // These two expressions are (up to constant factors) equivalent
134     // to computing the hyperbolic sine and cosine of a respectively
135     // The series consists of adding up these terms in an alternating fashion
136     PSinhRosenA = exp(alpha) - exp(-alpha);
137     PCoshRosenA = -exp(alpha) - exp(-alpha);
138     TheSum = PSinhRosenA;
139     for (unsigned nu=1; nu<=n; nu++)
140     {
141     if (IsPositive)
142     {
143     PHyperRosenA = PCoshRosenA;
144     IsPositive = false;
145     }
146     else // term to add should be negative
147     {
148     PHyperRosenA = PSinhRosenA;
149     IsPositive = true;
150     }
151     Term *= alpha / nu;
152     TheSum += Term * PHyperRosenA;
153     }
154     RosenB_ = TheSum / (alpha*Term);
155     }
156     else // pathological case of a=0
157     {
158     printf("WARNING, a = 0 in RosenB\n");
159     RosenB_ = (1. - pow(-1., n)) / (n + 1.);
160     }
161     return RosenB_;
162     }
163    
164     /** @brief Computes Rosen's D combinatorial factor
165     * Computes Rosen's D factor, an auxiliary quantity needed to
166     * compute integrals involving Slater-type orbitals of s symmetry.
167     * \f[
168     * RosenD^{mn}_p = \sum_k (-1)^k \frac{m! n!}
169     * {(p-k)!(m-p+k)!(n-k)!k!}
170     * \f]
171     * @return the value of Rosen's D factor
172     * @note N. Rosen, Phys. Rev., 38 (1931), 255
173     */
174     inline RealType RosenD(int m, int n, int p)
175     // [?ERROR?] can't return int
176     {
177     printf("RosenD\n");
178     exit(0);
179     if (m+n+p > maxFact)
180     {
181     printf("Error, arguments exceed maximum factorial computed %d > %d\n", m+n+p, maxFact);
182     ::exit(0);
183     }
184    
185     RealType RosenD_ = 0;
186     for (int k=max(p-m,0); k<=min(n,p); k++)
187     {
188     if (mod(k,2) == 0)
189     RosenD_ += (fact[m] / (fact[p-k] * fact[m-p+k])) * (fact[n] / (fact[n-k] * fact[k]));
190     else
191     RosenD_ -= (fact[m] / ( fact[p-k] * fact[m-p+k])) * (fact[n] / (fact[n-k] * fact[k]));
192     }
193     return RosenD_;
194     }
195    
196     /** @brief Computes Coulomb integral analytically over s-type STOs
197     * Computes the two-center Coulomb integral over Slater-type orbitals of s symmetry.
198     * @param a : Slater zeta exponent of first atom in inverse Bohr (au)
199     * @param b : Slater zeta exponent of second atom in inverse Bohr (au)
200     * @param m : principal quantum number of first atom
201     * @param n : principal quantum number of second atom
202     * @param R : internuclear distance in atomic units (bohr)
203     * @return value of the Coulomb potential energy integral
204     * @note N. Rosen, Phys. Rev., 38 (1931), 255
205     * @note In Rosen's paper, this integral is known as K2.
206     */
207     inline RealType sSTOCoulInt(RealType a, RealType b, int m, int n, RealType R)
208     {
209     int nu, p;
210     RealType x, K2;
211     RealType Factor1, Factor2, Term, OneElectronTerm;
212     RealType eps, epsi;
213    
214     // To speed up calculation, we terminate loop once contributions
215     // to integral fall below the bound, epsilon
216     RealType epsilon = 0.;
217    
218     // x is the argument of the auxiliary RosenA and RosenB functions
219     x = 2. * a * R;
220    
221     // First compute the two-electron component
222     RealType sSTOCoulInt_ = 0.;
223     if (x == 0.) // Pathological case
224     {
225     if ((a==b) && (m==n))
226     {
227     for (int nu=0; nu<=2*n-1; nu++)
228     {
229     K2 = 0.;
230     for (unsigned p=0; p<=2*n+m; p++) K2 += 1. / fact[p];
231     sSTOCoulInt_ += K2 * fact[2*n+m] / fact[m];
232     }
233     sSTOCoulInt_ = 2 * a / (n * fact[2*n]) * sSTOCoulInt_;
234     }
235     else
236     {
237     // Not implemented
238     printf("ERROR, sSTOCoulInt cannot compute from arguments\n");
239     printf("a = %lf b = %lf m = %d n = %d R = %lf\n",a, b, m, n, R);
240     exit(0);
241     }
242     }
243     else
244     {
245     OneElectronTerm = 1./R + pow(x, 2*m)/(fact[2*m]*R)*
246     ((x-2*m)*RosenA(2*m-1,x)-exp(-x)) + sSTOCoulInt_;
247     eps = epsilon / OneElectronTerm;
248     if (a == b)
249     {
250     // Apply Rosen (48)
251     Factor1 = -a*pow(a*R, 2*m)/(n*fact[2*m]);
252     for (int nu=0; nu<=2*n-1; nu++)
253     {
254     Factor2 = (2.*n-nu)/fact[nu]*pow(a*R,nu);
255     epsi = eps / fabs(Factor1 * Factor2);
256     K2 = 0.;
257     for (int p=0; p<=m+(nu-1)/2; p++)
258     {
259     Term = RosenD(2*m-1, nu, 2*p)/(2.*p+1.) *RosenA(2*m+nu-1-2*p,x);
260     K2 += Term;
261     if ((Term > 0) && (Term < epsi)) goto label1;
262     }
263     sSTOCoulInt_ += K2 * Factor2;
264     }
265     label1:
266     sSTOCoulInt_ *= Factor1;
267     }
268     else
269     {
270     Factor1 = -a*pow(a*R,2*m)/(2.*n*fact[2*m]);
271     epsi = eps/fabs(Factor1);
272     if (b == 0.)
273     printf("WARNING: b = 0 in sSTOCoulInt\n");
274     else
275     {
276     // Apply Rosen (54)
277     for (int nu=0; nu<=2*n-1; nu++)
278     {
279     K2 = 0;
280     for (int p=0; p<=2*m+nu-1; p++)
281     K2=K2+RosenD(2*m-1,nu,p)*RosenB(p,R*(a-b))
282     *RosenA(2*m+nu-1-p,R*(a+b));
283     Term = K2*(2*n-nu)/fact[nu]*pow(b*R, nu);
284     sSTOCoulInt_ += Term;
285     if (fabs(Term) < epsi) goto label2;
286     }
287     label2:
288     sSTOCoulInt_ *= Factor1;
289     }
290     }
291     // Now add the one-electron term from Rosen (47) = Rosen (53)
292     sSTOCoulInt_ += OneElectronTerm;
293     }
294     return sSTOCoulInt_;
295     }
296    
297     /**
298     * @brief Computes overlap integral analytically over s-type STOs
299     * Computes the overlap integral over two
300     * Slater-type orbitals of s symmetry.
301     * @param a : Slater zeta exponent of first atom in inverse Bohr (au)
302     * @param b : Slater zeta exponent of second atom in inverse Bohr (au)
303     * @param m : principal quantum number of first atom
304     * @param n : principal quantum number of second atom
305     * @param R : internuclear distance in atomic units (bohr)
306     * @return the value of the sSTOOvInt integral
307     * @note N. Rosen, Phys. Rev., 38 (1931), 255
308     * @note In the Rosen paper, this integral is known as I.
309     */
310     inline RealType sSTOOvInt(RealType a, RealType b, int m, int n, RealType R)
311     {
312     RealType Factor, Term, eps;
313    
314     // To speed up calculation, we terminate loop once contributions
315     // to integral fall below the bound, epsilon
316     RealType epsilon = 0.;
317     RealType sSTOOvInt_ = 0.;
318    
319     if (a == b)
320     {
321     Factor = pow(a*R, m+n+1)/sqrt(fact[2*m]*fact[2*n]);
322     eps = epsilon / fabs(Factor);
323     for (int q=0; q<=(m+n)/2; q++)
324     {
325     Term = RosenD(m,n,2*q)/(2.*q+1.)*RosenA(m+n-2*q,a*R);
326     sSTOOvInt_ += Term;
327     if (fabs(Term) < eps) exit(0);
328     }
329     sSTOOvInt_ *= Factor;
330     }
331     else
332     {
333     Factor = 0.5*pow(a*R, m+0.5)*pow(b*R,n+0.5)
334     /sqrt(fact[2*m]*fact[2*n]);
335     eps = epsilon / fabs(Factor);
336     for (int q=0; q<=m+n; q++)
337     {
338     Term = RosenD(m,n,q)*RosenB(q, R/2.*(a-b))
339     * RosenA(m+n-q,R/2.*(a+b));
340     sSTOOvInt_ += Term;
341     if (fabs(Term) < eps) exit(0);
342     }
343     sSTOOvInt_ *= Factor;
344     }
345     return sSTOOvInt_;
346     }
347    
348     /**
349     * @brief Computes kinetic energy integral analytically over s-type STOs
350     * Computes the overlap integral over two Slater-type orbitals of s symmetry.
351     * @param a : Slater zeta exponent of first atom in inverse Bohr (au)
352     * @param b : Slater zeta exponent of second atom in inverse Bohr (au)
353     * @param m : principal quantum number of first atom
354     * @param n : principal quantum number of second atom
355     * @param R : internuclear distance in atomic units (bohr)
356     * @return the value of the kinetic energy integral
357     * @note N. Rosen, Phys. Rev., 38 (1931), 255
358     * @note untested
359     */
360     inline RealType KinInt(RealType a, RealType b, int m, int n,RealType R)
361     {
362     RealType KinInt_ = -0.5*b*b*sSTOOvInt(a, b, m, n, R);
363     if (n > 0)
364     {
365     KinInt_ += b*b*pow(2*b/(2*b-1),0.5) * sSTOOvInt(a, b, m, n-1, R);
366     if (n > 1) KinInt_ += pow(n*(n-1)/((n-0.5)*(n-1.5)), 0.5)
367     * sSTOOvInt(a, b, m, n-2, R);
368     }
369     return KinInt_;
370     }
371    
372     /**
373     * @brief Computes derivative of Coulomb integral with respect to the interatomic distance
374     * Computes the two-center Coulomb integral over Slater-type orbitals of s symmetry.
375     * @param a: Slater zeta exponent of first atom in inverse Bohr (au)
376     * @param b: Slater zeta exponent of second atom in inverse Bohr (au)
377     * @param m: principal quantum number of first atom
378     * @param n: principal quantum number of second atom
379     * @param R: internuclear distance in atomic units (bohr)
380     * @return the derivative of the Coulomb potential energy integral
381     * @note Derived in QTPIE research notes, May 15 2007
382     */
383     inline RealType sSTOCoulIntGrad(RealType a, RealType b, int m, int n, RealType R)
384     {
385     RealType x, y, z, K2, TheSum;
386     // x is the argument of the auxiliary RosenA and RosenB functions
387     x = 2. * a * R;
388    
389     // First compute the two-electron component
390     RealType sSTOCoulIntGrad_ = 0.;
391     if (x==0) // Pathological case
392     {
393     printf("WARNING: argument given to sSTOCoulIntGrad is 0\n");
394     printf("a = %lf R= %lf\n", a, R);
395     }
396     else
397     {
398     if (a == b)
399     {
400     TheSum = 0.;
401     for (int nu=0; nu<=2*(n-1); nu++)
402     {
403     K2 = 0.;
404     for (int p=0; p<=(m+nu)/2; p++)
405     K2 += RosenD(2*m-1, nu+1, 2*p)/(2*p + 1.) * RosenA(2*m+nu-1-2*p, x);
406     TheSum += (2*n-nu-1)/fact[nu]*pow(a*R, nu) * K2;
407     }
408     sSTOCoulIntGrad_ = -pow(a, 2*m+2)*pow(R, 2*m) /(n*fact[2*m])*TheSum;
409     TheSum = 0.;
410     for (int nu=0; nu<=2*n-1; nu++)
411     {
412     K2 = 0.;
413     for (int p=0; p<=(m+nu-1)/2; p++)
414     K2 += RosenD(2*m-1, nu, 2*p)/(2*p + 1.) * RosenA(2*m+nu-2*p, x);
415     TheSum += (2*n-nu)/fact[nu]*pow(a*R,nu) * K2;
416     }
417     sSTOCoulIntGrad_ += 2*pow(a, 2*m+2)*pow(R, 2*m) /(n*fact[2*m])*TheSum;
418     }
419     else
420     {
421     // Slater exponents are different
422     // First calculate some useful arguments
423     y = R*(a+b);
424     z = R*(a-b);
425     TheSum = 0.;
426     for (int nu=0; nu<=2*n-1; nu++)
427     {
428     K2 = 0.;
429     for (int p=0; p<=2*m+nu; p++)
430     K2 += RosenD(2*m-1, nu+1, p)
431     * RosenB(p,z)*RosenA(2*m+nu-p, y);
432     TheSum += (2*n-nu-1)/fact[nu]*pow(b*R,nu) * K2;
433     }
434     sSTOCoulIntGrad_ = -b*pow(a,2*m+1)*pow(R,2*m)/
435     (2*n*fact[2*m])*TheSum;
436     TheSum = 0.;
437     for (int nu=0; nu<=2*n; nu++)
438     {
439     K2 = 0.;
440     for (int p=0; p<=2*m-1+nu; p++)
441     K2 += RosenD(2*m-1, nu, p)
442     * ((a-b)*RosenB(p+1,z)*RosenA(2*m+nu-p-1, y)
443     +(a+b)*RosenB(p ,z)*RosenA(2*m+nu-p , y));
444     TheSum += (2*n-nu)/fact[nu]*pow(b*R,nu) * K2;
445     }
446     sSTOCoulIntGrad_ += pow(a,2*m+1)*pow(R,2*m)/(2*n*fact[2*m])*TheSum;
447     }
448     // Now add one-electron terms and common term
449     sSTOCoulIntGrad_ = sSTOCoulIntGrad_ - (2.*m+1.)/sqr(R)
450     + 2.*m/R * sSTOCoulInt(a,b,m,n,R)
451     + pow(x,2*m)/(fact[2*m]*sqr(R)) * ((2.*m+1.)*exp(-x)
452     + 2.*m*(1.+2.*m-x)*RosenA(2*m-1,x));
453     }
454     return sSTOCoulIntGrad_;
455     }
456    
457     /**
458     * @brief Computes gradient of overlap integral with respect to the interatomic diatance
459     * Computes the derivative of the overlap integral over two Slater-type orbitals of s symmetry.
460     * @param a: Slater zeta exponent of first atom in inverse Bohr (au)
461     * @param b: Slater zeta exponent of second atom in inverse Bohr (au)
462     * @param m: principal quantum number of first atom
463     * @param n: principal quantum number of second atom
464     * @param R: internuclear distance in atomic units (bohr)
465     * @return the derivative of the sSTOOvInt integral
466     * @note Derived in QTPIE research notes, May 15 2007
467     */
468     inline RealType sSTOOvIntGrad(RealType a, RealType b, int m, int n, RealType R)
469     {
470     RealType w, x, y, z, TheSum;
471    
472     // Calculate first term
473     RealType sSTOOvIntGrad_ = (m+n+1.)/R * sSTOOvInt(a, b, m, n, R);
474    
475     // Calculate remaining terms; answers depend on exponents
476     TheSum = 0.;
477     x = a * R;
478     if (a == b)
479     {
480     for (int q=0; q<=(m+n)/2; q++)
481     TheSum += RosenD(m,n,2*q) / (2*q + 1.) * RosenA(m+n-2*q+1, x);
482     sSTOOvIntGrad_ -= a*pow(x,m+n+1)/ sqrt(fact[2*m]*fact[2*n])*TheSum;
483     }
484     else
485     {
486     w = b*R;
487     y = 0.5*R*(a+b);
488     z = 0.5*R*(a-b);
489     for (int q=0; q<m+n; q++)
490     TheSum = TheSum + RosenD(m,n,q) *
491     ((a-b)*RosenB(q+1,z)*RosenA(m+n-q ,y)
492     +(a+b)*RosenB(q ,z)*RosenA(m+n-q+1,y));
493     sSTOOvIntGrad_ -= 0.25*sqrt((pow(x, 2*m+1)*pow(w, 2*n+1))/(fact[2*m]*fact[2*n]))*TheSum;
494     }
495     return sSTOOvIntGrad_;
496     }
497    
498     /**
499     * @brief Calculates a Slater-type orbital exponent based on the hardness parameters
500     * @param Hardness: chemical hardness in atomic units
501     * @param n: principal quantum number
502     * @note Modified for use with OpenMD by Gezelter and Michalka.
503     */
504     inline RealType getSTOZeta(int n, RealType hardness)
505     {
506     // Approximate the exact value of the constant of proportionality
507     // by its value at a very small distance epsilon
508     // since the exact R = 0 case has not be programmed
509     RealType epsilon = 1.0e-8;
510    
511     // Assign orbital exponent
512     return pow(sSTOCoulInt(1., 1., n, n, epsilon) / hardness, -1./(3. + 2.*n));
513     }
514    
515     #endif

Properties

Name Value
svn:eol-style native