ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/Electrostatic.cpp
(Generate patch)

Comparing trunk/src/nonbonded/Electrostatic.cpp (file contents):
Revision 1895 by gezelter, Mon Jul 1 21:09:37 2013 UTC vs.
Revision 1929 by gezelter, Mon Aug 19 13:12:00 2013 UTC

# Line 44 | Line 44
44   #include <string.h>
45  
46   #include <cmath>
47 + #include <numeric>
48   #include "nonbonded/Electrostatic.hpp"
49   #include "utils/simError.h"
50   #include "types/NonBondedInteractionType.hpp"
# Line 55 | Line 56
56   #include "utils/PhysicalConstants.hpp"
57   #include "math/erfc.hpp"
58   #include "math/SquareMatrix.hpp"
59 + #include "primitives/Molecule.hpp"
60 + #ifdef IS_MPI
61 + #include <mpi.h>
62 + #endif
63  
64   namespace OpenMD {
65    
# Line 191 | Line 196 | namespace OpenMD {
196        simError();
197      }
198            
199 <    if (screeningMethod_ == DAMPED) {      
199 >    if (screeningMethod_ == DAMPED || summationMethod_ == esm_EWALD_FULL) {
200        if (!simParams_->haveDampingAlpha()) {
201          // first set a cutoff dependent alpha value
202          // we assume alpha depends linearly with rcut from 0 to 20.5 ang
203          dampingAlpha_ = 0.425 - cutoffRadius_* 0.02;
204 <        if (dampingAlpha_ < 0.0) dampingAlpha_ = 0.0;
200 <        
204 >        if (dampingAlpha_ < 0.0) dampingAlpha_ = 0.0;        
205          // throw warning
206          sprintf( painCave.errMsg,
207                   "Electrostatic::initialize: dampingAlpha was not specified in the\n"
# Line 213 | Line 217 | namespace OpenMD {
217        haveDampingAlpha_ = true;
218      }
219  
220 +
221      Etypes.clear();
222      Etids.clear();
223      FQtypes.clear();
# Line 262 | Line 267 | namespace OpenMD {
267        b3c = (5.0 * b2c + pow(2.0*a2, 3) * expTerm * invArootPi) / r2;
268        b4c = (7.0 * b3c + pow(2.0*a2, 4) * expTerm * invArootPi) / r2;
269        b5c = (9.0 * b4c + pow(2.0*a2, 5) * expTerm * invArootPi) / r2;
270 <      selfMult_ = b0c + a2 * invArootPi;
270 >      // Half the Smith self piece:
271 >      selfMult1_ = - a2 * invArootPi;
272 >      selfMult2_ = - 2.0 * a2 * a2 * invArootPi / 3.0;
273 >      selfMult4_ = - 4.0 * a2 * a2 * a2 * invArootPi / 5.0;
274      } else {
275        a2 = 0.0;
276        b0c = 1.0 / r;
# Line 271 | Line 279 | namespace OpenMD {
279        b3c = (5.0 * b2c) / r2;
280        b4c = (7.0 * b3c) / r2;
281        b5c = (9.0 * b4c) / r2;
282 <      selfMult_ = b0c;
282 >      selfMult1_ = 0.0;
283 >      selfMult2_ = 0.0;
284 >      selfMult4_ = 0.0;
285      }
286  
287      // higher derivatives of B_0 at the cutoff radius:
# Line 279 | Line 289 | namespace OpenMD {
289      db0c_2 =     -b1c + r2 * b2c;
290      db0c_3 =          3.0*r*b2c  - r2*r*b3c;
291      db0c_4 =          3.0*b2c  - 6.0*r2*b3c     + r2*r2*b4c;
292 <    db0c_5 =                    -15.0*r*b3c + 10.0*r2*r*b4c - r2*r2*r*b5c;
283 <    
292 >    db0c_5 =                    -15.0*r*b3c + 10.0*r2*r*b4c - r2*r2*r*b5c;  
293  
294 +    if (summationMethod_ != esm_EWALD_FULL) {
295 +      selfMult1_ -= b0c;
296 +      selfMult2_ += (db0c_2 + 2.0*db0c_1*ric) /  3.0;
297 +      selfMult4_ -= (db0c_4 + 4.0*db0c_3*ric) / 15.0;
298 +    }
299 +
300      // working variables for the splines:
301      RealType ri, ri2;
302      RealType b0, b1, b2, b3, b4, b5;
# Line 317 | Line 332 | namespace OpenMD {
332      vector<RealType> v31v, v32v;
333      vector<RealType> v41v, v42v, v43v;
334  
320    /*
321    vector<RealType> dv01v;
322    vector<RealType> dv11v;
323    vector<RealType> dv21v, dv22v;
324    vector<RealType> dv31v, dv32v;
325    vector<RealType> dv41v, dv42v, dv43v;
326    */
327
335      for (int i = 1; i < np_ + 1; i++) {
336        r = RealType(i) * dx;
337        rv.push_back(r);
# Line 488 | Line 495 | namespace OpenMD {
495  
496        case esm_SWITCHING_FUNCTION:
497        case esm_HARD:
498 +      case esm_EWALD_FULL:
499  
500          v01 = f;
501          v11 = g;
# Line 546 | Line 554 | namespace OpenMD {
554  
555          break;
556                  
549      case esm_EWALD_FULL:
557        case esm_EWALD_PME:
558        case esm_EWALD_SPME:
559        default :
# Line 575 | Line 582 | namespace OpenMD {
582        v41v.push_back(v41);
583        v42v.push_back(v42);
584        v43v.push_back(v43);
578      /*
579      dv01v.push_back(dv01);
580      dv11v.push_back(dv11);
581      dv21v.push_back(dv21);
582      dv22v.push_back(dv22);
583      dv31v.push_back(dv31);
584      dv32v.push_back(dv32);      
585      dv41v.push_back(dv41);
586      dv42v.push_back(dv42);
587      dv43v.push_back(dv43);
588      */
585      }
586  
587      // construct the spline structures and fill them with the values we've
# Line 609 | Line 605 | namespace OpenMD {
605      v42s->addPoints(rv, v42v);
606      v43s = new CubicSpline();
607      v43s->addPoints(rv, v43v);
612
613    /*
614    dv01s = new CubicSpline();
615    dv01s->addPoints(rv, dv01v);
616    dv11s = new CubicSpline();
617    dv11s->addPoints(rv, dv11v);
618    dv21s = new CubicSpline();
619    dv21s->addPoints(rv, dv21v);
620    dv22s = new CubicSpline();
621    dv22s->addPoints(rv, dv22v);
622    dv31s = new CubicSpline();
623    dv31s->addPoints(rv, dv31v);
624    dv32s = new CubicSpline();
625    dv32s->addPoints(rv, dv32v);
626    dv41s = new CubicSpline();
627    dv41s->addPoints(rv, dv41v);
628    dv42s = new CubicSpline();
629    dv42s->addPoints(rv, dv42v);
630    dv43s = new CubicSpline();
631    dv43s->addPoints(rv, dv43v);
632    */
608  
609      haveElectroSplines_ = true;
610  
# Line 704 | Line 679 | namespace OpenMD {
679        FQtids[atid] = fqtid;
680        Jij[fqtid].resize(nFlucq_);
681  
682 <      // Now, iterate over all known fluctuating and add to the coulomb integral map:
682 >      // Now, iterate over all known fluctuating and add to the
683 >      // coulomb integral map:
684        
685        std::set<int>::iterator it;
686        for( it = FQtypes.begin(); it != FQtypes.end(); ++it) {    
# Line 912 | Line 888 | namespace OpenMD {
888                          + rdQbr * rhat * (dv22 - 2.0*v22or));
889      }
890      
891 <    if ((a_is_Fluctuating || b_is_Fluctuating) && idat.excluded) {
891 >    if ((a_is_Fluctuating || b_is_Fluctuating)
892 >        && (idat.excluded || idat.sameRegion)) {
893        J = Jij[FQtids[idat.atid1]][FQtids[idat.atid2]];
894      }    
895      
# Line 933 | Line 910 | namespace OpenMD {
910          }
911          
912          // Fluctuating charge forces are handled via Coulomb integrals
913 <        // for excluded pairs (i.e. those connected via bonds) and
914 <        // with the standard charge-charge interaction otherwise.
913 >        // for excluded pairs (i.e. those connected via bonds), atoms
914 >        // within the same region (for metals) and with the standard
915 >        // charge-charge interaction otherwise.
916  
917 <        if (idat.excluded) {          
917 >        if (idat.excluded || idat.sameRegion) {          
918            if (a_is_Fluctuating || b_is_Fluctuating) {
919              coulInt = J->getValueAt( *(idat.rij) );
920              if (a_is_Fluctuating)  dUdCa += coulInt * C_b;
921              if (b_is_Fluctuating)  dUdCb += coulInt * C_a;
922 <            excluded_Pot += C_a * C_b * coulInt;
922 >            if (idat.excluded)
923 >              excluded_Pot += C_a * C_b * coulInt;
924            }          
925          } else {
926            if (a_is_Fluctuating) dUdCa += C_b * pref * v01;
927 <          if (a_is_Fluctuating) dUdCb += C_a * pref * v01;
927 >          if (b_is_Fluctuating) dUdCb += C_a * pref * v01;
928          }
929        }
930  
# Line 1104 | Line 1083 | namespace OpenMD {
1083  
1084          Tb += pref * 2.0 * cross(rhat,Qbr) * rdQar * v43;
1085  
1107        //  cerr << " tsum = " << Ta + Tb - cross(  *(idat.d) , F ) << "\n";
1086        }
1087      }
1088  
# Line 1156 | Line 1134 | namespace OpenMD {
1134      // logicals
1135      bool i_is_Charge = data.is_Charge;
1136      bool i_is_Dipole = data.is_Dipole;
1137 +    bool i_is_Quadrupole = data.is_Quadrupole;
1138      bool i_is_Fluctuating = data.is_Fluctuating;
1139      RealType C_a = data.fixedCharge;  
1140 <    RealType self, preVal, DadDa;
1141 <    
1140 >    RealType self(0.0), preVal, DdD, trQ, trQQ;
1141 >
1142 >    if (i_is_Dipole) {
1143 >      DdD = data.dipole.lengthSquare();
1144 >    }
1145 >        
1146      if (i_is_Fluctuating) {
1147        C_a += *(sdat.flucQ);
1148 <      // dVdFQ is really a force, so this is negative the derivative
1166 <      *(sdat.dVdFQ) -=  *(sdat.flucQ) * data.hardness + data.electronegativity;
1148 >      *(sdat.flucQfrc) -=  *(sdat.flucQ) * data.hardness + data.electronegativity;
1149        (*(sdat.excludedPot))[ELECTROSTATIC_FAMILY] += (*sdat.flucQ) *
1150          (*(sdat.flucQ) * data.hardness * 0.5 + data.electronegativity);
1151      }
# Line 1180 | Line 1162 | namespace OpenMD {
1162        }
1163  
1164        if (i_is_Dipole) {
1165 <        DadDa = data.dipole.lengthSquare();
1184 <        (*(sdat.pot))[ELECTROSTATIC_FAMILY] -= pre22_ * preRF_ * DadDa;
1165 >        (*(sdat.pot))[ELECTROSTATIC_FAMILY] -= pre22_ * preRF_ * DdD;
1166        }
1167        
1168        break;
1169        
1170      case esm_SHIFTED_FORCE:
1171      case esm_SHIFTED_POTENTIAL:
1172 <      if (i_is_Charge) {
1173 <        self = - selfMult_ * C_a * (C_a + *(sdat.skippedCharge)) * pre11_;
1174 <        (*(sdat.pot))[ELECTROSTATIC_FAMILY] += self;
1172 >    case esm_TAYLOR_SHIFTED:
1173 >    case esm_EWALD_FULL:
1174 >      if (i_is_Charge)
1175 >        self += selfMult1_ * pre11_ * C_a * (C_a + *(sdat.skippedCharge));      
1176 >      if (i_is_Dipole)
1177 >        self += selfMult2_ * pre22_ * DdD;      
1178 >      if (i_is_Quadrupole) {
1179 >        trQ = data.quadrupole.trace();
1180 >        trQQ = (data.quadrupole * data.quadrupole).trace();
1181 >        self += selfMult4_ * pre44_ * (2.0*trQQ + trQ*trQ);
1182 >        if (i_is_Charge)
1183 >          self -= selfMult2_ * pre14_ * 2.0 * C_a * trQ;
1184        }
1185 +      (*(sdat.pot))[ELECTROSTATIC_FAMILY] += self;      
1186        break;
1187      default:
1188        break;
# Line 1205 | Line 1196 | namespace OpenMD {
1196      // cases.
1197      return 12.0;
1198    }
1199 +
1200 +
1201 +  void Electrostatic::ReciprocalSpaceSum(RealType& pot) {
1202 +    
1203 +    RealType kPot = 0.0;
1204 +    RealType kVir = 0.0;
1205 +    
1206 +    const RealType mPoleConverter = 0.20819434; // converts from the
1207 +                                                // internal units of
1208 +                                                // Debye (for dipoles)
1209 +                                                // or Debye-angstroms
1210 +                                                // (for quadrupoles) to
1211 +                                                // electron angstroms or
1212 +                                                // electron-angstroms^2
1213 +    
1214 +    const RealType eConverter = 332.0637778; // convert the
1215 +                                             // Charge-Charge
1216 +                                             // electrostatic
1217 +                                             // interactions into kcal /
1218 +                                             // mol assuming distances
1219 +                                             // are measured in
1220 +                                             // angstroms.
1221 +
1222 +    Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat();
1223 +    Vector3d box = hmat.diagonals();
1224 +    RealType boxMax = box.max();
1225 +    
1226 +    //int kMax = int(2.0 * M_PI / (pow(dampingAlpha_,2)*cutoffRadius_ * boxMax) );
1227 +    int kMax = 7;
1228 +    int kSqMax = kMax*kMax + 2;
1229 +    
1230 +    int kLimit = kMax+1;
1231 +    int kLim2 = 2*kMax+1;
1232 +    int kSqLim = kSqMax;
1233 +    
1234 +    vector<RealType> AK(kSqLim+1, 0.0);
1235 +    RealType xcl = 2.0 * M_PI / box.x();
1236 +    RealType ycl = 2.0 * M_PI / box.y();
1237 +    RealType zcl = 2.0 * M_PI / box.z();
1238 +    RealType rcl = 2.0 * M_PI / boxMax;
1239 +    RealType rvol = 2.0 * M_PI /(box.x() * box.y() * box.z());
1240 +    
1241 +    if(dampingAlpha_ < 1.0e-12) return;
1242 +    
1243 +    RealType ralph = -0.25/pow(dampingAlpha_,2);
1244 +    
1245 +    // Calculate and store exponential factors  
1246 +    
1247 +    vector<vector<RealType> > elc;
1248 +    vector<vector<RealType> > emc;
1249 +    vector<vector<RealType> > enc;
1250 +    vector<vector<RealType> > els;
1251 +    vector<vector<RealType> > ems;
1252 +    vector<vector<RealType> > ens;
1253 +
1254 +    
1255 +    int nMax = info_->getNAtoms();
1256 +    
1257 +    elc.resize(kLimit+1);
1258 +    emc.resize(kLimit+1);
1259 +    enc.resize(kLimit+1);
1260 +    els.resize(kLimit+1);
1261 +    ems.resize(kLimit+1);
1262 +    ens.resize(kLimit+1);
1263 +
1264 +    for (int j = 0; j < kLimit+1; j++) {
1265 +      elc[j].resize(nMax);
1266 +      emc[j].resize(nMax);
1267 +      enc[j].resize(nMax);
1268 +      els[j].resize(nMax);
1269 +      ems[j].resize(nMax);
1270 +      ens[j].resize(nMax);
1271 +    }
1272 +    
1273 +    Vector3d t( 2.0 * M_PI );
1274 +    t.Vdiv(t, box);
1275 +
1276 +    
1277 +    SimInfo::MoleculeIterator mi;
1278 +    Molecule::AtomIterator ai;
1279 +    int i;
1280 +    Vector3d r;
1281 +    Vector3d tt;
1282 +    
1283 +    for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1284 +         mol = info_->nextMolecule(mi)) {
1285 +      for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1286 +          atom = mol->nextAtom(ai)) {  
1287 +        
1288 +        i = atom->getLocalIndex();
1289 +        r = atom->getPos();
1290 +        info_->getSnapshotManager()->getCurrentSnapshot()->wrapVector(r);
1291 +        
1292 +        tt.Vmul(t, r);
1293 +
1294 +        elc[1][i] = 1.0;
1295 +        emc[1][i] = 1.0;
1296 +        enc[1][i] = 1.0;
1297 +        els[1][i] = 0.0;
1298 +        ems[1][i] = 0.0;
1299 +        ens[1][i] = 0.0;
1300 +
1301 +        elc[2][i] = cos(tt.x());
1302 +        emc[2][i] = cos(tt.y());
1303 +        enc[2][i] = cos(tt.z());
1304 +        els[2][i] = sin(tt.x());
1305 +        ems[2][i] = sin(tt.y());
1306 +        ens[2][i] = sin(tt.z());
1307 +        
1308 +        for(int l = 3; l <= kLimit; l++) {
1309 +          elc[l][i]=elc[l-1][i]*elc[2][i]-els[l-1][i]*els[2][i];
1310 +          emc[l][i]=emc[l-1][i]*emc[2][i]-ems[l-1][i]*ems[2][i];
1311 +          enc[l][i]=enc[l-1][i]*enc[2][i]-ens[l-1][i]*ens[2][i];
1312 +          els[l][i]=els[l-1][i]*elc[2][i]+elc[l-1][i]*els[2][i];
1313 +          ems[l][i]=ems[l-1][i]*emc[2][i]+emc[l-1][i]*ems[2][i];
1314 +          ens[l][i]=ens[l-1][i]*enc[2][i]+enc[l-1][i]*ens[2][i];
1315 +        }
1316 +      }
1317 +    }
1318 +    
1319 +    // Calculate and store AK coefficients:
1320 +    
1321 +    RealType eksq = 1.0;
1322 +    RealType expf = 0.0;
1323 +    if (ralph < 0.0) expf = exp(ralph*rcl*rcl);
1324 +    for (i = 1; i <= kSqLim; i++) {
1325 +      RealType rksq = float(i)*rcl*rcl;
1326 +      eksq = expf*eksq;
1327 +      AK[i] = eConverter * eksq/rksq;
1328 +    }
1329 +    
1330 +    /*
1331 +     * Loop over all k vectors k = 2 pi (ll/Lx, mm/Ly, nn/Lz)
1332 +     * the values of ll, mm and nn are selected so that the symmetry of
1333 +     * reciprocal lattice is taken into account i.e. the following
1334 +     * rules apply.
1335 +     *
1336 +     * ll ranges over the values 0 to kMax only.
1337 +     *
1338 +     * mm ranges over 0 to kMax when ll=0 and over
1339 +     *            -kMax to kMax otherwise.
1340 +     * nn ranges over 1 to kMax when ll=mm=0 and over
1341 +     *            -kMax to kMax otherwise.
1342 +     *
1343 +     * Hence the result of the summation must be doubled at the end.    
1344 +     */
1345 +    
1346 +    std::vector<RealType> clm(nMax, 0.0);
1347 +    std::vector<RealType> slm(nMax, 0.0);
1348 +    std::vector<RealType> ckr(nMax, 0.0);
1349 +    std::vector<RealType> skr(nMax, 0.0);
1350 +    std::vector<RealType> ckc(nMax, 0.0);
1351 +    std::vector<RealType> cks(nMax, 0.0);
1352 +    std::vector<RealType> dkc(nMax, 0.0);
1353 +    std::vector<RealType> dks(nMax, 0.0);
1354 +    std::vector<RealType> qkc(nMax, 0.0);
1355 +    std::vector<RealType> qks(nMax, 0.0);
1356 +    std::vector<Vector3d> dxk(nMax, V3Zero);
1357 +    std::vector<Vector3d> qxk(nMax, V3Zero);
1358 +    RealType rl, rm, rn;
1359 +    Vector3d kVec;
1360 +    Vector3d Qk;
1361 +    Mat3x3d k2;
1362 +    RealType ckcs, ckss, dkcs, dkss, qkcs, qkss;
1363 +    int atid;
1364 +    ElectrostaticAtomData data;
1365 +    RealType C, dk, qk;
1366 +    Vector3d D;
1367 +    Mat3x3d  Q;
1368 +
1369 +    int mMin = kLimit;
1370 +    int nMin = kLimit + 1;
1371 +    for (int l = 1; l <= kLimit; l++) {
1372 +      int ll = l - 1;
1373 +      rl = xcl * float(ll);
1374 +      for (int mmm = mMin; mmm <= kLim2; mmm++) {
1375 +        int mm = mmm - kLimit;
1376 +        int m = abs(mm) + 1;
1377 +        rm = ycl * float(mm);
1378 +        // Set temporary products of exponential terms
1379 +        for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1380 +             mol = info_->nextMolecule(mi)) {
1381 +          for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1382 +              atom = mol->nextAtom(ai)) {
1383 +            
1384 +            i = atom->getLocalIndex();
1385 +            if(mm < 0) {
1386 +              clm[i]=elc[l][i]*emc[m][i]+els[l][i]*ems[m][i];
1387 +              slm[i]=els[l][i]*emc[m][i]-ems[m][i]*elc[l][i];
1388 +            } else {
1389 +              clm[i]=elc[l][i]*emc[m][i]-els[l][i]*ems[m][i];
1390 +              slm[i]=els[l][i]*emc[m][i]+ems[m][i]*elc[l][i];
1391 +            }
1392 +          }
1393 +        }
1394 +        for (int nnn = nMin; nnn <= kLim2; nnn++) {
1395 +          int nn = nnn - kLimit;          
1396 +          int n = abs(nn) + 1;
1397 +          rn = zcl * float(nn);
1398 +          // Test on magnitude of k vector:
1399 +          int kk=ll*ll + mm*mm + nn*nn;
1400 +          if(kk <= kSqLim) {
1401 +            kVec = Vector3d(rl, rm, rn);
1402 +            k2 = outProduct(kVec, kVec);
1403 +            // Calculate exp(ikr) terms
1404 +            for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1405 +                 mol = info_->nextMolecule(mi)) {
1406 +              for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1407 +                  atom = mol->nextAtom(ai)) {
1408 +                i = atom->getLocalIndex();
1409 +                
1410 +                if (nn < 0) {
1411 +                  ckr[i]=clm[i]*enc[n][i]+slm[i]*ens[n][i];
1412 +                  skr[i]=slm[i]*enc[n][i]-clm[i]*ens[n][i];
1413 +
1414 +                } else {
1415 +                  ckr[i]=clm[i]*enc[n][i]-slm[i]*ens[n][i];
1416 +                  skr[i]=slm[i]*enc[n][i]+clm[i]*ens[n][i];
1417 +                }
1418 +              }
1419 +            }
1420 +            
1421 +            // Calculate scalar and vector products for each site:
1422 +            
1423 +            for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1424 +                 mol = info_->nextMolecule(mi)) {
1425 +              for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1426 +                  atom = mol->nextAtom(ai)) {
1427 +                i = atom->getLocalIndex();
1428 +                int atid = atom->getAtomType()->getIdent();
1429 +                data = ElectrostaticMap[Etids[atid]];
1430 +                              
1431 +                if (data.is_Charge) {
1432 +                  C = data.fixedCharge;
1433 +                  if (atom->isFluctuatingCharge()) C += atom->getFlucQPos();
1434 +                  ckc[i] = C * ckr[i];
1435 +                  cks[i] = C * skr[i];
1436 +                }
1437 +                
1438 +                if (data.is_Dipole) {
1439 +                  D = atom->getDipole() * mPoleConverter;
1440 +                  dk = dot(D, kVec);
1441 +                  dxk[i] = cross(D, kVec);
1442 +                  dkc[i] = dk * ckr[i];
1443 +                  dks[i] = dk * skr[i];
1444 +                }
1445 +                if (data.is_Quadrupole) {
1446 +                  Q = atom->getQuadrupole();
1447 +                  Q *= mPoleConverter;
1448 +                  Qk = Q * kVec;
1449 +                  qk = dot(kVec, Qk);
1450 +                  qxk[i] = cross(kVec, Qk);
1451 +                  qkc[i] = qk * ckr[i];
1452 +                  qks[i] = qk * skr[i];
1453 +                }              
1454 +              }
1455 +            }
1456 +
1457 +            // calculate vector sums
1458 +            
1459 +            ckcs = std::accumulate(ckc.begin(),ckc.end(),0.0);
1460 +            ckss = std::accumulate(cks.begin(),cks.end(),0.0);
1461 +            dkcs = std::accumulate(dkc.begin(),dkc.end(),0.0);
1462 +            dkss = std::accumulate(dks.begin(),dks.end(),0.0);
1463 +            qkcs = std::accumulate(qkc.begin(),qkc.end(),0.0);
1464 +            qkss = std::accumulate(qks.begin(),qks.end(),0.0);
1465 +            
1466 + #ifdef IS_MPI
1467 +            MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &ckcs, 1, MPI::REALTYPE,
1468 +                                      MPI::SUM);
1469 +            MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &ckss, 1, MPI::REALTYPE,
1470 +                                      MPI::SUM);
1471 +            MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &dkcs, 1, MPI::REALTYPE,
1472 +                                      MPI::SUM);
1473 +            MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &dkss, 1, MPI::REALTYPE,
1474 +                                      MPI::SUM);
1475 +            MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &qkcs, 1, MPI::REALTYPE,
1476 +                                      MPI::SUM);
1477 +            MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &qkss, 1, MPI::REALTYPE,
1478 +                                      MPI::SUM);
1479 + #endif        
1480 +            
1481 +            // Accumulate potential energy and virial contribution:
1482 +
1483 +            kPot += 2.0 * rvol * AK[kk]*((ckss+dkcs-qkss)*(ckss+dkcs-qkss)
1484 +                                         + (ckcs-dkss-qkcs)*(ckcs-dkss-qkcs));
1485 +
1486 +            kVir += 2.0 * rvol  * AK[kk]*(ckcs*ckcs+ckss*ckss
1487 +                                          +4.0*(ckss*dkcs-ckcs*dkss)
1488 +                                          +3.0*(dkcs*dkcs+dkss*dkss)
1489 +                                          -6.0*(ckss*qkss+ckcs*qkcs)
1490 +                                          +8.0*(dkss*qkcs-dkcs*qkss)
1491 +                                          +5.0*(qkss*qkss+qkcs*qkcs));
1492 +            
1493 +            // Calculate force and torque for each site:
1494 +            
1495 +            for (Molecule* mol = info_->beginMolecule(mi); mol != NULL;
1496 +                 mol = info_->nextMolecule(mi)) {
1497 +              for(Atom* atom = mol->beginAtom(ai); atom != NULL;
1498 +                  atom = mol->nextAtom(ai)) {
1499 +                
1500 +                i = atom->getLocalIndex();
1501 +                atid = atom->getAtomType()->getIdent();
1502 +                data = ElectrostaticMap[Etids[atid]];
1503 +
1504 +                RealType qfrc = AK[kk]*((cks[i]+dkc[i]-qks[i])*(ckcs-dkss-qkcs)
1505 +                                     - (ckc[i]-dks[i]-qkc[i])*(ckss+dkcs-qkss));
1506 +                RealType qtrq1 = AK[kk]*(skr[i]*(ckcs-dkss-qkcs)
1507 +                                         -ckr[i]*(ckss+dkcs-qkss));
1508 +                RealType qtrq2 = 2.0*AK[kk]*(ckr[i]*(ckcs-dkss-qkcs)
1509 +                                             +skr[i]*(ckss+dkcs-qkss));
1510 +              
1511 +                atom->addFrc( 4.0 * rvol * qfrc * kVec );
1512 +                
1513 +                if (data.is_Dipole) {
1514 +                  atom->addTrq( 4.0 * rvol * qtrq1 * dxk[i] );
1515 +                }
1516 +                if (data.is_Quadrupole) {
1517 +                  atom->addTrq( 4.0 * rvol * qtrq2 * qxk[i] );
1518 +                }
1519 +              }
1520 +            }
1521 +          }
1522 +        }
1523 +        nMin = 1;
1524 +      }
1525 +      mMin = 1;
1526 +    }
1527 +    pot += kPot;  
1528 +  }
1529   }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines