ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/Electrostatic.cpp
Revision: 1922
Committed: Mon Aug 5 13:41:15 2013 UTC (11 years, 11 months ago) by gezelter
File size: 54531 byte(s)
Log Message:
Fixed a few ewald bugs

File Contents

# User Rev Content
1 gezelter 1502 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9     * 1. Redistributions of source code must retain the above copyright
10     * notice, this list of conditions and the following disclaimer.
11     *
12     * 2. Redistributions in binary form must reproduce the above copyright
13     * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31     *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 gezelter 1587 * [2] Fennell & Gezelter, J. Chem. Phys. 124 234104 (2006).
38 gezelter 1879 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 1502 */
42    
43     #include <stdio.h>
44     #include <string.h>
45    
46     #include <cmath>
47 gezelter 1915 #include <numeric>
48 gezelter 1502 #include "nonbonded/Electrostatic.hpp"
49     #include "utils/simError.h"
50     #include "types/NonBondedInteractionType.hpp"
51 gezelter 1710 #include "types/FixedChargeAdapter.hpp"
52 gezelter 1720 #include "types/FluctuatingChargeAdapter.hpp"
53 gezelter 1710 #include "types/MultipoleAdapter.hpp"
54 gezelter 1535 #include "io/Globals.hpp"
55 gezelter 1718 #include "nonbonded/SlaterIntegrals.hpp"
56     #include "utils/PhysicalConstants.hpp"
57 gezelter 1767 #include "math/erfc.hpp"
58 gezelter 1879 #include "math/SquareMatrix.hpp"
59 gezelter 1915 #include "primitives/Molecule.hpp"
60 gezelter 1502
61 gezelter 1915
62 gezelter 1502 namespace OpenMD {
63    
64     Electrostatic::Electrostatic(): name_("Electrostatic"), initialized_(false),
65 gezelter 1587 forceField_(NULL), info_(NULL),
66     haveCutoffRadius_(false),
67     haveDampingAlpha_(false),
68     haveDielectric_(false),
69 gezelter 1879 haveElectroSplines_(false)
70 gezelter 1587 {}
71 gezelter 1502
72     void Electrostatic::initialize() {
73 gezelter 1587
74 gezelter 1584 Globals* simParams_ = info_->getSimParams();
75 gezelter 1535
76 gezelter 1528 summationMap_["HARD"] = esm_HARD;
77 gezelter 1616 summationMap_["NONE"] = esm_HARD;
78 gezelter 1528 summationMap_["SWITCHING_FUNCTION"] = esm_SWITCHING_FUNCTION;
79     summationMap_["SHIFTED_POTENTIAL"] = esm_SHIFTED_POTENTIAL;
80     summationMap_["SHIFTED_FORCE"] = esm_SHIFTED_FORCE;
81 gezelter 1879 summationMap_["TAYLOR_SHIFTED"] = esm_TAYLOR_SHIFTED;
82 gezelter 1528 summationMap_["REACTION_FIELD"] = esm_REACTION_FIELD;
83     summationMap_["EWALD_FULL"] = esm_EWALD_FULL;
84     summationMap_["EWALD_PME"] = esm_EWALD_PME;
85     summationMap_["EWALD_SPME"] = esm_EWALD_SPME;
86     screeningMap_["DAMPED"] = DAMPED;
87     screeningMap_["UNDAMPED"] = UNDAMPED;
88    
89 gezelter 1502 // these prefactors convert the multipole interactions into kcal / mol
90     // all were computed assuming distances are measured in angstroms
91     // Charge-Charge, assuming charges are measured in electrons
92     pre11_ = 332.0637778;
93     // Charge-Dipole, assuming charges are measured in electrons, and
94     // dipoles are measured in debyes
95     pre12_ = 69.13373;
96 gezelter 1879 // Dipole-Dipole, assuming dipoles are measured in Debye
97 gezelter 1502 pre22_ = 14.39325;
98     // Charge-Quadrupole, assuming charges are measured in electrons, and
99     // quadrupoles are measured in 10^-26 esu cm^2
100 gezelter 1879 // This unit is also known affectionately as an esu centibarn.
101 gezelter 1502 pre14_ = 69.13373;
102 gezelter 1879 // Dipole-Quadrupole, assuming dipoles are measured in debyes and
103     // quadrupoles in esu centibarns:
104     pre24_ = 14.39325;
105     // Quadrupole-Quadrupole, assuming esu centibarns:
106     pre44_ = 14.39325;
107    
108 gezelter 1502 // conversions for the simulation box dipole moment
109     chargeToC_ = 1.60217733e-19;
110     angstromToM_ = 1.0e-10;
111     debyeToCm_ = 3.33564095198e-30;
112    
113 gezelter 1879 // Default number of points for electrostatic splines
114 gezelter 1502 np_ = 100;
115    
116     // variables to handle different summation methods for long-range
117     // electrostatics:
118 gezelter 1528 summationMethod_ = esm_HARD;
119 gezelter 1502 screeningMethod_ = UNDAMPED;
120     dielectric_ = 1.0;
121    
122 gezelter 1528 // check the summation method:
123     if (simParams_->haveElectrostaticSummationMethod()) {
124     string myMethod = simParams_->getElectrostaticSummationMethod();
125     toUpper(myMethod);
126     map<string, ElectrostaticSummationMethod>::iterator i;
127     i = summationMap_.find(myMethod);
128     if ( i != summationMap_.end() ) {
129     summationMethod_ = (*i).second;
130     } else {
131     // throw error
132     sprintf( painCave.errMsg,
133 gezelter 1536 "Electrostatic::initialize: Unknown electrostaticSummationMethod.\n"
134 gezelter 1528 "\t(Input file specified %s .)\n"
135 gezelter 1616 "\telectrostaticSummationMethod must be one of: \"hard\",\n"
136 gezelter 1879 "\t\"shifted_potential\", \"shifted_force\",\n"
137     "\t\"taylor_shifted\", or \"reaction_field\".\n",
138     myMethod.c_str() );
139 gezelter 1528 painCave.isFatal = 1;
140     simError();
141     }
142     } else {
143     // set ElectrostaticSummationMethod to the cutoffMethod:
144     if (simParams_->haveCutoffMethod()){
145     string myMethod = simParams_->getCutoffMethod();
146     toUpper(myMethod);
147     map<string, ElectrostaticSummationMethod>::iterator i;
148     i = summationMap_.find(myMethod);
149     if ( i != summationMap_.end() ) {
150     summationMethod_ = (*i).second;
151     }
152     }
153     }
154    
155     if (summationMethod_ == esm_REACTION_FIELD) {
156     if (!simParams_->haveDielectric()) {
157     // throw warning
158     sprintf( painCave.errMsg,
159     "SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n"
160     "\tA default value of %f will be used for the dielectric.\n", dielectric_);
161     painCave.isFatal = 0;
162     painCave.severity = OPENMD_INFO;
163     simError();
164     } else {
165     dielectric_ = simParams_->getDielectric();
166     }
167     haveDielectric_ = true;
168     }
169    
170     if (simParams_->haveElectrostaticScreeningMethod()) {
171     string myScreen = simParams_->getElectrostaticScreeningMethod();
172     toUpper(myScreen);
173     map<string, ElectrostaticScreeningMethod>::iterator i;
174     i = screeningMap_.find(myScreen);
175     if ( i != screeningMap_.end()) {
176     screeningMethod_ = (*i).second;
177     } else {
178     sprintf( painCave.errMsg,
179     "SimInfo error: Unknown electrostaticScreeningMethod.\n"
180     "\t(Input file specified %s .)\n"
181     "\telectrostaticScreeningMethod must be one of: \"undamped\"\n"
182     "or \"damped\".\n", myScreen.c_str() );
183     painCave.isFatal = 1;
184     simError();
185     }
186     }
187    
188     // check to make sure a cutoff value has been set:
189     if (!haveCutoffRadius_) {
190     sprintf( painCave.errMsg, "Electrostatic::initialize has no Default "
191     "Cutoff value!\n");
192     painCave.severity = OPENMD_ERROR;
193     painCave.isFatal = 1;
194     simError();
195     }
196    
197 gezelter 1915 if (screeningMethod_ == DAMPED || summationMethod_ == esm_EWALD_FULL) {
198 gezelter 1528 if (!simParams_->haveDampingAlpha()) {
199     // first set a cutoff dependent alpha value
200     // we assume alpha depends linearly with rcut from 0 to 20.5 ang
201     dampingAlpha_ = 0.425 - cutoffRadius_* 0.02;
202 gezelter 1915 if (dampingAlpha_ < 0.0) dampingAlpha_ = 0.0;
203 gezelter 1528 // throw warning
204     sprintf( painCave.errMsg,
205 gezelter 1750 "Electrostatic::initialize: dampingAlpha was not specified in the\n"
206     "\tinput file. A default value of %f (1/ang) will be used for the\n"
207     "\tcutoff of %f (ang).\n",
208 gezelter 1528 dampingAlpha_, cutoffRadius_);
209     painCave.severity = OPENMD_INFO;
210     painCave.isFatal = 0;
211     simError();
212     } else {
213     dampingAlpha_ = simParams_->getDampingAlpha();
214     }
215     haveDampingAlpha_ = true;
216     }
217    
218 gezelter 1915
219 gezelter 1895 Etypes.clear();
220     Etids.clear();
221     FQtypes.clear();
222     FQtids.clear();
223     ElectrostaticMap.clear();
224     Jij.clear();
225     nElectro_ = 0;
226     nFlucq_ = 0;
227    
228     Etids.resize( forceField_->getNAtomType(), -1);
229     FQtids.resize( forceField_->getNAtomType(), -1);
230    
231     set<AtomType*>::iterator at;
232     for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
233     if ((*at)->isElectrostatic()) nElectro_++;
234     if ((*at)->isFluctuatingCharge()) nFlucq_++;
235     }
236 gezelter 1528
237 gezelter 1895 Jij.resize(nFlucq_);
238    
239     for (at = simTypes_.begin(); at != simTypes_.end(); ++at) {
240     if ((*at)->isElectrostatic()) addType(*at);
241 gezelter 1879 }
242    
243     if (summationMethod_ == esm_REACTION_FIELD) {
244     preRF_ = (dielectric_ - 1.0) /
245     ((2.0 * dielectric_ + 1.0) * pow(cutoffRadius_,3) );
246 gezelter 1502 }
247    
248 gezelter 1879 RealType b0c, b1c, b2c, b3c, b4c, b5c;
249     RealType db0c_1, db0c_2, db0c_3, db0c_4, db0c_5;
250     RealType a2, expTerm, invArootPi;
251    
252     RealType r = cutoffRadius_;
253     RealType r2 = r * r;
254     RealType ric = 1.0 / r;
255     RealType ric2 = ric * ric;
256 gezelter 1502
257 gezelter 1879 if (screeningMethod_ == DAMPED) {
258     a2 = dampingAlpha_ * dampingAlpha_;
259     invArootPi = 1.0 / (dampingAlpha_ * sqrt(M_PI));
260     expTerm = exp(-a2 * r2);
261     // values of Smith's B_l functions at the cutoff radius:
262     b0c = erfc(dampingAlpha_ * r) / r;
263     b1c = ( b0c + 2.0*a2 * expTerm * invArootPi) / r2;
264     b2c = (3.0 * b1c + pow(2.0*a2, 2) * expTerm * invArootPi) / r2;
265     b3c = (5.0 * b2c + pow(2.0*a2, 3) * expTerm * invArootPi) / r2;
266     b4c = (7.0 * b3c + pow(2.0*a2, 4) * expTerm * invArootPi) / r2;
267     b5c = (9.0 * b4c + pow(2.0*a2, 5) * expTerm * invArootPi) / r2;
268 gezelter 1900 // Half the Smith self piece:
269     selfMult1_ = - a2 * invArootPi;
270     selfMult2_ = - 2.0 * a2 * a2 * invArootPi / 3.0;
271     selfMult4_ = - 4.0 * a2 * a2 * a2 * invArootPi / 5.0;
272 gezelter 1502 } else {
273 gezelter 1879 a2 = 0.0;
274     b0c = 1.0 / r;
275     b1c = ( b0c) / r2;
276     b2c = (3.0 * b1c) / r2;
277     b3c = (5.0 * b2c) / r2;
278     b4c = (7.0 * b3c) / r2;
279     b5c = (9.0 * b4c) / r2;
280 gezelter 1900 selfMult1_ = 0.0;
281     selfMult2_ = 0.0;
282     selfMult4_ = 0.0;
283 gezelter 1502 }
284 gezelter 1879
285     // higher derivatives of B_0 at the cutoff radius:
286     db0c_1 = -r * b1c;
287     db0c_2 = -b1c + r2 * b2c;
288     db0c_3 = 3.0*r*b2c - r2*r*b3c;
289     db0c_4 = 3.0*b2c - 6.0*r2*b3c + r2*r2*b4c;
290 gezelter 1900 db0c_5 = -15.0*r*b3c + 10.0*r2*r*b4c - r2*r2*r*b5c;
291 gezelter 1879
292 gezelter 1922 if (summationMethod_ != esm_EWALD_FULL) {
293 gezelter 1915 selfMult1_ -= b0c;
294     selfMult2_ += (db0c_2 + 2.0*db0c_1*ric) / 3.0;
295     selfMult4_ -= (db0c_4 + 4.0*db0c_3*ric) / 15.0;
296     }
297    
298 gezelter 1879 // working variables for the splines:
299     RealType ri, ri2;
300     RealType b0, b1, b2, b3, b4, b5;
301     RealType db0_1, db0_2, db0_3, db0_4, db0_5;
302     RealType f, fc, f0;
303     RealType g, gc, g0, g1, g2, g3, g4;
304     RealType h, hc, h1, h2, h3, h4;
305     RealType s, sc, s2, s3, s4;
306     RealType t, tc, t3, t4;
307     RealType u, uc, u4;
308    
309     // working variables for Taylor expansion:
310     RealType rmRc, rmRc2, rmRc3, rmRc4;
311    
312     // Approximate using splines using a maximum of 0.1 Angstroms
313     // between points.
314     int nptest = int((cutoffRadius_ + 2.0) / 0.1);
315     np_ = (np_ > nptest) ? np_ : nptest;
316    
317 gezelter 1613 // Add a 2 angstrom safety window to deal with cutoffGroups that
318     // have charged atoms longer than the cutoffRadius away from each
319 gezelter 1879 // other. Splining is almost certainly the best choice here.
320     // Direct calls to erfc would be preferrable if it is a very fast
321     // implementation.
322 gezelter 1613
323 gezelter 1879 RealType dx = (cutoffRadius_ + 2.0) / RealType(np_);
324    
325     // Storage vectors for the computed functions
326     vector<RealType> rv;
327     vector<RealType> v01v;
328     vector<RealType> v11v;
329     vector<RealType> v21v, v22v;
330     vector<RealType> v31v, v32v;
331     vector<RealType> v41v, v42v, v43v;
332    
333     for (int i = 1; i < np_ + 1; i++) {
334     r = RealType(i) * dx;
335     rv.push_back(r);
336    
337     ri = 1.0 / r;
338     ri2 = ri * ri;
339    
340     r2 = r * r;
341     expTerm = exp(-a2 * r2);
342    
343     // Taylor expansion factors (no need for factorials this way):
344     rmRc = r - cutoffRadius_;
345     rmRc2 = rmRc * rmRc / 2.0;
346     rmRc3 = rmRc2 * rmRc / 3.0;
347     rmRc4 = rmRc3 * rmRc / 4.0;
348    
349     // values of Smith's B_l functions at r:
350     if (screeningMethod_ == DAMPED) {
351     b0 = erfc(dampingAlpha_ * r) * ri;
352     b1 = ( b0 + 2.0*a2 * expTerm * invArootPi) * ri2;
353     b2 = (3.0 * b1 + pow(2.0*a2, 2) * expTerm * invArootPi) * ri2;
354     b3 = (5.0 * b2 + pow(2.0*a2, 3) * expTerm * invArootPi) * ri2;
355     b4 = (7.0 * b3 + pow(2.0*a2, 4) * expTerm * invArootPi) * ri2;
356     b5 = (9.0 * b4 + pow(2.0*a2, 5) * expTerm * invArootPi) * ri2;
357     } else {
358     b0 = ri;
359     b1 = ( b0) * ri2;
360     b2 = (3.0 * b1) * ri2;
361     b3 = (5.0 * b2) * ri2;
362     b4 = (7.0 * b3) * ri2;
363     b5 = (9.0 * b4) * ri2;
364     }
365    
366     // higher derivatives of B_0 at r:
367     db0_1 = -r * b1;
368     db0_2 = -b1 + r2 * b2;
369     db0_3 = 3.0*r*b2 - r2*r*b3;
370     db0_4 = 3.0*b2 - 6.0*r2*b3 + r2*r2*b4;
371     db0_5 = -15.0*r*b3 + 10.0*r2*r*b4 - r2*r2*r*b5;
372    
373     f = b0;
374     fc = b0c;
375     f0 = f - fc - rmRc*db0c_1;
376    
377     g = db0_1;
378     gc = db0c_1;
379     g0 = g - gc;
380     g1 = g0 - rmRc *db0c_2;
381     g2 = g1 - rmRc2*db0c_3;
382     g3 = g2 - rmRc3*db0c_4;
383     g4 = g3 - rmRc4*db0c_5;
384    
385     h = db0_2;
386     hc = db0c_2;
387     h1 = h - hc;
388     h2 = h1 - rmRc *db0c_3;
389     h3 = h2 - rmRc2*db0c_4;
390     h4 = h3 - rmRc3*db0c_5;
391    
392     s = db0_3;
393     sc = db0c_3;
394     s2 = s - sc;
395     s3 = s2 - rmRc *db0c_4;
396     s4 = s3 - rmRc2*db0c_5;
397    
398     t = db0_4;
399     tc = db0c_4;
400     t3 = t - tc;
401     t4 = t3 - rmRc *db0c_5;
402    
403     u = db0_5;
404     uc = db0c_5;
405     u4 = u - uc;
406    
407     // in what follows below, the various v functions are used for
408     // potentials and torques, while the w functions show up in the
409     // forces.
410    
411     switch (summationMethod_) {
412     case esm_SHIFTED_FORCE:
413    
414     v01 = f - fc - rmRc*gc;
415     v11 = g - gc - rmRc*hc;
416     v21 = g*ri - gc*ric - rmRc*(hc - gc*ric)*ric;
417     v22 = h - g*ri - (hc - gc*ric) - rmRc*(sc - (hc - gc*ric)*ric);
418 gezelter 1894 v31 = (h-g*ri)*ri - (hc-gc*ric)*ric - rmRc*(sc-2.0*(hc-gc*ric)*ric)*ric;
419 gezelter 1879 v32 = (s - 3.0*(h-g*ri)*ri) - (sc - 3.0*(hc-gc*ric)*ric)
420     - rmRc*(tc - 3.0*(sc-2.0*(hc-gc*ric)*ric)*ric);
421     v41 = (h - g*ri)*ri2 - (hc - gc*ric)*ric2
422     - rmRc*(sc - 3.0*(hc-gc*ric)*ric)*ric2;
423     v42 = (s-3.0*(h-g*ri)*ri)*ri - (sc-3.0*(hc-gc*ric)*ric)*ric
424     - rmRc*(tc - (4.0*sc - 9.0*(hc - gc*ric)*ric)*ric)*ric;
425    
426     v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri)
427     - (tc - 3.0*(2.0*sc - 5.0*(hc - gc*ric)*ric)*ric)
428     - rmRc*(uc-3.0*(2.0*tc - (7.0*sc - 15.0*(hc - gc*ric)*ric)*ric)*ric);
429    
430     dv01 = g - gc;
431     dv11 = h - hc;
432     dv21 = (h - g*ri)*ri - (hc - gc*ric)*ric;
433     dv22 = (s - (h - g*ri)*ri) - (sc - (hc - gc*ric)*ric);
434     dv31 = (s - 2.0*(h-g*ri)*ri)*ri - (sc - 2.0*(hc-gc*ric)*ric)*ric;
435     dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri)
436     - (tc - 3.0*(sc-2.0*(hc-gc*ric)*ric)*ric);
437     dv41 = (s - 3.0*(h - g*ri)*ri)*ri2 - (sc - 3.0*(hc - gc*ric)*ric)*ric2;
438     dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri
439     - (tc - (4.0*sc - 9.0*(hc-gc*ric)*ric)*ric)*ric;
440     dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri)
441     - (uc - 3.0*(2.0*tc - (7.0*sc - 15.0*(hc - gc*ric)*ric)*ric)*ric);
442    
443     break;
444    
445     case esm_TAYLOR_SHIFTED:
446    
447     v01 = f0;
448     v11 = g1;
449     v21 = g2 * ri;
450     v22 = h2 - v21;
451     v31 = (h3 - g3 * ri) * ri;
452     v32 = s3 - 3.0*v31;
453     v41 = (h4 - g4 * ri) * ri2;
454     v42 = s4 * ri - 3.0*v41;
455     v43 = t4 - 6.0*v42 - 3.0*v41;
456    
457     dv01 = g0;
458     dv11 = h1;
459     dv21 = (h2 - g2*ri)*ri;
460     dv22 = (s2 - (h2 - g2*ri)*ri);
461     dv31 = (s3 - 2.0*(h3-g3*ri)*ri)*ri;
462     dv32 = (t3 - 3.0*(s3-2.0*(h3-g3*ri)*ri)*ri);
463     dv41 = (s4 - 3.0*(h4 - g4*ri)*ri)*ri2;
464     dv42 = (t4 - (4.0*s4 - 9.0*(h4-g4*ri)*ri)*ri)*ri;
465     dv43 = (u4 - 3.0*(2.0*t4 - (7.0*s4 - 15.0*(h4 - g4*ri)*ri)*ri)*ri);
466    
467     break;
468    
469     case esm_SHIFTED_POTENTIAL:
470    
471     v01 = f - fc;
472     v11 = g - gc;
473     v21 = g*ri - gc*ric;
474     v22 = h - g*ri - (hc - gc*ric);
475 gezelter 1894 v31 = (h-g*ri)*ri - (hc-gc*ric)*ric;
476 gezelter 1879 v32 = (s - 3.0*(h-g*ri)*ri) - (sc - 3.0*(hc-gc*ric)*ric);
477     v41 = (h - g*ri)*ri2 - (hc - gc*ric)*ric2;
478     v42 = (s-3.0*(h-g*ri)*ri)*ri - (sc-3.0*(hc-gc*ric)*ric)*ric;
479     v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri)
480     - (tc - 3.0*(2.0*sc - 5.0*(hc - gc*ric)*ric)*ric);
481    
482     dv01 = g;
483     dv11 = h;
484     dv21 = (h - g*ri)*ri;
485     dv22 = (s - (h - g*ri)*ri);
486     dv31 = (s - 2.0*(h-g*ri)*ri)*ri;
487     dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri);
488     dv41 = (s - 3.0*(h - g*ri)*ri)*ri2;
489     dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri;
490     dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri);
491    
492     break;
493    
494     case esm_SWITCHING_FUNCTION:
495     case esm_HARD:
496 gezelter 1915 case esm_EWALD_FULL:
497 gezelter 1879
498     v01 = f;
499     v11 = g;
500     v21 = g*ri;
501     v22 = h - g*ri;
502     v31 = (h-g*ri)*ri;
503     v32 = (s - 3.0*(h-g*ri)*ri);
504     v41 = (h - g*ri)*ri2;
505     v42 = (s-3.0*(h-g*ri)*ri)*ri;
506     v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri);
507    
508     dv01 = g;
509     dv11 = h;
510     dv21 = (h - g*ri)*ri;
511     dv22 = (s - (h - g*ri)*ri);
512     dv31 = (s - 2.0*(h-g*ri)*ri)*ri;
513     dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri);
514     dv41 = (s - 3.0*(h - g*ri)*ri)*ri2;
515     dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri;
516     dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri);
517    
518     break;
519    
520     case esm_REACTION_FIELD:
521    
522     // following DL_POLY's lead for shifting the image charge potential:
523     f = b0 + preRF_ * r2;
524     fc = b0c + preRF_ * cutoffRadius_ * cutoffRadius_;
525    
526     g = db0_1 + preRF_ * 2.0 * r;
527     gc = db0c_1 + preRF_ * 2.0 * cutoffRadius_;
528    
529     h = db0_2 + preRF_ * 2.0;
530     hc = db0c_2 + preRF_ * 2.0;
531    
532     v01 = f - fc;
533     v11 = g - gc;
534     v21 = g*ri - gc*ric;
535     v22 = h - g*ri - (hc - gc*ric);
536 gezelter 1894 v31 = (h-g*ri)*ri - (hc-gc*ric)*ric;
537 gezelter 1879 v32 = (s - 3.0*(h-g*ri)*ri) - (sc - 3.0*(hc-gc*ric)*ric);
538     v41 = (h - g*ri)*ri2 - (hc - gc*ric)*ric2;
539     v42 = (s-3.0*(h-g*ri)*ri)*ri - (sc-3.0*(hc-gc*ric)*ric)*ric;
540     v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri)
541     - (tc - 3.0*(2.0*sc - 5.0*(hc - gc*ric)*ric)*ric);
542    
543     dv01 = g;
544     dv11 = h;
545     dv21 = (h - g*ri)*ri;
546     dv22 = (s - (h - g*ri)*ri);
547     dv31 = (s - 2.0*(h-g*ri)*ri)*ri;
548     dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri);
549     dv41 = (s - 3.0*(h - g*ri)*ri)*ri2;
550     dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri;
551     dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri);
552    
553     break;
554    
555     case esm_EWALD_PME:
556     case esm_EWALD_SPME:
557     default :
558     map<string, ElectrostaticSummationMethod>::iterator i;
559     std::string meth;
560     for (i = summationMap_.begin(); i != summationMap_.end(); ++i) {
561     if ((*i).second == summationMethod_) meth = (*i).first;
562     }
563     sprintf( painCave.errMsg,
564     "Electrostatic::initialize: electrostaticSummationMethod %s \n"
565     "\thas not been implemented yet. Please select one of:\n"
566     "\t\"hard\", \"shifted_potential\", or \"shifted_force\"\n",
567     meth.c_str() );
568     painCave.isFatal = 1;
569     simError();
570     break;
571     }
572    
573     // Add these computed values to the storage vectors for spline creation:
574     v01v.push_back(v01);
575     v11v.push_back(v11);
576     v21v.push_back(v21);
577     v22v.push_back(v22);
578     v31v.push_back(v31);
579     v32v.push_back(v32);
580     v41v.push_back(v41);
581     v42v.push_back(v42);
582     v43v.push_back(v43);
583 gezelter 1502 }
584    
585 gezelter 1879 // construct the spline structures and fill them with the values we've
586     // computed:
587    
588     v01s = new CubicSpline();
589     v01s->addPoints(rv, v01v);
590     v11s = new CubicSpline();
591     v11s->addPoints(rv, v11v);
592     v21s = new CubicSpline();
593     v21s->addPoints(rv, v21v);
594     v22s = new CubicSpline();
595     v22s->addPoints(rv, v22v);
596     v31s = new CubicSpline();
597     v31s->addPoints(rv, v31v);
598     v32s = new CubicSpline();
599     v32s->addPoints(rv, v32v);
600     v41s = new CubicSpline();
601     v41s->addPoints(rv, v41v);
602     v42s = new CubicSpline();
603     v42s->addPoints(rv, v42v);
604     v43s = new CubicSpline();
605     v43s->addPoints(rv, v43v);
606    
607     haveElectroSplines_ = true;
608    
609 gezelter 1502 initialized_ = true;
610     }
611    
612     void Electrostatic::addType(AtomType* atomType){
613 gezelter 1895
614 gezelter 1502 ElectrostaticAtomData electrostaticAtomData;
615     electrostaticAtomData.is_Charge = false;
616     electrostaticAtomData.is_Dipole = false;
617     electrostaticAtomData.is_Quadrupole = false;
618 gezelter 1723 electrostaticAtomData.is_Fluctuating = false;
619 gezelter 1502
620 gezelter 1710 FixedChargeAdapter fca = FixedChargeAdapter(atomType);
621 gezelter 1502
622 gezelter 1710 if (fca.isFixedCharge()) {
623 gezelter 1502 electrostaticAtomData.is_Charge = true;
624 gezelter 1720 electrostaticAtomData.fixedCharge = fca.getCharge();
625 gezelter 1502 }
626    
627 gezelter 1710 MultipoleAdapter ma = MultipoleAdapter(atomType);
628     if (ma.isMultipole()) {
629     if (ma.isDipole()) {
630 gezelter 1502 electrostaticAtomData.is_Dipole = true;
631 gezelter 1879 electrostaticAtomData.dipole = ma.getDipole();
632 gezelter 1502 }
633 gezelter 1710 if (ma.isQuadrupole()) {
634 gezelter 1502 electrostaticAtomData.is_Quadrupole = true;
635 gezelter 1879 electrostaticAtomData.quadrupole = ma.getQuadrupole();
636 gezelter 1502 }
637     }
638    
639 gezelter 1718 FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter(atomType);
640 gezelter 1502
641 gezelter 1718 if (fqa.isFluctuatingCharge()) {
642 gezelter 1720 electrostaticAtomData.is_Fluctuating = true;
643     electrostaticAtomData.electronegativity = fqa.getElectronegativity();
644     electrostaticAtomData.hardness = fqa.getHardness();
645     electrostaticAtomData.slaterN = fqa.getSlaterN();
646     electrostaticAtomData.slaterZeta = fqa.getSlaterZeta();
647 gezelter 1718 }
648    
649 gezelter 1895 int atid = atomType->getIdent();
650     int etid = Etypes.size();
651     int fqtid = FQtypes.size();
652    
653     pair<set<int>::iterator,bool> ret;
654     ret = Etypes.insert( atid );
655 gezelter 1502 if (ret.second == false) {
656     sprintf( painCave.errMsg,
657     "Electrostatic already had a previous entry with ident %d\n",
658 gezelter 1895 atid);
659 gezelter 1502 painCave.severity = OPENMD_INFO;
660     painCave.isFatal = 0;
661     simError();
662     }
663    
664 gezelter 1895 Etids[ atid ] = etid;
665     ElectrostaticMap.push_back(electrostaticAtomData);
666 gezelter 1718
667 gezelter 1895 if (electrostaticAtomData.is_Fluctuating) {
668     ret = FQtypes.insert( atid );
669     if (ret.second == false) {
670     sprintf( painCave.errMsg,
671     "Electrostatic already had a previous fluctuating charge entry with ident %d\n",
672     atid );
673     painCave.severity = OPENMD_INFO;
674     painCave.isFatal = 0;
675     simError();
676     }
677     FQtids[atid] = fqtid;
678     Jij[fqtid].resize(nFlucq_);
679    
680 gezelter 1921 // Now, iterate over all known fluctuating and add to the
681     // coulomb integral map:
682 gezelter 1895
683     std::set<int>::iterator it;
684     for( it = FQtypes.begin(); it != FQtypes.end(); ++it) {
685     int etid2 = Etids[ (*it) ];
686     int fqtid2 = FQtids[ (*it) ];
687     ElectrostaticAtomData eaData2 = ElectrostaticMap[ etid2 ];
688 gezelter 1718 RealType a = electrostaticAtomData.slaterZeta;
689 gezelter 1720 RealType b = eaData2.slaterZeta;
690 gezelter 1718 int m = electrostaticAtomData.slaterN;
691 gezelter 1720 int n = eaData2.slaterN;
692 gezelter 1895
693 gezelter 1718 // Create the spline of the coulombic integral for s-type
694     // Slater orbitals. Add a 2 angstrom safety window to deal
695     // with cutoffGroups that have charged atoms longer than the
696     // cutoffRadius away from each other.
697 gezelter 1895
698 gezelter 1720 RealType rval;
699 gezelter 1718 RealType dr = (cutoffRadius_ + 2.0) / RealType(np_ - 1);
700     vector<RealType> rvals;
701 gezelter 1879 vector<RealType> Jvals;
702     // don't start at i = 0, as rval = 0 is undefined for the
703     // slater overlap integrals.
704 gezelter 1761 for (int i = 1; i < np_+1; i++) {
705 gezelter 1718 rval = RealType(i) * dr;
706     rvals.push_back(rval);
707 gezelter 1879 Jvals.push_back(sSTOCoulInt( a, b, m, n, rval *
708     PhysicalConstants::angstromToBohr ) *
709     PhysicalConstants::hartreeToKcal );
710 gezelter 1718 }
711    
712 gezelter 1879 CubicSpline* J = new CubicSpline();
713     J->addPoints(rvals, Jvals);
714 gezelter 1895 Jij[fqtid][fqtid2] = J;
715     Jij[fqtid2].resize( nFlucq_ );
716     Jij[fqtid2][fqtid] = J;
717     }
718     }
719 gezelter 1502 return;
720     }
721    
722 gezelter 1584 void Electrostatic::setCutoffRadius( RealType rCut ) {
723     cutoffRadius_ = rCut;
724 gezelter 1528 haveCutoffRadius_ = true;
725 gezelter 1502 }
726 gezelter 1584
727 gezelter 1502 void Electrostatic::setElectrostaticSummationMethod( ElectrostaticSummationMethod esm ) {
728     summationMethod_ = esm;
729     }
730     void Electrostatic::setElectrostaticScreeningMethod( ElectrostaticScreeningMethod sm ) {
731     screeningMethod_ = sm;
732     }
733     void Electrostatic::setDampingAlpha( RealType alpha ) {
734     dampingAlpha_ = alpha;
735     haveDampingAlpha_ = true;
736     }
737     void Electrostatic::setReactionFieldDielectric( RealType dielectric ){
738     dielectric_ = dielectric;
739     haveDielectric_ = true;
740     }
741    
742 gezelter 1536 void Electrostatic::calcForce(InteractionData &idat) {
743 gezelter 1502
744 gezelter 1893 if (!initialized_) initialize();
745    
746 gezelter 1895 data1 = ElectrostaticMap[Etids[idat.atid1]];
747     data2 = ElectrostaticMap[Etids[idat.atid2]];
748 gezelter 1502
749 gezelter 1893 U = 0.0; // Potential
750     F.zero(); // Force
751     Ta.zero(); // Torque on site a
752     Tb.zero(); // Torque on site b
753     Ea.zero(); // Electric field at site a
754     Eb.zero(); // Electric field at site b
755     dUdCa = 0.0; // fluctuating charge force at site a
756     dUdCb = 0.0; // fluctuating charge force at site a
757 gezelter 1879
758     // Indirect interactions mediated by the reaction field.
759 gezelter 1893 indirect_Pot = 0.0; // Potential
760     indirect_F.zero(); // Force
761     indirect_Ta.zero(); // Torque on site a
762     indirect_Tb.zero(); // Torque on site b
763 gezelter 1587
764 gezelter 1879 // Excluded potential that is still computed for fluctuating charges
765 gezelter 1893 excluded_Pot= 0.0;
766 gezelter 1879
767    
768 gezelter 1502 // some variables we'll need independent of electrostatic type:
769    
770 gezelter 1879 ri = 1.0 / *(idat.rij);
771 gezelter 1893 rhat = *(idat.d) * ri;
772 gezelter 1879
773 gezelter 1502 // logicals
774    
775 gezelter 1893 a_is_Charge = data1.is_Charge;
776     a_is_Dipole = data1.is_Dipole;
777     a_is_Quadrupole = data1.is_Quadrupole;
778     a_is_Fluctuating = data1.is_Fluctuating;
779 gezelter 1502
780 gezelter 1893 b_is_Charge = data2.is_Charge;
781     b_is_Dipole = data2.is_Dipole;
782     b_is_Quadrupole = data2.is_Quadrupole;
783     b_is_Fluctuating = data2.is_Fluctuating;
784 gezelter 1879
785     // Obtain all of the required radial function values from the
786     // spline structures:
787 gezelter 1502
788 gezelter 1879 // needed for fields (and forces):
789     if (a_is_Charge || b_is_Charge) {
790     v01s->getValueAndDerivativeAt( *(idat.rij), v01, dv01);
791     }
792     if (a_is_Dipole || b_is_Dipole) {
793     v11s->getValueAndDerivativeAt( *(idat.rij), v11, dv11);
794     v11or = ri * v11;
795     }
796     if (a_is_Quadrupole || b_is_Quadrupole || (a_is_Dipole && b_is_Dipole)) {
797     v21s->getValueAndDerivativeAt( *(idat.rij), v21, dv21);
798     v22s->getValueAndDerivativeAt( *(idat.rij), v22, dv22);
799     v22or = ri * v22;
800     }
801 gezelter 1721
802 gezelter 1879 // needed for potentials (and forces and torques):
803     if ((a_is_Dipole && b_is_Quadrupole) ||
804     (b_is_Dipole && a_is_Quadrupole)) {
805     v31s->getValueAndDerivativeAt( *(idat.rij), v31, dv31);
806     v32s->getValueAndDerivativeAt( *(idat.rij), v32, dv32);
807     v31or = v31 * ri;
808     v32or = v32 * ri;
809     }
810     if (a_is_Quadrupole && b_is_Quadrupole) {
811     v41s->getValueAndDerivativeAt( *(idat.rij), v41, dv41);
812     v42s->getValueAndDerivativeAt( *(idat.rij), v42, dv42);
813     v43s->getValueAndDerivativeAt( *(idat.rij), v43, dv43);
814     v42or = v42 * ri;
815     v43or = v43 * ri;
816     }
817    
818     // calculate the single-site contributions (fields, etc).
819    
820     if (a_is_Charge) {
821     C_a = data1.fixedCharge;
822    
823     if (a_is_Fluctuating) {
824     C_a += *(idat.flucQ1);
825 gezelter 1721 }
826    
827 gezelter 1587 if (idat.excluded) {
828 gezelter 1879 *(idat.skippedCharge2) += C_a;
829     } else {
830     // only do the field if we're not excluded:
831     Eb -= C_a * pre11_ * dv01 * rhat;
832 gezelter 1587 }
833     }
834 gezelter 1879
835     if (a_is_Dipole) {
836     D_a = *(idat.dipole1);
837     rdDa = dot(rhat, D_a);
838     rxDa = cross(rhat, D_a);
839     if (!idat.excluded)
840     Eb -= pre12_ * ((dv11-v11or) * rdDa * rhat + v11or * D_a);
841 gezelter 1502 }
842    
843 gezelter 1879 if (a_is_Quadrupole) {
844     Q_a = *(idat.quadrupole1);
845     trQa = Q_a.trace();
846     Qar = Q_a * rhat;
847     rQa = rhat * Q_a;
848     rdQar = dot(rhat, Qar);
849     rxQar = cross(rhat, Qar);
850     if (!idat.excluded)
851     Eb -= pre14_ * (trQa * rhat * dv21 + 2.0 * Qar * v22or
852     + rdQar * rhat * (dv22 - 2.0*v22or));
853     }
854    
855     if (b_is_Charge) {
856     C_b = data2.fixedCharge;
857 gezelter 1502
858 gezelter 1879 if (b_is_Fluctuating)
859     C_b += *(idat.flucQ2);
860    
861 gezelter 1587 if (idat.excluded) {
862 gezelter 1879 *(idat.skippedCharge1) += C_b;
863     } else {
864     // only do the field if we're not excluded:
865     Ea += C_b * pre11_ * dv01 * rhat;
866 gezelter 1587 }
867     }
868 gezelter 1502
869 gezelter 1879 if (b_is_Dipole) {
870     D_b = *(idat.dipole2);
871     rdDb = dot(rhat, D_b);
872     rxDb = cross(rhat, D_b);
873     if (!idat.excluded)
874     Ea += pre12_ * ((dv11-v11or) * rdDb * rhat + v11or * D_b);
875 gezelter 1502 }
876    
877 gezelter 1879 if (b_is_Quadrupole) {
878     Q_b = *(idat.quadrupole2);
879     trQb = Q_b.trace();
880     Qbr = Q_b * rhat;
881     rQb = rhat * Q_b;
882     rdQbr = dot(rhat, Qbr);
883     rxQbr = cross(rhat, Qbr);
884     if (!idat.excluded)
885     Ea += pre14_ * (trQb * rhat * dv21 + 2.0 * Qbr * v22or
886     + rdQbr * rhat * (dv22 - 2.0*v22or));
887 gezelter 1721 }
888 gezelter 1502
889 gezelter 1879 if ((a_is_Fluctuating || b_is_Fluctuating) && idat.excluded) {
890 gezelter 1895 J = Jij[FQtids[idat.atid1]][FQtids[idat.atid2]];
891 gezelter 1879 }
892    
893     if (a_is_Charge) {
894 gezelter 1502
895 gezelter 1879 if (b_is_Charge) {
896     pref = pre11_ * *(idat.electroMult);
897     U += C_a * C_b * pref * v01;
898     F += C_a * C_b * pref * dv01 * rhat;
899    
900     // If this is an excluded pair, there are still indirect
901     // interactions via the reaction field we must worry about:
902 gezelter 1616
903 gezelter 1879 if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) {
904     rfContrib = preRF_ * pref * C_a * C_b * *(idat.r2);
905     indirect_Pot += rfContrib;
906     indirect_F += rfContrib * 2.0 * ri * rhat;
907 gezelter 1502 }
908    
909 gezelter 1879 // Fluctuating charge forces are handled via Coulomb integrals
910     // for excluded pairs (i.e. those connected via bonds) and
911     // with the standard charge-charge interaction otherwise.
912 gezelter 1502
913 gezelter 1879 if (idat.excluded) {
914     if (a_is_Fluctuating || b_is_Fluctuating) {
915     coulInt = J->getValueAt( *(idat.rij) );
916     if (a_is_Fluctuating) dUdCa += coulInt * C_b;
917     if (b_is_Fluctuating) dUdCb += coulInt * C_a;
918     excluded_Pot += C_a * C_b * coulInt;
919     }
920 gezelter 1502 } else {
921 gezelter 1879 if (a_is_Fluctuating) dUdCa += C_b * pref * v01;
922     if (a_is_Fluctuating) dUdCb += C_a * pref * v01;
923 gezelter 1721 }
924 gezelter 1502 }
925    
926 gezelter 1879 if (b_is_Dipole) {
927     pref = pre12_ * *(idat.electroMult);
928     U += C_a * pref * v11 * rdDb;
929     F += C_a * pref * ((dv11 - v11or) * rdDb * rhat + v11or * D_b);
930     Tb += C_a * pref * v11 * rxDb;
931 gezelter 1502
932 gezelter 1879 if (a_is_Fluctuating) dUdCa += pref * v11 * rdDb;
933 gezelter 1502
934 gezelter 1879 // Even if we excluded this pair from direct interactions, we
935     // still have the reaction-field-mediated charge-dipole
936     // interaction:
937 gezelter 1502
938 gezelter 1879 if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) {
939     rfContrib = C_a * pref * preRF_ * 2.0 * *(idat.rij);
940     indirect_Pot += rfContrib * rdDb;
941     indirect_F += rfContrib * D_b / (*idat.rij);
942     indirect_Tb += C_a * pref * preRF_ * rxDb;
943 gezelter 1502 }
944     }
945    
946 gezelter 1879 if (b_is_Quadrupole) {
947     pref = pre14_ * *(idat.electroMult);
948     U += C_a * pref * (v21 * trQb + v22 * rdQbr);
949     F += C_a * pref * (trQb * dv21 * rhat + 2.0 * Qbr * v22or);
950     F += C_a * pref * rdQbr * rhat * (dv22 - 2.0*v22or);
951     Tb += C_a * pref * 2.0 * rxQbr * v22;
952 gezelter 1502
953 gezelter 1879 if (a_is_Fluctuating) dUdCa += pref * (v21 * trQb + v22 * rdQbr);
954 gezelter 1502 }
955     }
956    
957 gezelter 1879 if (a_is_Dipole) {
958 gezelter 1502
959 gezelter 1879 if (b_is_Charge) {
960     pref = pre12_ * *(idat.electroMult);
961 gezelter 1502
962 gezelter 1879 U -= C_b * pref * v11 * rdDa;
963     F -= C_b * pref * ((dv11-v11or) * rdDa * rhat + v11or * D_a);
964     Ta -= C_b * pref * v11 * rxDa;
965 gezelter 1502
966 gezelter 1879 if (b_is_Fluctuating) dUdCb -= pref * v11 * rdDa;
967 gezelter 1587
968 gezelter 1879 // Even if we excluded this pair from direct interactions,
969     // we still have the reaction-field-mediated charge-dipole
970     // interaction:
971     if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) {
972     rfContrib = C_b * pref * preRF_ * 2.0 * *(idat.rij);
973     indirect_Pot -= rfContrib * rdDa;
974     indirect_F -= rfContrib * D_a / (*idat.rij);
975     indirect_Ta -= C_b * pref * preRF_ * rxDa;
976     }
977     }
978 gezelter 1587
979 gezelter 1879 if (b_is_Dipole) {
980     pref = pre22_ * *(idat.electroMult);
981     DadDb = dot(D_a, D_b);
982     DaxDb = cross(D_a, D_b);
983 gezelter 1502
984 gezelter 1879 U -= pref * (DadDb * v21 + rdDa * rdDb * v22);
985     F -= pref * (dv21 * DadDb * rhat + v22or * (rdDb * D_a + rdDa * D_b));
986     F -= pref * (rdDa * rdDb) * (dv22 - 2.0*v22or) * rhat;
987     Ta += pref * ( v21 * DaxDb - v22 * rdDb * rxDa);
988     Tb += pref * (-v21 * DaxDb - v22 * rdDa * rxDb);
989 gezelter 1723
990 gezelter 1879 // Even if we excluded this pair from direct interactions, we
991     // still have the reaction-field-mediated dipole-dipole
992     // interaction:
993     if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) {
994     rfContrib = -pref * preRF_ * 2.0;
995     indirect_Pot += rfContrib * DadDb;
996     indirect_Ta += rfContrib * DaxDb;
997     indirect_Tb -= rfContrib * DaxDb;
998 gezelter 1723 }
999 gezelter 1502 }
1000    
1001 gezelter 1879 if (b_is_Quadrupole) {
1002     pref = pre24_ * *(idat.electroMult);
1003     DadQb = D_a * Q_b;
1004     DadQbr = dot(D_a, Qbr);
1005     DaxQbr = cross(D_a, Qbr);
1006 gezelter 1502
1007 gezelter 1879 U -= pref * ((trQb*rdDa + 2.0*DadQbr)*v31 + rdDa*rdQbr*v32);
1008     F -= pref * (trQb*D_a + 2.0*DadQb) * v31or;
1009     F -= pref * (trQb*rdDa + 2.0*DadQbr) * (dv31-v31or) * rhat;
1010     F -= pref * (D_a*rdQbr + 2.0*rdDa*rQb) * v32or;
1011     F -= pref * (rdDa * rdQbr * rhat * (dv32-3.0*v32or));
1012     Ta += pref * ((-trQb*rxDa + 2.0 * DaxQbr)*v31 - rxDa*rdQbr*v32);
1013     Tb += pref * ((2.0*cross(DadQb, rhat) - 2.0*DaxQbr)*v31
1014     - 2.0*rdDa*rxQbr*v32);
1015     }
1016     }
1017 gezelter 1502
1018 gezelter 1879 if (a_is_Quadrupole) {
1019     if (b_is_Charge) {
1020     pref = pre14_ * *(idat.electroMult);
1021     U += C_b * pref * (v21 * trQa + v22 * rdQar);
1022     F += C_b * pref * (trQa * rhat * dv21 + 2.0 * Qar * v22or);
1023     F += C_b * pref * rdQar * rhat * (dv22 - 2.0*v22or);
1024     Ta += C_b * pref * 2.0 * rxQar * v22;
1025 gezelter 1502
1026 gezelter 1879 if (b_is_Fluctuating) dUdCb += pref * (v21 * trQa + v22 * rdQar);
1027     }
1028     if (b_is_Dipole) {
1029     pref = pre24_ * *(idat.electroMult);
1030     DbdQa = D_b * Q_a;
1031     DbdQar = dot(D_b, Qar);
1032     DbxQar = cross(D_b, Qar);
1033 gezelter 1502
1034 gezelter 1879 U += pref * ((trQa*rdDb + 2.0*DbdQar)*v31 + rdDb*rdQar*v32);
1035     F += pref * (trQa*D_b + 2.0*DbdQa) * v31or;
1036     F += pref * (trQa*rdDb + 2.0*DbdQar) * (dv31-v31or) * rhat;
1037     F += pref * (D_b*rdQar + 2.0*rdDb*rQa) * v32or;
1038     F += pref * (rdDb * rdQar * rhat * (dv32-3.0*v32or));
1039     Ta += pref * ((-2.0*cross(DbdQa, rhat) + 2.0*DbxQar)*v31
1040     + 2.0*rdDb*rxQar*v32);
1041     Tb += pref * ((trQa*rxDb - 2.0 * DbxQar)*v31 + rxDb*rdQar*v32);
1042 gezelter 1502 }
1043 gezelter 1879 if (b_is_Quadrupole) {
1044     pref = pre44_ * *(idat.electroMult); // yes
1045     QaQb = Q_a * Q_b;
1046     trQaQb = QaQb.trace();
1047     rQaQb = rhat * QaQb;
1048     QaQbr = QaQb * rhat;
1049     QaxQb = cross(Q_a, Q_b);
1050     rQaQbr = dot(rQa, Qbr);
1051     rQaxQbr = cross(rQa, Qbr);
1052    
1053     U += pref * (trQa * trQb + 2.0 * trQaQb) * v41;
1054     U += pref * (trQa * rdQbr + trQb * rdQar + 4.0 * rQaQbr) * v42;
1055     U += pref * (rdQar * rdQbr) * v43;
1056 gezelter 1502
1057 gezelter 1879 F += pref * rhat * (trQa * trQb + 2.0 * trQaQb)*dv41;
1058     F += pref*rhat*(trQa*rdQbr + trQb*rdQar + 4.0*rQaQbr)*(dv42-2.0*v42or);
1059     F += pref * rhat * (rdQar * rdQbr)*(dv43 - 4.0*v43or);
1060 gezelter 1502
1061 gezelter 1879 F += pref * 2.0 * (trQb*rQa + trQa*rQb) * v42or;
1062     F += pref * 4.0 * (rQaQb + QaQbr) * v42or;
1063     F += pref * 2.0 * (rQa*rdQbr + rdQar*rQb) * v43or;
1064 gezelter 1502
1065 gezelter 1879 Ta += pref * (- 4.0 * QaxQb * v41);
1066     Ta += pref * (- 2.0 * trQb * cross(rQa, rhat)
1067     + 4.0 * cross(rhat, QaQbr)
1068     - 4.0 * rQaxQbr) * v42;
1069     Ta += pref * 2.0 * cross(rhat,Qar) * rdQbr * v43;
1070 gezelter 1502
1071    
1072 gezelter 1879 Tb += pref * (+ 4.0 * QaxQb * v41);
1073     Tb += pref * (- 2.0 * trQa * cross(rQb, rhat)
1074     - 4.0 * cross(rQaQb, rhat)
1075     + 4.0 * rQaxQbr) * v42;
1076     // Possible replacement for line 2 above:
1077     // + 4.0 * cross(rhat, QbQar)
1078 gezelter 1502
1079 gezelter 1879 Tb += pref * 2.0 * cross(rhat,Qbr) * rdQar * v43;
1080 gezelter 1502
1081     }
1082     }
1083    
1084 gezelter 1879 if (idat.doElectricField) {
1085     *(idat.eField1) += Ea * *(idat.electroMult);
1086     *(idat.eField2) += Eb * *(idat.electroMult);
1087     }
1088 gezelter 1502
1089 gezelter 1879 if (a_is_Fluctuating) *(idat.dVdFQ1) += dUdCa * *(idat.sw);
1090     if (b_is_Fluctuating) *(idat.dVdFQ2) += dUdCb * *(idat.sw);
1091    
1092 gezelter 1587 if (!idat.excluded) {
1093    
1094 gezelter 1879 *(idat.vpair) += U;
1095     (*(idat.pot))[ELECTROSTATIC_FAMILY] += U * *(idat.sw);
1096     *(idat.f1) += F * *(idat.sw);
1097 gezelter 1587
1098 gezelter 1879 if (a_is_Dipole || a_is_Quadrupole)
1099     *(idat.t1) += Ta * *(idat.sw);
1100 gezelter 1587
1101 gezelter 1879 if (b_is_Dipole || b_is_Quadrupole)
1102     *(idat.t2) += Tb * *(idat.sw);
1103    
1104 gezelter 1587 } else {
1105    
1106     // only accumulate the forces and torques resulting from the
1107     // indirect reaction field terms.
1108 gezelter 1616
1109 gezelter 1879 *(idat.vpair) += indirect_Pot;
1110     (*(idat.excludedPot))[ELECTROSTATIC_FAMILY] += excluded_Pot;
1111     (*(idat.pot))[ELECTROSTATIC_FAMILY] += *(idat.sw) * indirect_Pot;
1112     *(idat.f1) += *(idat.sw) * indirect_F;
1113 gezelter 1761
1114 gezelter 1879 if (a_is_Dipole || a_is_Quadrupole)
1115     *(idat.t1) += *(idat.sw) * indirect_Ta;
1116    
1117     if (b_is_Dipole || b_is_Quadrupole)
1118     *(idat.t2) += *(idat.sw) * indirect_Tb;
1119 gezelter 1502 }
1120     return;
1121 gezelter 1879 }
1122 gezelter 1502
1123 gezelter 1545 void Electrostatic::calcSelfCorrection(SelfData &sdat) {
1124 gezelter 1879
1125 gezelter 1502 if (!initialized_) initialize();
1126 gezelter 1586
1127 gezelter 1895 ElectrostaticAtomData data = ElectrostaticMap[Etids[sdat.atid]];
1128 gezelter 1879
1129 gezelter 1502 // logicals
1130     bool i_is_Charge = data.is_Charge;
1131     bool i_is_Dipole = data.is_Dipole;
1132 gezelter 1900 bool i_is_Quadrupole = data.is_Quadrupole;
1133 jmichalk 1734 bool i_is_Fluctuating = data.is_Fluctuating;
1134 gezelter 1879 RealType C_a = data.fixedCharge;
1135 gezelter 1900 RealType self(0.0), preVal, DdD, trQ, trQQ;
1136    
1137     if (i_is_Dipole) {
1138     DdD = data.dipole.lengthSquare();
1139     }
1140    
1141 jmichalk 1734 if (i_is_Fluctuating) {
1142 gezelter 1879 C_a += *(sdat.flucQ);
1143 jmichalk 1734 // dVdFQ is really a force, so this is negative the derivative
1144     *(sdat.dVdFQ) -= *(sdat.flucQ) * data.hardness + data.electronegativity;
1145 gezelter 1761 (*(sdat.excludedPot))[ELECTROSTATIC_FAMILY] += (*sdat.flucQ) *
1146     (*(sdat.flucQ) * data.hardness * 0.5 + data.electronegativity);
1147 jmichalk 1734 }
1148 gezelter 1502
1149 gezelter 1879 switch (summationMethod_) {
1150     case esm_REACTION_FIELD:
1151    
1152     if (i_is_Charge) {
1153     // Self potential [see Wang and Hermans, "Reaction Field
1154     // Molecular Dynamics Simulation with Friedman’s Image Charge
1155     // Method," J. Phys. Chem. 99, 12001-12007 (1995).]
1156     preVal = pre11_ * preRF_ * C_a * C_a;
1157     (*(sdat.pot))[ELECTROSTATIC_FAMILY] -= 0.5 * preVal / cutoffRadius_;
1158     }
1159    
1160 gezelter 1502 if (i_is_Dipole) {
1161 gezelter 1900 (*(sdat.pot))[ELECTROSTATIC_FAMILY] -= pre22_ * preRF_ * DdD;
1162 gezelter 1502 }
1163 gezelter 1879
1164     break;
1165    
1166     case esm_SHIFTED_FORCE:
1167     case esm_SHIFTED_POTENTIAL:
1168 gezelter 1907 case esm_TAYLOR_SHIFTED:
1169 gezelter 1921 case esm_EWALD_FULL:
1170 gezelter 1900 if (i_is_Charge)
1171     self += selfMult1_ * pre11_ * C_a * (C_a + *(sdat.skippedCharge));
1172     if (i_is_Dipole)
1173     self += selfMult2_ * pre22_ * DdD;
1174     if (i_is_Quadrupole) {
1175     trQ = data.quadrupole.trace();
1176     trQQ = (data.quadrupole * data.quadrupole).trace();
1177     self += selfMult4_ * pre44_ * (2.0*trQQ + trQ*trQ);
1178     if (i_is_Charge)
1179     self -= selfMult2_ * pre14_ * 2.0 * C_a * trQ;
1180 gezelter 1502 }
1181 gezelter 1900 (*(sdat.pot))[ELECTROSTATIC_FAMILY] += self;
1182 gezelter 1879 break;
1183     default:
1184     break;
1185 gezelter 1502 }
1186     }
1187 gezelter 1879
1188 gezelter 1545 RealType Electrostatic::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) {
1189 gezelter 1505 // This seems to work moderately well as a default. There's no
1190     // inherent scale for 1/r interactions that we can standardize.
1191     // 12 angstroms seems to be a reasonably good guess for most
1192     // cases.
1193     return 12.0;
1194     }
1195 gezelter 1915
1196    
1197 gezelter 1921 void Electrostatic::ReciprocalSpaceSum(potVec& pot) {
1198 gezelter 1915
1199     RealType kPot = 0.0;
1200     RealType kVir = 0.0;
1201    
1202     const RealType mPoleConverter = 0.20819434; // converts from the
1203     // internal units of
1204     // Debye (for dipoles)
1205     // or Debye-angstroms
1206     // (for quadrupoles) to
1207     // electron angstroms or
1208     // electron-angstroms^2
1209    
1210     const RealType eConverter = 332.0637778; // convert the
1211     // Charge-Charge
1212     // electrostatic
1213     // interactions into kcal /
1214     // mol assuming distances
1215     // are measured in
1216     // angstroms.
1217    
1218     Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat();
1219     Vector3d box = hmat.diagonals();
1220     RealType boxMax = box.max();
1221    
1222 gezelter 1921 //int kMax = int(2.0 * M_PI / (pow(dampingAlpha_,2)*cutoffRadius_ * boxMax) );
1223 gezelter 1922 int kMax = 7;
1224 gezelter 1915 int kSqMax = kMax*kMax + 2;
1225    
1226     int kLimit = kMax+1;
1227     int kLim2 = 2*kMax+1;
1228     int kSqLim = kSqMax;
1229    
1230     vector<RealType> AK(kSqLim+1, 0.0);
1231     RealType xcl = 2.0 * M_PI / box.x();
1232     RealType ycl = 2.0 * M_PI / box.y();
1233     RealType zcl = 2.0 * M_PI / box.z();
1234     RealType rcl = 2.0 * M_PI / boxMax;
1235     RealType rvol = 2.0 * M_PI /(box.x() * box.y() * box.z());
1236    
1237     if(dampingAlpha_ < 1.0e-12) return;
1238    
1239     RealType ralph = -0.25/pow(dampingAlpha_,2);
1240    
1241     // Calculate and store exponential factors
1242    
1243     vector<vector<Vector3d> > eCos;
1244     vector<vector<Vector3d> > eSin;
1245    
1246     int nMax = info_->getNAtoms();
1247    
1248     eCos.resize(kLimit+1);
1249     eSin.resize(kLimit+1);
1250     for (int j = 0; j < kLimit+1; j++) {
1251     eCos[j].resize(nMax);
1252     eSin[j].resize(nMax);
1253     }
1254    
1255     Vector3d t( 2.0 * M_PI );
1256     t.Vdiv(t, box);
1257 gezelter 1922
1258 gezelter 1915
1259     SimInfo::MoleculeIterator mi;
1260     Molecule::AtomIterator ai;
1261     int i;
1262     Vector3d r;
1263     Vector3d tt;
1264     Vector3d w;
1265     Vector3d u;
1266 gezelter 1921 Vector3d a;
1267     Vector3d b;
1268 gezelter 1915
1269     for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1270     mol = info_->nextMolecule(mi)) {
1271     for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1272     atom = mol->nextAtom(ai)) {
1273    
1274     i = atom->getLocalIndex();
1275     r = atom->getPos();
1276     info_->getSnapshotManager()->getCurrentSnapshot()->wrapVector(r);
1277    
1278     tt.Vmul(t, r);
1279 gezelter 1921
1280 gezelter 1915
1281     eCos[1][i] = Vector3d(1.0, 1.0, 1.0);
1282     eSin[1][i] = Vector3d(0.0, 0.0, 0.0);
1283     eCos[2][i] = Vector3d(cos(tt.x()), cos(tt.y()), cos(tt.z()));
1284     eSin[2][i] = Vector3d(sin(tt.x()), sin(tt.y()), sin(tt.z()));
1285 gezelter 1922
1286 gezelter 1921 u = eCos[2][i];
1287     w = eSin[2][i];
1288 gezelter 1915
1289     for(int l = 3; l <= kLimit; l++) {
1290 gezelter 1922 eCos[l][i].x() = eCos[l-1][i].x()*eCos[2][i].x() - eSin[l-1][i].x()*eSin[2][i].x();
1291     eCos[l][i].y() = eCos[l-1][i].y()*eCos[2][i].y() - eSin[l-1][i].y()*eSin[2][i].y();
1292     eCos[l][i].z() = eCos[l-1][i].z()*eCos[2][i].z() - eSin[l-1][i].z()*eSin[2][i].z();
1293    
1294     eSin[l][i].x() = eSin[l-1][i].x()*eCos[2][i].x() + eCos[l-1][i].x()*eSin[2][i].x();
1295     eSin[l][i].y() = eSin[l-1][i].y()*eCos[2][i].y() + eCos[l-1][i].y()*eSin[2][i].z();
1296     eSin[l][i].z() = eSin[l-1][i].z()*eCos[2][i].z() + eCos[l-1][i].z()*eSin[2][i].y();
1297    
1298    
1299     // a.Vmul(eCos[l-1][i], u);
1300     // b.Vmul(eSin[l-1][i], w);
1301     // eCos[l][i] = a - b;
1302     // a.Vmul(eSin[l-1][i], u);
1303     // b.Vmul(eCos[l-1][i], w);
1304     // eSin[l][i] = a + b;
1305    
1306 gezelter 1915 }
1307     }
1308     }
1309    
1310     // Calculate and store AK coefficients:
1311    
1312     RealType eksq = 1.0;
1313     RealType expf = 0.0;
1314     if (ralph < 0.0) expf = exp(ralph*rcl*rcl);
1315     for (i = 1; i <= kSqLim; i++) {
1316     RealType rksq = float(i)*rcl*rcl;
1317     eksq = expf*eksq;
1318     AK[i] = eConverter * eksq/rksq;
1319     }
1320    
1321     /*
1322     * Loop over all k vectors k = 2 pi (ll/Lx, mm/Ly, nn/Lz)
1323     * the values of ll, mm and nn are selected so that the symmetry of
1324     * reciprocal lattice is taken into account i.e. the following
1325     * rules apply.
1326     *
1327     * ll ranges over the values 0 to kMax only.
1328     *
1329     * mm ranges over 0 to kMax when ll=0 and over
1330     * -kMax to kMax otherwise.
1331     * nn ranges over 1 to kMax when ll=mm=0 and over
1332     * -kMax to kMax otherwise.
1333     *
1334     * Hence the result of the summation must be doubled at the end.
1335     */
1336    
1337     std::vector<RealType> clm(nMax, 0.0);
1338     std::vector<RealType> slm(nMax, 0.0);
1339     std::vector<RealType> ckr(nMax, 0.0);
1340     std::vector<RealType> skr(nMax, 0.0);
1341     std::vector<RealType> ckc(nMax, 0.0);
1342     std::vector<RealType> cks(nMax, 0.0);
1343     std::vector<RealType> dkc(nMax, 0.0);
1344     std::vector<RealType> dks(nMax, 0.0);
1345     std::vector<RealType> qkc(nMax, 0.0);
1346     std::vector<RealType> qks(nMax, 0.0);
1347     std::vector<Vector3d> dxk(nMax, V3Zero);
1348     std::vector<Vector3d> qxk(nMax, V3Zero);
1349    
1350     int mMin = kLimit;
1351     int nMin = kLimit + 1;
1352     for (int l = 1; l <= kLimit; l++) {
1353 gezelter 1922 int ll = l - 1;
1354 gezelter 1915 RealType rl = xcl * float(ll);
1355     for (int mmm = mMin; mmm <= kLim2; mmm++) {
1356     int mm = mmm - kLimit;
1357     int m = abs(mm) + 1;
1358     RealType rm = ycl * float(mm);
1359     // Set temporary products of exponential terms
1360     for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1361     mol = info_->nextMolecule(mi)) {
1362     for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1363     atom = mol->nextAtom(ai)) {
1364    
1365     i = atom->getLocalIndex();
1366     if(mm < 0) {
1367     clm[i] = eCos[l][i].x()*eCos[m][i].y()
1368     + eSin[l][i].x()*eSin[m][i].y();
1369 gezelter 1921 slm[i] = eSin[l][i].x()*eCos[m][i].y()
1370 gezelter 1915 - eSin[m][i].y()*eCos[l][i].x();
1371     } else {
1372     clm[i] = eCos[l][i].x()*eCos[m][i].y()
1373     - eSin[l][i].x()*eSin[m][i].y();
1374     slm[i] = eSin[l][i].x()*eCos[m][i].y()
1375 gezelter 1922 + eSin[m][i].y()*eCos[l][i].x();
1376 gezelter 1915 }
1377     }
1378     }
1379     for (int nnn = nMin; nnn <= kLim2; nnn++) {
1380     int nn = nnn - kLimit;
1381     int n = abs(nn) + 1;
1382     RealType rn = zcl * float(nn);
1383     // Test on magnitude of k vector:
1384     int kk=ll*ll + mm*mm + nn*nn;
1385     if(kk <= kSqLim) {
1386     Vector3d kVec = Vector3d(rl, rm, rn);
1387     Mat3x3d k2 = outProduct(kVec, kVec);
1388     // Calculate exp(ikr) terms
1389     for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1390     mol = info_->nextMolecule(mi)) {
1391     for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1392     atom = mol->nextAtom(ai)) {
1393     i = atom->getLocalIndex();
1394    
1395     if (nn < 0) {
1396     ckr[i]=clm[i]*eCos[n][i].z()+slm[i]*eSin[n][i].z();
1397     skr[i]=slm[i]*eCos[n][i].z()-clm[i]*eSin[n][i].z();
1398     } else {
1399     ckr[i]=clm[i]*eCos[n][i].z()-slm[i]*eSin[n][i].z();
1400     skr[i]=slm[i]*eCos[n][i].z()+clm[i]*eSin[n][i].z();
1401     }
1402     }
1403     }
1404    
1405     // Calculate scalar and vector products for each site:
1406    
1407     for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1408     mol = info_->nextMolecule(mi)) {
1409     for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1410     atom = mol->nextAtom(ai)) {
1411 gezelter 1921 i = atom->getLocalIndex();
1412 gezelter 1915 int atid = atom->getAtomType()->getIdent();
1413     ElectrostaticAtomData data = ElectrostaticMap[Etids[atid]];
1414    
1415     if (data.is_Charge) {
1416     RealType C = data.fixedCharge;
1417     if (atom->isFluctuatingCharge()) C += atom->getFlucQPos();
1418     ckc[i] = C * ckr[i];
1419 gezelter 1921 cks[i] = C * skr[i];
1420 gezelter 1915 }
1421    
1422     if (data.is_Dipole) {
1423     Vector3d D = atom->getDipole() * mPoleConverter;
1424     RealType dk = dot(kVec, D);
1425     dxk[i] = cross(kVec, D);
1426     dkc[i] = dk * ckr[i];
1427     dks[i] = dk * skr[i];
1428     }
1429     if (data.is_Quadrupole) {
1430     Mat3x3d Q = atom->getQuadrupole();
1431     Q *= mPoleConverter;
1432     RealType qk = -( Q * k2 ).trace();
1433     qxk[i] = -2.0 * cross(k2, Q);
1434     qkc[i] = qk * ckr[i];
1435     qks[i] = qk * skr[i];
1436     }
1437     }
1438     }
1439 gezelter 1921
1440 gezelter 1915 // calculate vector sums
1441    
1442     RealType ckcs = std::accumulate(ckc.begin(),ckc.end(),0.0);
1443     RealType ckss = std::accumulate(cks.begin(),cks.end(),0.0);
1444     RealType dkcs = std::accumulate(dkc.begin(),dkc.end(),0.0);
1445     RealType dkss = std::accumulate(dks.begin(),dks.end(),0.0);
1446     RealType qkcs = std::accumulate(qkc.begin(),qkc.end(),0.0);
1447     RealType qkss = std::accumulate(qks.begin(),qks.end(),0.0);
1448 gezelter 1921
1449 gezelter 1915
1450     #ifdef IS_MPI
1451     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &ckcs, 1, MPI::REALTYPE,
1452     MPI::SUM);
1453     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &ckss, 1, MPI::REALTYPE,
1454     MPI::SUM);
1455     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &dkcs, 1, MPI::REALTYPE,
1456     MPI::SUM);
1457     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &dkss, 1, MPI::REALTYPE,
1458     MPI::SUM);
1459     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &qkcs, 1, MPI::REALTYPE,
1460     MPI::SUM);
1461     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &qkss, 1, MPI::REALTYPE,
1462     MPI::SUM);
1463     #endif
1464    
1465 gezelter 1921 // Accumulate potential energy and virial contribution:
1466    
1467 gezelter 1915 kPot += 2.0 * rvol * AK[kk]*((ckss+dkcs-qkss)*(ckss+dkcs-qkss)
1468     + (ckcs-dkss-qkcs)*(ckcs-dkss-qkss));
1469 gezelter 1921
1470 gezelter 1915 kVir -= 2.0 * rvol * AK[kk]*(ckcs*ckcs+ckss*ckss
1471     +4.0*(ckss*dkcs-ckcs*dkss)
1472     +3.0*(dkcs*dkcs+dkss*dkss)
1473     -6.0*(ckss*qkss+ckcs*qkcs)
1474     +8.0*(dkss*qkcs-dkcs*qkss)
1475     +5.0*(qkss*qkss+qkcs*qkcs));
1476    
1477     // Calculate force and torque for each site:
1478    
1479     for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1480     mol = info_->nextMolecule(mi)) {
1481     for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1482     atom = mol->nextAtom(ai)) {
1483    
1484     i = atom->getLocalIndex();
1485     int atid = atom->getAtomType()->getIdent();
1486     ElectrostaticAtomData data = ElectrostaticMap[Etids[atid]];
1487    
1488     RealType qfrc = AK[kk]*((cks[i]+dkc[i]-qks[i])*(ckcs-dkss-qkcs)
1489 gezelter 1922 - (ckc[i]-dks[i]-qkc[i])*(ckss+dkcs-qkss));
1490 gezelter 1915 RealType qtrq1 = AK[kk]*(skr[i]*(ckcs-dkss-qkcs)
1491     -ckr[i]*(ckss+dkcs-qkss));
1492     RealType qtrq2 = 2.0*AK[kk]*(ckr[i]*(ckcs-dkss-qkcs)+
1493     skr[i]*(ckss+dkcs-qkss));
1494    
1495     atom->addFrc( 4.0 * rvol * qfrc * kVec );
1496    
1497     if (data.is_Dipole) {
1498     atom->addTrq( 4.0 * rvol * qtrq1 * dxk[i] );
1499     }
1500     if (data.is_Quadrupole) {
1501     atom->addTrq( 4.0 * rvol * qtrq2 * qxk[i] );
1502     }
1503     }
1504     }
1505     }
1506     }
1507 gezelter 1922 nMin = 1;
1508 gezelter 1915 }
1509 gezelter 1922 mMin = 1;
1510 gezelter 1915 }
1511 gezelter 1921 cerr << "kPot = " << kPot << "\n";
1512     pot[ELECTROSTATIC_FAMILY] += kPot;
1513 gezelter 1915 }
1514 gezelter 1502 }

Properties

Name Value
svn:eol-style native