ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/Electrostatic.cpp
Revision: 1938
Committed: Thu Oct 31 15:32:17 2013 UTC (11 years, 8 months ago) by gezelter
File size: 54099 byte(s)
Log Message:
Some MPI include re-ordering to work with the Intel MPI implementation.

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

Properties

Name Value
svn:eol-style native