ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/SlaterIntegrals.hpp
Revision: 1730
Committed: Wed May 30 17:19:13 2012 UTC (13 years ago) by gezelter
File size: 18559 byte(s)
Log Message:
Fixed compilation issues

File Contents

# Content
1 /***************************************************************************
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 <cstdlib>
65 #include "math/Factorials.hpp"
66
67 #ifndef NONBONDED_SLATERINTEGRALS_HPP
68 #define NONBONDED_SLATERINTEGRALS_HPP
69
70 template <typename T> inline T sqr(T t) { return t*t; }
71 template <typename T> inline T mod(T x, T m)
72 { return x<0 ? m - 1 - ((-x) - 1)%m : x%m; }
73
74 // #include "parameters.h"
75
76 /**
77 * @brief Computes Rosen's Guillimer-Zener function A
78 * Computes Rosen's A integral, an auxiliary quantity needed to
79 * compute integrals involving Slater-type orbitals of s symmetry.
80 * \f[
81 * A_n(\alpha) = \int_1^\infty x^n e^{-\alpha x}dx
82 * = \frac{n! e^{-\alpha}}{\alpha^{n+1}}\sum_{\nu=0}^n
83 * \frac{\alpha^\nu}{\nu!}
84 * \f]
85 * @param n - principal quantum number
86 * @param alpha - Slater exponent
87 * @return the value of Rosen's A integral
88 * @note N. Rosen, Phys. Rev., 38 (1931), 255
89 */
90 inline RealType RosenA(int n, RealType a)
91 {
92 RealType RosenA_ = 0.;
93 if (a != 0.)
94 {
95 RealType Term = 1.;
96 RosenA_ = Term;
97 for (unsigned nu=1; nu<=n; nu++)
98 {
99 Term *= a / nu;
100 RosenA_ += Term;
101 }
102 RosenA_ = (RosenA_/Term) * (exp(-a)/a);
103 }
104 return RosenA_;
105 }
106
107 /**
108 * @brief Computes Rosen's Guillimer-Zener function B
109 * Computes Rosen's B integral, an auxiliary quantity needed to
110 * compute integrals involving Slater-type orbitals of s symmetry.
111 * \f[
112 * B_n(\alpha) = \int_{-1}^1 x^n e^{-\alpha x} dx
113 * = \frac{n!}{\alpha^{n+1}}
114 * \sum_{\nu=0}^n \frac{e^\alpha(-\alpha)^\nu
115 * - e^{-\alpha} \alpha^\nu}{\nu!}
116 * \f]
117 * @param n - principal quantum number
118 * @param alpha - Slater exponent
119 * @return the value of Rosen's B integral
120 * @note N. Rosen, Phys. Rev., 38 (1931), 255
121 */
122 inline RealType RosenB(int n, RealType alpha)
123 {
124 RealType TheSum, Term;
125 RealType RosenB_, PSinhRosenA, PCoshRosenA, PHyperRosenA;
126 int nu;
127 bool IsPositive;
128 if (alpha != 0.)
129 {
130 Term = 1.;
131 TheSum = 1.;
132 IsPositive = true;
133
134 // These two expressions are (up to constant factors) equivalent
135 // to computing the hyperbolic sine and cosine of a respectively
136 // The series consists of adding up these terms in an alternating fashion
137 PSinhRosenA = exp(alpha) - exp(-alpha);
138 PCoshRosenA = -exp(alpha) - exp(-alpha);
139 TheSum = PSinhRosenA;
140 for (unsigned nu=1; nu<=n; nu++)
141 {
142 if (IsPositive)
143 {
144 PHyperRosenA = PCoshRosenA;
145 IsPositive = false;
146 }
147 else // term to add should be negative
148 {
149 PHyperRosenA = PSinhRosenA;
150 IsPositive = true;
151 }
152 Term *= alpha / nu;
153 TheSum += Term * PHyperRosenA;
154 }
155 RosenB_ = TheSum / (alpha*Term);
156 }
157 else // pathological case of a=0
158 {
159 printf("WARNING, a = 0 in RosenB\n");
160 RosenB_ = (1. - pow(-1., n)) / (n + 1.);
161 }
162 return RosenB_;
163 }
164
165 /** @brief Computes Rosen's D combinatorial factor
166 * Computes Rosen's D factor, an auxiliary quantity needed to
167 * compute integrals involving Slater-type orbitals of s symmetry.
168 * \f[
169 * RosenD^{mn}_p = \sum_k (-1)^k \frac{m! n!}
170 * {(p-k)!(m-p+k)!(n-k)!k!}
171 * \f]
172 * @return the value of Rosen's D factor
173 * @note N. Rosen, Phys. Rev., 38 (1931), 255
174 */
175 inline RealType RosenD(int m, int n, int p)
176 // [?ERROR?] can't return int
177 {
178 printf("RosenD\n");
179 exit(0);
180 if (m+n+p > maxFact)
181 {
182 printf("Error, arguments exceed maximum factorial computed %d > %d\n", m+n+p, maxFact);
183 ::exit(0);
184 }
185
186 RealType RosenD_ = 0;
187 for (int k=max(p-m,0); k<=min(n,p); k++)
188 {
189 if (mod(k,2) == 0)
190 RosenD_ += (fact[m] / (fact[p-k] * fact[m-p+k])) * (fact[n] / (fact[n-k] * fact[k]));
191 else
192 RosenD_ -= (fact[m] / ( fact[p-k] * fact[m-p+k])) * (fact[n] / (fact[n-k] * fact[k]));
193 }
194 return RosenD_;
195 }
196
197 /** @brief Computes Coulomb integral analytically over s-type STOs
198 * Computes the two-center Coulomb integral over Slater-type orbitals of s symmetry.
199 * @param a : Slater zeta exponent of first atom in inverse Bohr (au)
200 * @param b : Slater zeta exponent of second atom in inverse Bohr (au)
201 * @param m : principal quantum number of first atom
202 * @param n : principal quantum number of second atom
203 * @param R : internuclear distance in atomic units (bohr)
204 * @return value of the Coulomb potential energy integral
205 * @note N. Rosen, Phys. Rev., 38 (1931), 255
206 * @note In Rosen's paper, this integral is known as K2.
207 */
208 inline RealType sSTOCoulInt(RealType a, RealType b, int m, int n, RealType R)
209 {
210 int nu, p;
211 RealType x, K2;
212 RealType Factor1, Factor2, Term, OneElectronTerm;
213 RealType eps, epsi;
214
215 // To speed up calculation, we terminate loop once contributions
216 // to integral fall below the bound, epsilon
217 RealType epsilon = 0.;
218
219 // x is the argument of the auxiliary RosenA and RosenB functions
220 x = 2. * a * R;
221
222 // First compute the two-electron component
223 RealType sSTOCoulInt_ = 0.;
224 if (x == 0.) // Pathological case
225 {
226 if ((a==b) && (m==n))
227 {
228 for (int nu=0; nu<=2*n-1; nu++)
229 {
230 K2 = 0.;
231 for (unsigned p=0; p<=2*n+m; p++) K2 += 1. / fact[p];
232 sSTOCoulInt_ += K2 * fact[2*n+m] / fact[m];
233 }
234 sSTOCoulInt_ = 2 * a / (n * fact[2*n]) * sSTOCoulInt_;
235 }
236 else
237 {
238 // Not implemented
239 printf("ERROR, sSTOCoulInt cannot compute from arguments\n");
240 printf("a = %lf b = %lf m = %d n = %d R = %lf\n",a, b, m, n, R);
241 exit(0);
242 }
243 }
244 else
245 {
246 OneElectronTerm = 1./R + pow(x, 2*m)/(fact[2*m]*R)*
247 ((x-2*m)*RosenA(2*m-1,x)-exp(-x)) + sSTOCoulInt_;
248 eps = epsilon / OneElectronTerm;
249 if (a == b)
250 {
251 // Apply Rosen (48)
252 Factor1 = -a*pow(a*R, 2*m)/(n*fact[2*m]);
253 for (int nu=0; nu<=2*n-1; nu++)
254 {
255 Factor2 = (2.*n-nu)/fact[nu]*pow(a*R,nu);
256 epsi = eps / fabs(Factor1 * Factor2);
257 K2 = 0.;
258 for (int p=0; p<=m+(nu-1)/2; p++)
259 {
260 Term = RosenD(2*m-1, nu, 2*p)/(2.*p+1.) *RosenA(2*m+nu-1-2*p,x);
261 K2 += Term;
262 if ((Term > 0) && (Term < epsi)) goto label1;
263 }
264 sSTOCoulInt_ += K2 * Factor2;
265 }
266 label1:
267 sSTOCoulInt_ *= Factor1;
268 }
269 else
270 {
271 Factor1 = -a*pow(a*R,2*m)/(2.*n*fact[2*m]);
272 epsi = eps/fabs(Factor1);
273 if (b == 0.)
274 printf("WARNING: b = 0 in sSTOCoulInt\n");
275 else
276 {
277 // Apply Rosen (54)
278 for (int nu=0; nu<=2*n-1; nu++)
279 {
280 K2 = 0;
281 for (int p=0; p<=2*m+nu-1; p++)
282 K2=K2+RosenD(2*m-1,nu,p)*RosenB(p,R*(a-b))
283 *RosenA(2*m+nu-1-p,R*(a+b));
284 Term = K2*(2*n-nu)/fact[nu]*pow(b*R, nu);
285 sSTOCoulInt_ += Term;
286 if (fabs(Term) < epsi) goto label2;
287 }
288 label2:
289 sSTOCoulInt_ *= Factor1;
290 }
291 }
292 // Now add the one-electron term from Rosen (47) = Rosen (53)
293 sSTOCoulInt_ += OneElectronTerm;
294 }
295 return sSTOCoulInt_;
296 }
297
298 /**
299 * @brief Computes overlap integral analytically over s-type STOs
300 * Computes the overlap integral over two
301 * Slater-type orbitals of s symmetry.
302 * @param a : Slater zeta exponent of first atom in inverse Bohr (au)
303 * @param b : Slater zeta exponent of second atom in inverse Bohr (au)
304 * @param m : principal quantum number of first atom
305 * @param n : principal quantum number of second atom
306 * @param R : internuclear distance in atomic units (bohr)
307 * @return the value of the sSTOOvInt integral
308 * @note N. Rosen, Phys. Rev., 38 (1931), 255
309 * @note In the Rosen paper, this integral is known as I.
310 */
311 inline RealType sSTOOvInt(RealType a, RealType b, int m, int n, RealType R)
312 {
313 RealType Factor, Term, eps;
314
315 // To speed up calculation, we terminate loop once contributions
316 // to integral fall below the bound, epsilon
317 RealType epsilon = 0.;
318 RealType sSTOOvInt_ = 0.;
319
320 if (a == b)
321 {
322 Factor = pow(a*R, m+n+1)/sqrt(fact[2*m]*fact[2*n]);
323 eps = epsilon / fabs(Factor);
324 for (int q=0; q<=(m+n)/2; q++)
325 {
326 Term = RosenD(m,n,2*q)/(2.*q+1.)*RosenA(m+n-2*q,a*R);
327 sSTOOvInt_ += Term;
328 if (fabs(Term) < eps) exit(0);
329 }
330 sSTOOvInt_ *= Factor;
331 }
332 else
333 {
334 Factor = 0.5*pow(a*R, m+0.5)*pow(b*R,n+0.5)
335 /sqrt(fact[2*m]*fact[2*n]);
336 eps = epsilon / fabs(Factor);
337 for (int q=0; q<=m+n; q++)
338 {
339 Term = RosenD(m,n,q)*RosenB(q, R/2.*(a-b))
340 * RosenA(m+n-q,R/2.*(a+b));
341 sSTOOvInt_ += Term;
342 if (fabs(Term) < eps) exit(0);
343 }
344 sSTOOvInt_ *= Factor;
345 }
346 return sSTOOvInt_;
347 }
348
349 /**
350 * @brief Computes kinetic energy integral analytically over s-type STOs
351 * Computes the overlap integral over two Slater-type orbitals of s symmetry.
352 * @param a : Slater zeta exponent of first atom in inverse Bohr (au)
353 * @param b : Slater zeta exponent of second atom in inverse Bohr (au)
354 * @param m : principal quantum number of first atom
355 * @param n : principal quantum number of second atom
356 * @param R : internuclear distance in atomic units (bohr)
357 * @return the value of the kinetic energy integral
358 * @note N. Rosen, Phys. Rev., 38 (1931), 255
359 * @note untested
360 */
361 inline RealType KinInt(RealType a, RealType b, int m, int n,RealType R)
362 {
363 RealType KinInt_ = -0.5*b*b*sSTOOvInt(a, b, m, n, R);
364 if (n > 0)
365 {
366 KinInt_ += b*b*pow(2*b/(2*b-1),0.5) * sSTOOvInt(a, b, m, n-1, R);
367 if (n > 1) KinInt_ += pow(n*(n-1)/((n-0.5)*(n-1.5)), 0.5)
368 * sSTOOvInt(a, b, m, n-2, R);
369 }
370 return KinInt_;
371 }
372
373 /**
374 * @brief Computes derivative of Coulomb integral with respect to the interatomic distance
375 * Computes the two-center Coulomb integral over Slater-type orbitals of s symmetry.
376 * @param a: Slater zeta exponent of first atom in inverse Bohr (au)
377 * @param b: Slater zeta exponent of second atom in inverse Bohr (au)
378 * @param m: principal quantum number of first atom
379 * @param n: principal quantum number of second atom
380 * @param R: internuclear distance in atomic units (bohr)
381 * @return the derivative of the Coulomb potential energy integral
382 * @note Derived in QTPIE research notes, May 15 2007
383 */
384 inline RealType sSTOCoulIntGrad(RealType a, RealType b, int m, int n, RealType R)
385 {
386 RealType x, y, z, K2, TheSum;
387 // x is the argument of the auxiliary RosenA and RosenB functions
388 x = 2. * a * R;
389
390 // First compute the two-electron component
391 RealType sSTOCoulIntGrad_ = 0.;
392 if (x==0) // Pathological case
393 {
394 printf("WARNING: argument given to sSTOCoulIntGrad is 0\n");
395 printf("a = %lf R= %lf\n", a, R);
396 }
397 else
398 {
399 if (a == b)
400 {
401 TheSum = 0.;
402 for (int nu=0; nu<=2*(n-1); nu++)
403 {
404 K2 = 0.;
405 for (int p=0; p<=(m+nu)/2; p++)
406 K2 += RosenD(2*m-1, nu+1, 2*p)/(2*p + 1.) * RosenA(2*m+nu-1-2*p, x);
407 TheSum += (2*n-nu-1)/fact[nu]*pow(a*R, nu) * K2;
408 }
409 sSTOCoulIntGrad_ = -pow(a, 2*m+2)*pow(R, 2*m) /(n*fact[2*m])*TheSum;
410 TheSum = 0.;
411 for (int nu=0; nu<=2*n-1; nu++)
412 {
413 K2 = 0.;
414 for (int p=0; p<=(m+nu-1)/2; p++)
415 K2 += RosenD(2*m-1, nu, 2*p)/(2*p + 1.) * RosenA(2*m+nu-2*p, x);
416 TheSum += (2*n-nu)/fact[nu]*pow(a*R,nu) * K2;
417 }
418 sSTOCoulIntGrad_ += 2*pow(a, 2*m+2)*pow(R, 2*m) /(n*fact[2*m])*TheSum;
419 }
420 else
421 {
422 // Slater exponents are different
423 // First calculate some useful arguments
424 y = R*(a+b);
425 z = R*(a-b);
426 TheSum = 0.;
427 for (int nu=0; nu<=2*n-1; nu++)
428 {
429 K2 = 0.;
430 for (int p=0; p<=2*m+nu; p++)
431 K2 += RosenD(2*m-1, nu+1, p)
432 * RosenB(p,z)*RosenA(2*m+nu-p, y);
433 TheSum += (2*n-nu-1)/fact[nu]*pow(b*R,nu) * K2;
434 }
435 sSTOCoulIntGrad_ = -b*pow(a,2*m+1)*pow(R,2*m)/
436 (2*n*fact[2*m])*TheSum;
437 TheSum = 0.;
438 for (int nu=0; nu<=2*n; nu++)
439 {
440 K2 = 0.;
441 for (int p=0; p<=2*m-1+nu; p++)
442 K2 += RosenD(2*m-1, nu, p)
443 * ((a-b)*RosenB(p+1,z)*RosenA(2*m+nu-p-1, y)
444 +(a+b)*RosenB(p ,z)*RosenA(2*m+nu-p , y));
445 TheSum += (2*n-nu)/fact[nu]*pow(b*R,nu) * K2;
446 }
447 sSTOCoulIntGrad_ += pow(a,2*m+1)*pow(R,2*m)/(2*n*fact[2*m])*TheSum;
448 }
449 // Now add one-electron terms and common term
450 sSTOCoulIntGrad_ = sSTOCoulIntGrad_ - (2.*m+1.)/sqr(R)
451 + 2.*m/R * sSTOCoulInt(a,b,m,n,R)
452 + pow(x,2*m)/(fact[2*m]*sqr(R)) * ((2.*m+1.)*exp(-x)
453 + 2.*m*(1.+2.*m-x)*RosenA(2*m-1,x));
454 }
455 return sSTOCoulIntGrad_;
456 }
457
458 /**
459 * @brief Computes gradient of overlap integral with respect to the interatomic diatance
460 * Computes the derivative of the overlap integral over two Slater-type orbitals of s symmetry.
461 * @param a: Slater zeta exponent of first atom in inverse Bohr (au)
462 * @param b: Slater zeta exponent of second atom in inverse Bohr (au)
463 * @param m: principal quantum number of first atom
464 * @param n: principal quantum number of second atom
465 * @param R: internuclear distance in atomic units (bohr)
466 * @return the derivative of the sSTOOvInt integral
467 * @note Derived in QTPIE research notes, May 15 2007
468 */
469 inline RealType sSTOOvIntGrad(RealType a, RealType b, int m, int n, RealType R)
470 {
471 RealType w, x, y, z, TheSum;
472
473 // Calculate first term
474 RealType sSTOOvIntGrad_ = (m+n+1.)/R * sSTOOvInt(a, b, m, n, R);
475
476 // Calculate remaining terms; answers depend on exponents
477 TheSum = 0.;
478 x = a * R;
479 if (a == b)
480 {
481 for (int q=0; q<=(m+n)/2; q++)
482 TheSum += RosenD(m,n,2*q) / (2*q + 1.) * RosenA(m+n-2*q+1, x);
483 sSTOOvIntGrad_ -= a*pow(x,m+n+1)/ sqrt(fact[2*m]*fact[2*n])*TheSum;
484 }
485 else
486 {
487 w = b*R;
488 y = 0.5*R*(a+b);
489 z = 0.5*R*(a-b);
490 for (int q=0; q<m+n; q++)
491 TheSum = TheSum + RosenD(m,n,q) *
492 ((a-b)*RosenB(q+1,z)*RosenA(m+n-q ,y)
493 +(a+b)*RosenB(q ,z)*RosenA(m+n-q+1,y));
494 sSTOOvIntGrad_ -= 0.25*sqrt((pow(x, 2*m+1)*pow(w, 2*n+1))/(fact[2*m]*fact[2*n]))*TheSum;
495 }
496 return sSTOOvIntGrad_;
497 }
498
499 /**
500 * @brief Calculates a Slater-type orbital exponent based on the hardness parameters
501 * @param Hardness: chemical hardness in atomic units
502 * @param n: principal quantum number
503 * @note Modified for use with OpenMD by Gezelter and Michalka.
504 */
505 inline RealType getSTOZeta(int n, RealType hardness)
506 {
507 // Approximate the exact value of the constant of proportionality
508 // by its value at a very small distance epsilon
509 // since the exact R = 0 case has not be programmed
510 RealType epsilon = 1.0e-8;
511
512 // Assign orbital exponent
513 return pow(sSTOCoulInt(1., 1., n, n, epsilon) / hardness, -1./(3. + 2.*n));
514 }
515
516 #endif

Properties

Name Value
svn:eol-style native