OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ForceMatrixDecomposition.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44#include "parallel/ForceMatrixDecomposition.hpp"
45
46#include "brains/PairList.hpp"
49#include "nonbonded/NonBondedInteraction.hpp"
50
51using namespace std;
52namespace OpenMD {
53
54 ForceMatrixDecomposition::ForceMatrixDecomposition(SimInfo* info,
55 InteractionManager* iMan) :
56 ForceDecomposition(info, iMan) {
57 // Row and colum scans must visit all surrounding cells
58 cellOffsets_.clear();
59 cellOffsets_.push_back(Vector3i(-1, -1, -1));
60 cellOffsets_.push_back(Vector3i(0, -1, -1));
61 cellOffsets_.push_back(Vector3i(1, -1, -1));
62 cellOffsets_.push_back(Vector3i(-1, 0, -1));
63 cellOffsets_.push_back(Vector3i(0, 0, -1));
64 cellOffsets_.push_back(Vector3i(1, 0, -1));
65 cellOffsets_.push_back(Vector3i(-1, 1, -1));
66 cellOffsets_.push_back(Vector3i(0, 1, -1));
67 cellOffsets_.push_back(Vector3i(1, 1, -1));
68 cellOffsets_.push_back(Vector3i(-1, -1, 0));
69 cellOffsets_.push_back(Vector3i(0, -1, 0));
70 cellOffsets_.push_back(Vector3i(1, -1, 0));
71 cellOffsets_.push_back(Vector3i(-1, 0, 0));
72 cellOffsets_.push_back(Vector3i(0, 0, 0));
73 cellOffsets_.push_back(Vector3i(1, 0, 0));
74 cellOffsets_.push_back(Vector3i(-1, 1, 0));
75 cellOffsets_.push_back(Vector3i(0, 1, 0));
76 cellOffsets_.push_back(Vector3i(1, 1, 0));
77 cellOffsets_.push_back(Vector3i(-1, -1, 1));
78 cellOffsets_.push_back(Vector3i(0, -1, 1));
79 cellOffsets_.push_back(Vector3i(1, -1, 1));
80 cellOffsets_.push_back(Vector3i(-1, 0, 1));
81 cellOffsets_.push_back(Vector3i(0, 0, 1));
82 cellOffsets_.push_back(Vector3i(1, 0, 1));
83 cellOffsets_.push_back(Vector3i(-1, 1, 1));
84 cellOffsets_.push_back(Vector3i(0, 1, 1));
85 cellOffsets_.push_back(Vector3i(1, 1, 1));
86 }
87
88 ForceMatrixDecomposition::~ForceMatrixDecomposition() {
89#ifdef IS_MPI
90 delete AtomPlanIntRow;
91 delete AtomPlanRealRow;
92 delete AtomPlanVectorRow;
93 delete AtomPlanMatrixRow;
94 delete AtomPlanPotRow;
95 delete AtomPlanIntColumn;
96 delete AtomPlanRealColumn;
97 delete AtomPlanVectorColumn;
98 delete AtomPlanMatrixColumn;
99 delete AtomPlanPotColumn;
100 delete cgPlanIntRow;
101 delete cgPlanVectorRow;
102 delete cgPlanIntColumn;
103 delete cgPlanVectorColumn;
104#endif
105 }
106
107 /**
108 * distributeInitialData is essentially a copy of the older fortran
109 * SimulationSetup
110 */
112 snap_ = sman_->getCurrentSnapshot();
113 atomStorageLayout_ = sman_->getAtomStorageLayout();
114 ff_ = info_->getForceField();
115 nLocal_ = snap_->getNumberOfAtoms();
116
117 nGroups_ = info_->getNLocalCutoffGroups();
118 // gather the information for atomtype IDs (atids):
119 idents = info_->getIdentArray();
120 regions = info_->getRegions();
121 AtomLocalToGlobal = info_->getGlobalAtomIndices();
122 cgLocalToGlobal = info_->getGlobalGroupIndices();
123 vector<int> globalGroupMembership = info_->getGlobalGroupMembership();
124
125 massFactors = info_->getMassFactors();
126
127 PairList* excludes = info_->getExcludedInteractions();
128 PairList* oneTwo = info_->getOneTwoInteractions();
129 PairList* oneThree = info_->getOneThreeInteractions();
130 PairList* oneFour = info_->getOneFourInteractions();
131
132 if (needVelocities_)
133 snap_->cgData.setStorageLayout(DataStorage::dslPosition |
134 DataStorage::dslVelocity);
135 else
136 snap_->cgData.setStorageLayout(DataStorage::dslPosition);
137
138#ifdef IS_MPI
139
140 MPI_Comm row = rowComm.getComm();
141 MPI_Comm col = colComm.getComm();
142
143 AtomPlanIntRow = new Plan<int>(row, nLocal_);
144 AtomPlanRealRow = new Plan<RealType>(row, nLocal_);
145 AtomPlanVectorRow = new Plan<Vector3d>(row, nLocal_);
146 AtomPlanMatrixRow = new Plan<Mat3x3d>(row, nLocal_);
147 AtomPlanPotRow = new Plan<potVec>(row, nLocal_);
148
149 AtomPlanIntColumn = new Plan<int>(col, nLocal_);
150 AtomPlanRealColumn = new Plan<RealType>(col, nLocal_);
151 AtomPlanVectorColumn = new Plan<Vector3d>(col, nLocal_);
152 AtomPlanMatrixColumn = new Plan<Mat3x3d>(col, nLocal_);
153 AtomPlanPotColumn = new Plan<potVec>(col, nLocal_);
154
155 cgPlanIntRow = new Plan<int>(row, nGroups_);
156 cgPlanVectorRow = new Plan<Vector3d>(row, nGroups_);
157 cgPlanIntColumn = new Plan<int>(col, nGroups_);
158 cgPlanVectorColumn = new Plan<Vector3d>(col, nGroups_);
159
160 nAtomsInRow_ = AtomPlanIntRow->getSize();
161 nAtomsInCol_ = AtomPlanIntColumn->getSize();
162 nGroupsInRow_ = cgPlanIntRow->getSize();
163 nGroupsInCol_ = cgPlanIntColumn->getSize();
164
165 // Modify the data storage objects with the correct layouts and sizes:
166 atomRowData.resize(nAtomsInRow_);
167 atomRowData.setStorageLayout(atomStorageLayout_);
168 atomColData.resize(nAtomsInCol_);
169 atomColData.setStorageLayout(atomStorageLayout_);
170 cgRowData.resize(nGroupsInRow_);
171 cgRowData.setStorageLayout(DataStorage::dslPosition);
172 cgColData.resize(nGroupsInCol_);
173 if (needVelocities_)
174 // we only need column velocities if we need them.
175 cgColData.setStorageLayout(DataStorage::dslPosition |
176 DataStorage::dslVelocity);
177 else
178 cgColData.setStorageLayout(DataStorage::dslPosition);
179
180 identsRow.resize(nAtomsInRow_);
181 identsCol.resize(nAtomsInCol_);
182
183 AtomPlanIntRow->gather(idents, identsRow);
184 AtomPlanIntColumn->gather(idents, identsCol);
185
186 regionsRow.resize(nAtomsInRow_);
187 regionsCol.resize(nAtomsInCol_);
188
189 AtomPlanIntRow->gather(regions, regionsRow);
190 AtomPlanIntColumn->gather(regions, regionsCol);
191
192 // allocate memory for the parallel objects
193 atypesRow.resize(nAtomsInRow_);
194 atypesCol.resize(nAtomsInCol_);
195
196 for (int i = 0; i < nAtomsInRow_; i++)
197 atypesRow[i] = ff_->getAtomType(identsRow[i]);
198 for (int i = 0; i < nAtomsInCol_; i++)
199 atypesCol[i] = ff_->getAtomType(identsCol[i]);
200
201 pot_row.resize(nAtomsInRow_);
202 pot_col.resize(nAtomsInCol_);
203
204 expot_row.resize(nAtomsInRow_);
205 expot_col.resize(nAtomsInCol_);
206
207 selepot_row.resize(nAtomsInRow_);
208 selepot_col.resize(nAtomsInCol_);
209
210 AtomRowToGlobal.resize(nAtomsInRow_);
211 AtomColToGlobal.resize(nAtomsInCol_);
212 AtomPlanIntRow->gather(AtomLocalToGlobal, AtomRowToGlobal);
213 AtomPlanIntColumn->gather(AtomLocalToGlobal, AtomColToGlobal);
214
215 cgRowToGlobal.resize(nGroupsInRow_);
216 cgColToGlobal.resize(nGroupsInCol_);
217 cgPlanIntRow->gather(cgLocalToGlobal, cgRowToGlobal);
218 cgPlanIntColumn->gather(cgLocalToGlobal, cgColToGlobal);
219
220 massFactorsRow.resize(nAtomsInRow_);
221 massFactorsCol.resize(nAtomsInCol_);
222 AtomPlanRealRow->gather(massFactors, massFactorsRow);
223 AtomPlanRealColumn->gather(massFactors, massFactorsCol);
224
225 groupListRow_.clear();
226 groupListRow_.resize(nGroupsInRow_);
227 for (int i = 0; i < nGroupsInRow_; i++) {
228 int gid = cgRowToGlobal[i];
229 for (int j = 0; j < nAtomsInRow_; j++) {
230 int aid = AtomRowToGlobal[j];
231 if (globalGroupMembership[aid] == gid) groupListRow_[i].push_back(j);
232 }
233 }
234
235 groupListCol_.clear();
236 groupListCol_.resize(nGroupsInCol_);
237 for (int i = 0; i < nGroupsInCol_; i++) {
238 int gid = cgColToGlobal[i];
239 for (int j = 0; j < nAtomsInCol_; j++) {
240 int aid = AtomColToGlobal[j];
241 if (globalGroupMembership[aid] == gid) groupListCol_[i].push_back(j);
242 }
243 }
244
245 excludesForAtom.clear();
246 excludesForAtom.resize(nAtomsInRow_);
247 toposForAtom.clear();
248 toposForAtom.resize(nAtomsInRow_);
249 topoDist.clear();
250 topoDist.resize(nAtomsInRow_);
251 for (int i = 0; i < nAtomsInRow_; i++) {
252 int iglob = AtomRowToGlobal[i];
253
254 for (int j = 0; j < nAtomsInCol_; j++) {
255 int jglob = AtomColToGlobal[j];
256
257 if (excludes->hasPair(iglob, jglob)) excludesForAtom[i].push_back(j);
258
259 if (oneTwo->hasPair(iglob, jglob)) {
260 toposForAtom[i].push_back(j);
261 topoDist[i].push_back(1);
262 } else {
263 if (oneThree->hasPair(iglob, jglob)) {
264 toposForAtom[i].push_back(j);
265 topoDist[i].push_back(2);
266 } else {
267 if (oneFour->hasPair(iglob, jglob)) {
268 toposForAtom[i].push_back(j);
269 topoDist[i].push_back(3);
270 }
271 }
272 }
273 }
274 }
275
276#else
277 excludesForAtom.clear();
278 excludesForAtom.resize(nLocal_);
279 toposForAtom.clear();
280 toposForAtom.resize(nLocal_);
281 topoDist.clear();
282 topoDist.resize(nLocal_);
283
284 for (int i = 0; i < nLocal_; i++) {
285 int iglob = AtomLocalToGlobal[i];
286
287 for (int j = 0; j < nLocal_; j++) {
288 int jglob = AtomLocalToGlobal[j];
289
290 if (excludes->hasPair(iglob, jglob)) excludesForAtom[i].push_back(j);
291
292 if (oneTwo->hasPair(iglob, jglob)) {
293 toposForAtom[i].push_back(j);
294 topoDist[i].push_back(1);
295 } else {
296 if (oneThree->hasPair(iglob, jglob)) {
297 toposForAtom[i].push_back(j);
298 topoDist[i].push_back(2);
299 } else {
300 if (oneFour->hasPair(iglob, jglob)) {
301 toposForAtom[i].push_back(j);
302 topoDist[i].push_back(3);
303 }
304 }
305 }
306 }
307 }
308#endif
309
310 // allocate memory for the parallel objects
311 atypesLocal.resize(nLocal_);
312
313 for (int i = 0; i < nLocal_; i++)
314 atypesLocal[i] = ff_->getAtomType(idents[i]);
315
316 groupList_.clear();
317 groupList_.resize(nGroups_);
318 for (int i = 0; i < nGroups_; i++) {
319 int gid = cgLocalToGlobal[i];
320 for (int j = 0; j < nLocal_; j++) {
321 int aid = AtomLocalToGlobal[j];
322 if (globalGroupMembership[aid] == gid) { groupList_[i].push_back(j); }
323 }
324 }
325 }
326
327 int ForceMatrixDecomposition::getTopologicalDistance(int atom1, int atom2) {
328 for (unsigned int j = 0; j < toposForAtom[atom1].size(); j++) {
329 if (toposForAtom[atom1][j] == atom2) return topoDist[atom1][j];
330 }
331 return 0;
332 }
333
334 void ForceMatrixDecomposition::zeroWorkArrays() {
335 pairwisePot = 0.0;
336 selfPot = 0.0;
337 excludedPot = 0.0;
338 excludedSelfPot = 0.0;
339 selectedPot = 0.0;
340 selectedSelfPot = 0.0;
341
342#ifdef IS_MPI
343 if (atomStorageLayout_ & DataStorage::dslForce) {
344 fill(atomRowData.force.begin(), atomRowData.force.end(), V3Zero);
345 fill(atomColData.force.begin(), atomColData.force.end(), V3Zero);
346 }
347
348 if (atomStorageLayout_ & DataStorage::dslTorque) {
349 fill(atomRowData.torque.begin(), atomRowData.torque.end(), V3Zero);
350 fill(atomColData.torque.begin(), atomColData.torque.end(), V3Zero);
351 }
352
353 fill(pot_row.begin(), pot_row.end(),
355
356 fill(pot_col.begin(), pot_col.end(),
358
359 fill(expot_row.begin(), expot_row.end(),
361
362 fill(expot_col.begin(), expot_col.end(),
364
365 fill(selepot_row.begin(), selepot_row.end(),
367
368 fill(selepot_col.begin(), selepot_col.end(),
370
371 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
372 fill(atomRowData.particlePot.begin(), atomRowData.particlePot.end(), 0.0);
373 fill(atomColData.particlePot.begin(), atomColData.particlePot.end(), 0.0);
374 }
375
376 if (atomStorageLayout_ & DataStorage::dslDensity) {
377 fill(atomRowData.density.begin(), atomRowData.density.end(), 0.0);
378 fill(atomColData.density.begin(), atomColData.density.end(), 0.0);
379 }
380
381 if (atomStorageLayout_ & DataStorage::dslFunctional) {
382 fill(atomRowData.functional.begin(), atomRowData.functional.end(), 0.0);
383 fill(atomColData.functional.begin(), atomColData.functional.end(), 0.0);
384 }
385
386 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
387 fill(atomRowData.functionalDerivative.begin(),
388 atomRowData.functionalDerivative.end(), 0.0);
389 fill(atomColData.functionalDerivative.begin(),
390 atomColData.functionalDerivative.end(), 0.0);
391 }
392
393 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
394 fill(atomRowData.skippedCharge.begin(), atomRowData.skippedCharge.end(),
395 0.0);
396 fill(atomColData.skippedCharge.begin(), atomColData.skippedCharge.end(),
397 0.0);
398 }
399
400 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
401 fill(atomRowData.flucQFrc.begin(), atomRowData.flucQFrc.end(), 0.0);
402 fill(atomColData.flucQFrc.begin(), atomColData.flucQFrc.end(), 0.0);
403 }
404
405 if (atomStorageLayout_ & DataStorage::dslElectricField) {
406 fill(atomRowData.electricField.begin(), atomRowData.electricField.end(),
407 V3Zero);
408 fill(atomColData.electricField.begin(), atomColData.electricField.end(),
409 V3Zero);
410 }
411
412 if (atomStorageLayout_ & DataStorage::dslSitePotential) {
413 fill(atomRowData.sitePotential.begin(), atomRowData.sitePotential.end(),
414 0.0);
415 fill(atomColData.sitePotential.begin(), atomColData.sitePotential.end(),
416 0.0);
417 }
418
419#endif
420 // even in parallel, we need to zero out the local arrays:
421
422 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
423 fill(snap_->atomData.particlePot.begin(),
424 snap_->atomData.particlePot.end(), 0.0);
425 }
426
427 if (atomStorageLayout_ & DataStorage::dslDensity) {
428 fill(snap_->atomData.density.begin(), snap_->atomData.density.end(), 0.0);
429 }
430
431 if (atomStorageLayout_ & DataStorage::dslFunctional) {
432 fill(snap_->atomData.functional.begin(), snap_->atomData.functional.end(),
433 0.0);
434 }
435
436 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
437 fill(snap_->atomData.functionalDerivative.begin(),
438 snap_->atomData.functionalDerivative.end(), 0.0);
439 }
440
441 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
442 fill(snap_->atomData.skippedCharge.begin(),
443 snap_->atomData.skippedCharge.end(), 0.0);
444 }
445
446 if (atomStorageLayout_ & DataStorage::dslElectricField) {
447 fill(snap_->atomData.electricField.begin(),
448 snap_->atomData.electricField.end(), V3Zero);
449 }
450 if (atomStorageLayout_ & DataStorage::dslSitePotential) {
451 fill(snap_->atomData.sitePotential.begin(),
452 snap_->atomData.sitePotential.end(), 0.0);
453 }
454 }
455
456 void ForceMatrixDecomposition::distributeData() {
457#ifdef IS_MPI
458
459 snap_ = sman_->getCurrentSnapshot();
460 atomStorageLayout_ = sman_->getAtomStorageLayout();
461
462 bool needsCG = true;
463 if (info_->getNCutoffGroups() != info_->getNAtoms()) needsCG = false;
464
465 // gather up the atomic positions
466 AtomPlanVectorRow->gather(snap_->atomData.position, atomRowData.position);
467 AtomPlanVectorColumn->gather(snap_->atomData.position,
468 atomColData.position);
469
470 // gather up the cutoff group positions
471
472 if (needsCG) {
473 cgPlanVectorRow->gather(snap_->cgData.position, cgRowData.position);
474
475 cgPlanVectorColumn->gather(snap_->cgData.position, cgColData.position);
476 }
477
478 if (needVelocities_) {
479 // gather up the atomic velocities
480 AtomPlanVectorColumn->gather(snap_->atomData.velocity,
481 atomColData.velocity);
482
483 if (needsCG) {
484 cgPlanVectorColumn->gather(snap_->cgData.velocity, cgColData.velocity);
485 }
486 }
487
488 // if needed, gather the atomic rotation matrices
489 if (atomStorageLayout_ & DataStorage::dslAmat) {
490 AtomPlanMatrixRow->gather(snap_->atomData.aMat, atomRowData.aMat);
491 AtomPlanMatrixColumn->gather(snap_->atomData.aMat, atomColData.aMat);
492 }
493
494 // if needed, gather the atomic eletrostatic information
495 if (atomStorageLayout_ & DataStorage::dslDipole) {
496 AtomPlanVectorRow->gather(snap_->atomData.dipole, atomRowData.dipole);
497 AtomPlanVectorColumn->gather(snap_->atomData.dipole, atomColData.dipole);
498 }
499
500 if (atomStorageLayout_ & DataStorage::dslQuadrupole) {
501 AtomPlanMatrixRow->gather(snap_->atomData.quadrupole,
502 atomRowData.quadrupole);
503 AtomPlanMatrixColumn->gather(snap_->atomData.quadrupole,
504 atomColData.quadrupole);
505 }
506
507 // if needed, gather the atomic fluctuating charge values
508 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
509 AtomPlanRealRow->gather(snap_->atomData.flucQPos, atomRowData.flucQPos);
510 AtomPlanRealColumn->gather(snap_->atomData.flucQPos,
511 atomColData.flucQPos);
512 }
513
514#endif
515 }
516
517 /* collects information obtained during the pre-pair loop onto local
518 * data structures.
519 */
520 void ForceMatrixDecomposition::collectIntermediateData() {
521#ifdef IS_MPI
522
523 snap_ = sman_->getCurrentSnapshot();
524 atomStorageLayout_ = sman_->getAtomStorageLayout();
525
526 if (atomStorageLayout_ & DataStorage::dslDensity) {
527 AtomPlanRealRow->scatter(atomRowData.density, snap_->atomData.density);
528
529 int n = snap_->atomData.density.size();
530 vector<RealType> rho_tmp(n, 0.0);
531 AtomPlanRealColumn->scatter(atomColData.density, rho_tmp);
532 for (int i = 0; i < n; i++)
533 snap_->atomData.density[i] += rho_tmp[i];
534 }
535
536 // this isn't necessary if we don't have polarizable atoms, but
537 // we'll leave it here for now.
538 if (atomStorageLayout_ & DataStorage::dslElectricField) {
539 AtomPlanVectorRow->scatter(atomRowData.electricField,
540 snap_->atomData.electricField);
541
542 int n = snap_->atomData.electricField.size();
543 vector<Vector3d> field_tmp(n, V3Zero);
544 AtomPlanVectorColumn->scatter(atomColData.electricField, field_tmp);
545 for (int i = 0; i < n; i++)
546 snap_->atomData.electricField[i] += field_tmp[i];
547 }
548#endif
549 }
550
551 /*
552 * redistributes information obtained during the pre-pair loop out to
553 * row and column-indexed data structures
554 */
555 void ForceMatrixDecomposition::distributeIntermediateData() {
556#ifdef IS_MPI
557 snap_ = sman_->getCurrentSnapshot();
558 atomStorageLayout_ = sman_->getAtomStorageLayout();
559
560 if (atomStorageLayout_ & DataStorage::dslFunctional) {
561 AtomPlanRealRow->gather(snap_->atomData.functional,
562 atomRowData.functional);
563 AtomPlanRealColumn->gather(snap_->atomData.functional,
564 atomColData.functional);
565 }
566
567 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
568 AtomPlanRealRow->gather(snap_->atomData.functionalDerivative,
569 atomRowData.functionalDerivative);
570 AtomPlanRealColumn->gather(snap_->atomData.functionalDerivative,
571 atomColData.functionalDerivative);
572 }
573#endif
574 }
575
576 void ForceMatrixDecomposition::collectData() {
577#ifdef IS_MPI
578 snap_ = sman_->getCurrentSnapshot();
579 atomStorageLayout_ = sman_->getAtomStorageLayout();
580
581 int n = snap_->atomData.force.size();
582 vector<Vector3d> frc_tmp(n, V3Zero);
583
584 AtomPlanVectorRow->scatter(atomRowData.force, frc_tmp);
585 for (int i = 0; i < n; i++) {
586 snap_->atomData.force[i] += frc_tmp[i];
587 frc_tmp[i] = 0.0;
588 }
589
590 AtomPlanVectorColumn->scatter(atomColData.force, frc_tmp);
591 for (int i = 0; i < n; i++) {
592 snap_->atomData.force[i] += frc_tmp[i];
593 }
594
595 if (atomStorageLayout_ & DataStorage::dslTorque) {
596 int nt = snap_->atomData.torque.size();
597 vector<Vector3d> trq_tmp(nt, V3Zero);
598
599 AtomPlanVectorRow->scatter(atomRowData.torque, trq_tmp);
600 for (int i = 0; i < nt; i++) {
601 snap_->atomData.torque[i] += trq_tmp[i];
602 trq_tmp[i] = 0.0;
603 }
604
605 AtomPlanVectorColumn->scatter(atomColData.torque, trq_tmp);
606 for (int i = 0; i < nt; i++)
607 snap_->atomData.torque[i] += trq_tmp[i];
608 }
609
610 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
611 int ns = snap_->atomData.skippedCharge.size();
612 vector<RealType> skch_tmp(ns, 0.0);
613
614 AtomPlanRealRow->scatter(atomRowData.skippedCharge, skch_tmp);
615 for (int i = 0; i < ns; i++) {
616 snap_->atomData.skippedCharge[i] += skch_tmp[i];
617 skch_tmp[i] = 0.0;
618 }
619
620 AtomPlanRealColumn->scatter(atomColData.skippedCharge, skch_tmp);
621 for (int i = 0; i < ns; i++)
622 snap_->atomData.skippedCharge[i] += skch_tmp[i];
623 }
624
625 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
626 int nq = snap_->atomData.flucQFrc.size();
627 vector<RealType> fqfrc_tmp(nq, 0.0);
628
629 AtomPlanRealRow->scatter(atomRowData.flucQFrc, fqfrc_tmp);
630 for (int i = 0; i < nq; i++) {
631 snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
632 fqfrc_tmp[i] = 0.0;
633 }
634
635 AtomPlanRealColumn->scatter(atomColData.flucQFrc, fqfrc_tmp);
636 for (int i = 0; i < nq; i++)
637 snap_->atomData.flucQFrc[i] += fqfrc_tmp[i];
638 }
639
640 if (atomStorageLayout_ & DataStorage::dslElectricField) {
641 int nef = snap_->atomData.electricField.size();
642 vector<Vector3d> efield_tmp(nef, V3Zero);
643
644 AtomPlanVectorRow->scatter(atomRowData.electricField, efield_tmp);
645 for (int i = 0; i < nef; i++) {
646 snap_->atomData.electricField[i] += efield_tmp[i];
647 efield_tmp[i] = 0.0;
648 }
649
650 AtomPlanVectorColumn->scatter(atomColData.electricField, efield_tmp);
651 for (int i = 0; i < nef; i++)
652 snap_->atomData.electricField[i] += efield_tmp[i];
653 }
654
655 if (atomStorageLayout_ & DataStorage::dslSitePotential) {
656 int nsp = snap_->atomData.sitePotential.size();
657 vector<RealType> sp_tmp(nsp, 0.0);
658
659 AtomPlanRealRow->scatter(atomRowData.sitePotential, sp_tmp);
660 for (int i = 0; i < nsp; i++) {
661 snap_->atomData.sitePotential[i] += sp_tmp[i];
662 sp_tmp[i] = 0.0;
663 }
664
665 AtomPlanRealColumn->scatter(atomColData.sitePotential, sp_tmp);
666 for (int i = 0; i < nsp; i++)
667 snap_->atomData.sitePotential[i] += sp_tmp[i];
668 }
669
670 nLocal_ = snap_->getNumberOfAtoms();
671
672 vector<potVec> pot_temp(nLocal_,
674 vector<potVec> expot_temp(nLocal_,
676 vector<potVec> selepot_temp(nLocal_,
678
679 // scatter/gather pot_row into the members of my column
680
681 AtomPlanPotRow->scatter(pot_row, pot_temp);
682 AtomPlanPotRow->scatter(expot_row, expot_temp);
683 AtomPlanPotRow->scatter(selepot_row, selepot_temp);
684
685 for (std::size_t ii = 0; ii < pot_temp.size(); ii++)
686 pairwisePot += pot_temp[ii];
687
688 for (std::size_t ii = 0; ii < expot_temp.size(); ii++)
689 excludedPot += expot_temp[ii];
690
691 for (std::size_t ii = 0; ii < selepot_temp.size(); ii++)
692 selectedPot += selepot_temp[ii];
693
694 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
695 // This is the pairwise contribution to the particle pot. The
696 // embedding contribution is added in each of the low level
697 // non-bonded routines. In single processor, this is done in
698 // unpackInteractionData, not in collectData.
699 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
700 for (int i = 0; i < nLocal_; i++) {
701 // factor of two is because the total potential terms are divided
702 // by 2 in parallel due to row/ column scatter
703 snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
704 }
705 }
706 }
707
708 fill(pot_temp.begin(), pot_temp.end(),
710 fill(expot_temp.begin(), expot_temp.end(),
712 fill(selepot_temp.begin(), selepot_temp.end(),
714
715 AtomPlanPotColumn->scatter(pot_col, pot_temp);
716 AtomPlanPotColumn->scatter(expot_col, expot_temp);
717 AtomPlanPotColumn->scatter(selepot_col, selepot_temp);
718
719 for (std::size_t ii = 0; ii < pot_temp.size(); ii++)
720 pairwisePot += pot_temp[ii];
721
722 for (std::size_t ii = 0; ii < expot_temp.size(); ii++)
723 excludedPot += expot_temp[ii];
724
725 for (std::size_t ii = 0; ii < selepot_temp.size(); ii++)
726 selectedPot += selepot_temp[ii];
727
728 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
729 // This is the pairwise contribution to the particle pot. The
730 // embedding contribution is added in each of the low level
731 // non-bonded routines. In single processor, this is done in
732 // unpackInteractionData, not in collectData.
733 for (int ii = 0; ii < N_INTERACTION_FAMILIES; ii++) {
734 for (int i = 0; i < nLocal_; i++) {
735 // factor of two is because the total potential terms are divided
736 // by 2 in parallel due to row/ column scatter
737 snap_->atomData.particlePot[i] += 2.0 * pot_temp[i](ii);
738 }
739 }
740 }
741
742 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
743 int npp = snap_->atomData.particlePot.size();
744 vector<RealType> ppot_temp(npp, 0.0);
745
746 // This is the direct or embedding contribution to the particle
747 // pot.
748
749 AtomPlanRealRow->scatter(atomRowData.particlePot, ppot_temp);
750 for (int i = 0; i < npp; i++) {
751 snap_->atomData.particlePot[i] += ppot_temp[i];
752 }
753
754 fill(ppot_temp.begin(), ppot_temp.end(), 0.0);
755
756 AtomPlanRealColumn->scatter(atomColData.particlePot, ppot_temp);
757 for (int i = 0; i < npp; i++) {
758 snap_->atomData.particlePot[i] += ppot_temp[i];
759 }
760 }
761
762 MPI_Allreduce(MPI_IN_PLACE, &pairwisePot[0], N_INTERACTION_FAMILIES,
763 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
764
765 MPI_Allreduce(MPI_IN_PLACE, &excludedPot[0], N_INTERACTION_FAMILIES,
766 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
767
768 MPI_Allreduce(MPI_IN_PLACE, &selectedPot[0], N_INTERACTION_FAMILIES,
769 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
770
771 // Here be dragons.
772 MPI_Comm col = colComm.getComm();
773
774 MPI_Allreduce(MPI_IN_PLACE, &snap_->frameData.conductiveHeatFlux[0], 3,
775 MPI_REALTYPE, MPI_SUM, col);
776#endif
777 }
778
779 /**
780 * Collects information obtained during the post-pair (and embedding
781 * functional) loops onto local data structures.
782 */
784#ifdef IS_MPI
785 MPI_Allreduce(MPI_IN_PLACE, &selfPot[0], N_INTERACTION_FAMILIES,
786 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
787
788 MPI_Allreduce(MPI_IN_PLACE, &excludedSelfPot[0], N_INTERACTION_FAMILIES,
789 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
790
791 MPI_Allreduce(MPI_IN_PLACE, &selectedSelfPot[0], N_INTERACTION_FAMILIES,
792 MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
793#endif
794 }
795
796 int& ForceMatrixDecomposition::getNAtomsInRow() {
797#ifdef IS_MPI
798 return nAtomsInRow_;
799#else
800 return nLocal_;
801#endif
802 }
803
804 /**
805 * returns the list of atoms belonging to this group.
806 */
808#ifdef IS_MPI
809 return groupListRow_[cg1];
810#else
811 return groupList_[cg1];
812#endif
813 }
814
815 vector<int>& ForceMatrixDecomposition::getAtomsInGroupColumn(int cg2) {
816#ifdef IS_MPI
817 return groupListCol_[cg2];
818#else
819 return groupList_[cg2];
820#endif
821 }
822
823 Vector3d ForceMatrixDecomposition::getIntergroupVector(int cg1, int cg2) {
824 Vector3d d;
825#ifdef IS_MPI
826 d = cgColData.position[cg2] - cgRowData.position[cg1];
827#else
828 d = snap_->cgData.position[cg2] - snap_->cgData.position[cg1];
829#endif
830
831 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(d); }
832 return d;
833 }
834
835 Vector3d& ForceMatrixDecomposition::getGroupVelocityColumn(int cg2) {
836#ifdef IS_MPI
837 return cgColData.velocity[cg2];
838#else
839 return snap_->cgData.velocity[cg2];
840#endif
841 }
842
843 Vector3d& ForceMatrixDecomposition::getAtomVelocityColumn(int atom2) {
844#ifdef IS_MPI
845 return atomColData.velocity[atom2];
846#else
847 return snap_->atomData.velocity[atom2];
848#endif
849 }
850
851 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorRow(int atom1,
852 int cg1) {
853 Vector3d d;
854
855#ifdef IS_MPI
856 d = cgRowData.position[cg1] - atomRowData.position[atom1];
857#else
858 d = snap_->cgData.position[cg1] - snap_->atomData.position[atom1];
859#endif
860 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(d); }
861 return d;
862 }
863
864 Vector3d ForceMatrixDecomposition::getAtomToGroupVectorColumn(int atom2,
865 int cg2) {
866 Vector3d d;
867
868#ifdef IS_MPI
869 d = cgColData.position[cg2] - atomColData.position[atom2];
870#else
871 d = snap_->cgData.position[cg2] - snap_->atomData.position[atom2];
872#endif
873 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(d); }
874 return d;
875 }
876
877 RealType& ForceMatrixDecomposition::getMassFactorRow(int atom1) {
878#ifdef IS_MPI
879 return massFactorsRow[atom1];
880#else
881 return massFactors[atom1];
882#endif
883 }
884
885 RealType& ForceMatrixDecomposition::getMassFactorColumn(int atom2) {
886#ifdef IS_MPI
887 return massFactorsCol[atom2];
888#else
889 return massFactors[atom2];
890#endif
891 }
892
893 Vector3d ForceMatrixDecomposition::getInteratomicVector(int atom1,
894 int atom2) {
895 Vector3d d;
896
897#ifdef IS_MPI
898 d = atomColData.position[atom2] - atomRowData.position[atom1];
899#else
900 d = snap_->atomData.position[atom2] - snap_->atomData.position[atom1];
901#endif
902 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(d); }
903 return d;
904 }
905
906 vector<int>& ForceMatrixDecomposition::getExcludesForAtom(int atom1) {
907 return excludesForAtom[atom1];
908 }
909
910 /**
911 * We need to exclude some overcounted interactions that result from
912 * the parallel decomposition.
913 */
914#ifdef IS_MPI
915 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2, int, int) {
916 // in MPI, we have to look up the unique IDs for each atom
917 int unique_id_1 = AtomRowToGlobal[atom1];
918 int unique_id_2 = AtomColToGlobal[atom2];
919
920 if (unique_id_1 == unique_id_2) return true;
921
922 // this prevents us from doing the pair on multiple processors
923 if (unique_id_1 < unique_id_2) {
924 if ((unique_id_1 + unique_id_2) % 2 == 0) return true;
925 } else {
926 if ((unique_id_1 + unique_id_2) % 2 == 1) return true;
927 }
928
929 return false;
930 }
931
932#else
933 bool ForceMatrixDecomposition::skipAtomPair(int atom1, int atom2, int cg1,
934 int cg2) {
935 int unique_id_1 = AtomLocalToGlobal[atom1];
936 int unique_id_2 = AtomLocalToGlobal[atom2];
937 int group1 = cgLocalToGlobal[cg1];
938 int group2 = cgLocalToGlobal[cg2];
939
940 if (unique_id_1 == unique_id_2) return true;
941
942 if (group1 == group2) {
943 if (unique_id_1 < unique_id_2) return true;
944 }
945
946 return false;
947 }
948#endif
949
950 /**
951 * We need to handle the interactions for atoms who are involved in
952 * the same rigid body as well as some short range interactions
953 * (bonds, bends, torsions) differently from other interactions.
954 * We'll still visit the pairwise routines, but with a flag that
955 * tells those routines to exclude the pair from direct long range
956 * interactions. Some indirect interactions (notably reaction
957 * field) must still be handled for these pairs.
958 */
959 bool ForceMatrixDecomposition::excludeAtomPair(int atom1, int atom2) {
960 // excludesForAtom was constructed to use row/column indices in the MPI
961 // version, and to use local IDs in the non-MPI version:
962
963 for (vector<int>::iterator i = excludesForAtom[atom1].begin();
964 i != excludesForAtom[atom1].end(); ++i) {
965 if ((*i) == atom2) return true;
966 }
967
968 return false;
969 }
970
971 void ForceMatrixDecomposition::addForceToAtomRow(int atom1, Vector3d fg) {
972#ifdef IS_MPI
973 atomRowData.force[atom1] += fg;
974#else
975 snap_->atomData.force[atom1] += fg;
976#endif
977 }
978
979 void ForceMatrixDecomposition::addForceToAtomColumn(int atom2, Vector3d fg) {
980#ifdef IS_MPI
981 atomColData.force[atom2] += fg;
982#else
983 snap_->atomData.force[atom2] += fg;
984#endif
985 }
986
987 // filling interaction blocks with pointers
988 void ForceMatrixDecomposition::fillInteractionData(InteractionData& idat,
989 int atom1, int atom2,
990 bool newAtom1) {
991 idat.excluded = excludeAtomPair(atom1, atom2);
992
993 if (newAtom1) {
994#ifdef IS_MPI
995 idat.atid1 = identsRow[atom1];
996 idat.atid2 = identsCol[atom2];
997
998 if (regionsRow[atom1] >= 0 && regionsCol[atom2] >= 0) {
999 idat.sameRegion = (regionsRow[atom1] == regionsCol[atom2]);
1000 } else {
1001 idat.sameRegion = false;
1002 }
1003
1004 if (atomStorageLayout_ & DataStorage::dslAmat) {
1005 idat.A1 = atomRowData.aMat[atom1];
1006 idat.A2 = atomColData.aMat[atom2];
1007 }
1008
1009 if (atomStorageLayout_ & DataStorage::dslTorque) {
1010 idat.t1 = atomRowData.torque[atom1];
1011 idat.t2 = atomColData.torque[atom2];
1012 }
1013
1014 if (atomStorageLayout_ & DataStorage::dslDipole) {
1015 idat.D_1 = atomRowData.dipole[atom1];
1016 idat.D_2 = atomColData.dipole[atom2];
1017 }
1018
1019 if (atomStorageLayout_ & DataStorage::dslQuadrupole) {
1020 idat.Q_1 = atomRowData.quadrupole[atom1];
1021 idat.Q_2 = atomColData.quadrupole[atom2];
1022 }
1023
1024 if (atomStorageLayout_ & DataStorage::dslDensity) {
1025 idat.rho1 = atomRowData.density[atom1];
1026 idat.rho2 = atomColData.density[atom2];
1027 }
1028
1029 if (atomStorageLayout_ & DataStorage::dslFunctional) {
1030 idat.frho1 = atomRowData.functional[atom1];
1031 idat.frho2 = atomColData.functional[atom2];
1032 }
1033
1034 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
1035 idat.dfrho1 = atomRowData.functionalDerivative[atom1];
1036 idat.dfrho2 = atomColData.functionalDerivative[atom2];
1037 }
1038
1039 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
1040 idat.particlePot1 = atomRowData.particlePot[atom1];
1041 idat.particlePot2 = atomColData.particlePot[atom2];
1042 }
1043
1044 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
1045 idat.skippedCharge1 = atomRowData.skippedCharge[atom1];
1046 idat.skippedCharge2 = atomColData.skippedCharge[atom2];
1047 }
1048
1049 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
1050 idat.flucQ1 = atomRowData.flucQPos[atom1];
1051 idat.flucQ2 = atomColData.flucQPos[atom2];
1052 }
1053
1054#else
1055
1056 idat.atid1 = idents[atom1];
1057 idat.atid2 = idents[atom2];
1058
1059 if (regions[atom1] >= 0 && regions[atom2] >= 0) {
1060 idat.sameRegion = (regions[atom1] == regions[atom2]);
1061 } else {
1062 idat.sameRegion = false;
1063 }
1064
1065 if (atomStorageLayout_ & DataStorage::dslAmat) {
1066 idat.A1 = snap_->atomData.aMat[atom1];
1067 idat.A2 = snap_->atomData.aMat[atom2];
1068 }
1069
1070 if (atomStorageLayout_ & DataStorage::dslTorque) {
1071 idat.t1 = snap_->atomData.torque[atom1];
1072 idat.t2 = snap_->atomData.torque[atom2];
1073 }
1074
1075 if (atomStorageLayout_ & DataStorage::dslDipole) {
1076 idat.D_1 = snap_->atomData.dipole[atom1];
1077 idat.D_2 = snap_->atomData.dipole[atom2];
1078 }
1079
1080 if (atomStorageLayout_ & DataStorage::dslQuadrupole) {
1081 idat.Q_1 = snap_->atomData.quadrupole[atom1];
1082 idat.Q_2 = snap_->atomData.quadrupole[atom2];
1083 }
1084
1085 if (atomStorageLayout_ & DataStorage::dslDensity) {
1086 idat.rho1 = snap_->atomData.density[atom1];
1087 idat.rho2 = snap_->atomData.density[atom2];
1088 }
1089
1090 if (atomStorageLayout_ & DataStorage::dslFunctional) {
1091 idat.frho1 = snap_->atomData.functional[atom1];
1092 idat.frho2 = snap_->atomData.functional[atom2];
1093 }
1094
1095 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
1096 idat.dfrho1 = snap_->atomData.functionalDerivative[atom1];
1097 idat.dfrho2 = snap_->atomData.functionalDerivative[atom2];
1098 }
1099
1100 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
1101 idat.particlePot1 = snap_->atomData.particlePot[atom1];
1102 idat.particlePot2 = snap_->atomData.particlePot[atom2];
1103 }
1104
1105 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
1106 idat.skippedCharge1 = snap_->atomData.skippedCharge[atom1];
1107 idat.skippedCharge2 = snap_->atomData.skippedCharge[atom2];
1108 }
1109
1110 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
1111 idat.flucQ1 = snap_->atomData.flucQPos[atom1];
1112 idat.flucQ2 = snap_->atomData.flucQPos[atom2];
1113 }
1114#endif
1115
1116 } else {
1117 // atom1 is not new, so don't bother updating properties of that atom:
1118#ifdef IS_MPI
1119 idat.atid2 = identsCol[atom2];
1120
1121 if (regionsRow[atom1] >= 0 && regionsCol[atom2] >= 0) {
1122 idat.sameRegion = (regionsRow[atom1] == regionsCol[atom2]);
1123 } else {
1124 idat.sameRegion = false;
1125 }
1126
1127 if (atomStorageLayout_ & DataStorage::dslAmat) {
1128 idat.A2 = atomColData.aMat[atom2];
1129 }
1130
1131 if (atomStorageLayout_ & DataStorage::dslTorque) {
1132 idat.t2 = atomColData.torque[atom2];
1133 }
1134
1135 if (atomStorageLayout_ & DataStorage::dslDipole) {
1136 idat.D_2 = atomColData.dipole[atom2];
1137 }
1138
1139 if (atomStorageLayout_ & DataStorage::dslQuadrupole) {
1140 idat.Q_2 = atomColData.quadrupole[atom2];
1141 }
1142
1143 if (atomStorageLayout_ & DataStorage::dslDensity) {
1144 idat.rho2 = atomColData.density[atom2];
1145 }
1146
1147 if (atomStorageLayout_ & DataStorage::dslFunctional) {
1148 idat.frho2 = atomColData.functional[atom2];
1149 }
1150
1151 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
1152 idat.dfrho2 = atomColData.functionalDerivative[atom2];
1153 }
1154
1155 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
1156 idat.particlePot2 = atomColData.particlePot[atom2];
1157 }
1158
1159 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
1160 idat.skippedCharge2 = atomColData.skippedCharge[atom2];
1161 }
1162
1163 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
1164 idat.flucQ2 = atomColData.flucQPos[atom2];
1165 }
1166
1167#else
1168 idat.atid2 = idents[atom2];
1169
1170 if (regions[atom1] >= 0 && regions[atom2] >= 0) {
1171 idat.sameRegion = (regions[atom1] == regions[atom2]);
1172 } else {
1173 idat.sameRegion = false;
1174 }
1175
1176 if (atomStorageLayout_ & DataStorage::dslAmat) {
1177 idat.A2 = snap_->atomData.aMat[atom2];
1178 }
1179
1180 if (atomStorageLayout_ & DataStorage::dslTorque) {
1181 idat.t2 = snap_->atomData.torque[atom2];
1182 }
1183
1184 if (atomStorageLayout_ & DataStorage::dslDipole) {
1185 idat.D_2 = snap_->atomData.dipole[atom2];
1186 }
1187
1188 if (atomStorageLayout_ & DataStorage::dslQuadrupole) {
1189 idat.Q_2 = snap_->atomData.quadrupole[atom2];
1190 }
1191
1192 if (atomStorageLayout_ & DataStorage::dslDensity) {
1193 idat.rho2 = snap_->atomData.density[atom2];
1194 }
1195
1196 if (atomStorageLayout_ & DataStorage::dslFunctional) {
1197 idat.frho2 = snap_->atomData.functional[atom2];
1198 }
1199
1200 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
1201 idat.dfrho2 = snap_->atomData.functionalDerivative[atom2];
1202 }
1203
1204 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
1205 idat.particlePot2 = snap_->atomData.particlePot[atom2];
1206 }
1207
1208 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
1209 idat.skippedCharge2 = snap_->atomData.skippedCharge[atom2];
1210 }
1211
1212 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
1213 idat.flucQ2 = snap_->atomData.flucQPos[atom2];
1214 }
1215
1216#endif
1217 }
1218 }
1219
1220 void ForceMatrixDecomposition::unpackInteractionData(InteractionData& idat,
1221 int atom1, int atom2) {
1222#ifdef IS_MPI
1223 pot_row[atom1] += 0.5 * idat.pot;
1224 pot_col[atom2] += 0.5 * idat.pot;
1225 expot_row[atom1] += 0.5 * idat.excludedPot;
1226 expot_col[atom2] += 0.5 * idat.excludedPot;
1227 selepot_row[atom1] += 0.5 * idat.selePot;
1228 selepot_col[atom2] += 0.5 * idat.selePot;
1229
1230 atomRowData.force[atom1] += idat.f1;
1231 atomColData.force[atom2] -= idat.f1;
1232
1233 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
1234 atomRowData.flucQFrc[atom1] -= idat.dVdFQ1;
1235 atomColData.flucQFrc[atom2] -= idat.dVdFQ2;
1236 }
1237
1238 if (atomStorageLayout_ & DataStorage::dslElectricField) {
1239 atomRowData.electricField[atom1] += idat.eField1;
1240 atomColData.electricField[atom2] += idat.eField2;
1241 }
1242
1243 if (atomStorageLayout_ & DataStorage::dslSitePotential) {
1244 atomRowData.sitePotential[atom1] += idat.sPot1;
1245 atomColData.sitePotential[atom2] += idat.sPot2;
1246 }
1247
1248 if (atomStorageLayout_ & DataStorage::dslTorque) {
1249 atomRowData.torque[atom1] = idat.t1;
1250 atomColData.torque[atom2] = idat.t2;
1251 }
1252
1253 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
1254 atomRowData.skippedCharge[atom1] = idat.skippedCharge1;
1255 atomColData.skippedCharge[atom2] = idat.skippedCharge2;
1256 }
1257
1258#else
1259 pairwisePot += idat.pot;
1260 excludedPot += idat.excludedPot;
1261 selectedPot += idat.selePot;
1262
1263 snap_->atomData.force[atom1] += idat.f1;
1264 snap_->atomData.force[atom2] -= idat.f1;
1265
1266 if (idat.doParticlePot) {
1267 // This is the pairwise contribution to the particle pot. The
1268 // self and embedding contribution is added in each of the low
1269 // level non-bonded routines. In parallel, this calculation is
1270 // done in collectData, not in unpackInteractionData.
1271 snap_->atomData.particlePot[atom1] += idat.vpair * idat.sw;
1272 snap_->atomData.particlePot[atom2] += idat.vpair * idat.sw;
1273 }
1274
1275 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
1276 snap_->atomData.flucQFrc[atom1] -= idat.dVdFQ1;
1277 snap_->atomData.flucQFrc[atom2] -= idat.dVdFQ2;
1278 }
1279
1280 if (atomStorageLayout_ & DataStorage::dslElectricField) {
1281 snap_->atomData.electricField[atom1] += idat.eField1;
1282 snap_->atomData.electricField[atom2] += idat.eField2;
1283 }
1284
1285 if (atomStorageLayout_ & DataStorage::dslSitePotential) {
1286 snap_->atomData.sitePotential[atom1] += idat.sPot1;
1287 snap_->atomData.sitePotential[atom2] += idat.sPot2;
1288 }
1289
1290 if (atomStorageLayout_ & DataStorage::dslTorque) {
1291 snap_->atomData.torque[atom1] = idat.t1;
1292 snap_->atomData.torque[atom2] = idat.t2;
1293 }
1294
1295 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
1296 snap_->atomData.skippedCharge[atom1] = idat.skippedCharge1;
1297 snap_->atomData.skippedCharge[atom2] = idat.skippedCharge2;
1298 }
1299
1300#endif
1301 }
1302 void ForceMatrixDecomposition::unpackPrePairData(InteractionData& idat,
1303 int atom1, int atom2) {
1304#ifdef IS_MPI
1305
1306 if (atomStorageLayout_ & DataStorage::dslDensity) {
1307 atomRowData.density[atom1] = idat.rho1;
1308 atomColData.density[atom2] = idat.rho2;
1309 }
1310
1311#else
1312
1313 if (atomStorageLayout_ & DataStorage::dslDensity) {
1314 snap_->atomData.density[atom1] = idat.rho1;
1315 snap_->atomData.density[atom2] = idat.rho2;
1316 }
1317
1318#endif
1319 }
1320
1321 /*
1322 * buildNeighborList
1323 *
1324 * Constructs the Verlet neighbor list for a force-matrix
1325 * decomposition. In this case, each processor is responsible for
1326 * row-site interactions with column-sites.
1327 *
1328 * neighborList is returned as a packed array of neighboring
1329 * column-ordered CutoffGroups. The starting position in
1330 * neighborList for each row-ordered CutoffGroup is given by the
1331 * returned vector point.
1332 */
1333 void ForceMatrixDecomposition::buildNeighborList(
1334 vector<int>& neighborList, vector<int>& point,
1335 vector<Vector3d>& savedPositions) {
1336 neighborList.clear();
1337 point.clear();
1338 int len = 0;
1339
1340 bool doAllPairs = false;
1341
1342 Snapshot* snap_ = sman_->getCurrentSnapshot();
1343 Mat3x3d box;
1344 Mat3x3d invBox;
1345
1346 Vector3d rs, scaled, dr;
1347 Vector3i whichCell;
1348 int cellIndex;
1349
1350#ifdef IS_MPI
1351 cellListRow_.clear();
1352 cellListCol_.clear();
1353 point.resize(nGroupsInRow_ + 1);
1354#else
1355 cellList_.clear();
1356 point.resize(nGroups_ + 1);
1357#endif
1358
1359 if (!usePeriodicBoundaryConditions_) {
1360 box = snap_->getBoundingBox();
1361 invBox = snap_->getInvBoundingBox();
1362 } else {
1363 box = snap_->getHmat();
1364 invBox = snap_->getInvHmat();
1365 }
1366
1367 Vector3d A = box.getColumn(0);
1368 Vector3d B = box.getColumn(1);
1369 Vector3d C = box.getColumn(2);
1370
1371 // Required for triclinic cells
1372 Vector3d AxB = cross(A, B);
1373 Vector3d BxC = cross(B, C);
1374 Vector3d CxA = cross(C, A);
1375
1376 // unit vectors perpendicular to the faces of the triclinic cell:
1377 AxB.normalize();
1378 BxC.normalize();
1379 CxA.normalize();
1380
1381 // A set of perpendicular lengths in triclinic cells:
1382 RealType Wa = abs(dot(A, BxC));
1383 RealType Wb = abs(dot(B, CxA));
1384 RealType Wc = abs(dot(C, AxB));
1385
1386 nCells_.x() = int(Wa / rList_);
1387 nCells_.y() = int(Wb / rList_);
1388 nCells_.z() = int(Wc / rList_);
1389
1390 // handle small boxes where the cell offsets can end up repeating cells
1391 if (nCells_.x() < 3) doAllPairs = true;
1392 if (nCells_.y() < 3) doAllPairs = true;
1393 if (nCells_.z() < 3) doAllPairs = true;
1394
1395 int nCtot = nCells_.x() * nCells_.y() * nCells_.z();
1396
1397#ifdef IS_MPI
1398 cellListRow_.resize(nCtot);
1399 cellListCol_.resize(nCtot);
1400#else
1401 cellList_.resize(nCtot);
1402#endif
1403
1404 if (!doAllPairs) {
1405#ifdef IS_MPI
1406
1407 for (int i = 0; i < nGroupsInRow_; i++) {
1408 rs = cgRowData.position[i];
1409
1410 // scaled positions relative to the box vectors
1411 scaled = invBox * rs;
1412
1413 // wrap the vector back into the unit box by subtracting integer box
1414 // numbers
1415 for (int j = 0; j < 3; j++) {
1416 scaled[j] -= roundMe(scaled[j]);
1417 scaled[j] += 0.5;
1418 // Handle the special case when an object is exactly on the
1419 // boundary (a scaled coordinate of 1.0 is the same as
1420 // scaled coordinate of 0.0)
1421 if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1422 }
1423
1424 // find xyz-indices of cell that cutoffGroup is in.
1425 whichCell.x() = int(nCells_.x() * scaled.x());
1426 whichCell.y() = int(nCells_.y() * scaled.y());
1427 whichCell.z() = int(nCells_.z() * scaled.z());
1428
1429 // find single index of this cell:
1430 cellIndex = Vlinear(whichCell, nCells_);
1431
1432 // add this cutoff group to the list of groups in this cell;
1433 cellListRow_[cellIndex].push_back(i);
1434 }
1435 for (int i = 0; i < nGroupsInCol_; i++) {
1436 rs = cgColData.position[i];
1437
1438 // scaled positions relative to the box vectors
1439 scaled = invBox * rs;
1440
1441 // wrap the vector back into the unit box by subtracting integer box
1442 // numbers
1443 for (int j = 0; j < 3; j++) {
1444 scaled[j] -= roundMe(scaled[j]);
1445 scaled[j] += 0.5;
1446 // Handle the special case when an object is exactly on the
1447 // boundary (a scaled coordinate of 1.0 is the same as
1448 // scaled coordinate of 0.0)
1449 if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1450 }
1451
1452 // find xyz-indices of cell that cutoffGroup is in.
1453 whichCell.x() = int(nCells_.x() * scaled.x());
1454 whichCell.y() = int(nCells_.y() * scaled.y());
1455 whichCell.z() = int(nCells_.z() * scaled.z());
1456
1457 // find single index of this cell:
1458 cellIndex = Vlinear(whichCell, nCells_);
1459
1460 // add this cutoff group to the list of groups in this cell;
1461 cellListCol_[cellIndex].push_back(i);
1462 }
1463
1464#else
1465 for (int i = 0; i < nGroups_; i++) {
1466 rs = snap_->cgData.position[i];
1467
1468 // scaled positions relative to the box vectors
1469 scaled = invBox * rs;
1470
1471 // wrap the vector back into the unit box by subtracting integer box
1472 // numbers
1473 for (int j = 0; j < 3; j++) {
1474 scaled[j] -= roundMe(scaled[j]);
1475 scaled[j] += 0.5;
1476 // Handle the special case when an object is exactly on the
1477 // boundary (a scaled coordinate of 1.0 is the same as
1478 // scaled coordinate of 0.0)
1479 if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1480 }
1481
1482 // find xyz-indices of cell that cutoffGroup is in.
1483 whichCell.x() = int(nCells_.x() * scaled.x());
1484 whichCell.y() = int(nCells_.y() * scaled.y());
1485 whichCell.z() = int(nCells_.z() * scaled.z());
1486
1487 // find single index of this cell:
1488 cellIndex = Vlinear(whichCell, nCells_);
1489
1490 // add this cutoff group to the list of groups in this cell;
1491 cellList_[cellIndex].push_back(i);
1492 }
1493
1494#endif
1495
1496#ifdef IS_MPI
1497 for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1498 rs = cgRowData.position[j1];
1499#else
1500
1501 for (int j1 = 0; j1 < nGroups_; j1++) {
1502 rs = snap_->cgData.position[j1];
1503#endif
1504 point[j1] = len;
1505
1506 // scaled positions relative to the box vectors
1507 scaled = invBox * rs;
1508
1509 // wrap the vector back into the unit box by subtracting integer box
1510 // numbers
1511 for (int j = 0; j < 3; j++) {
1512 scaled[j] -= roundMe(scaled[j]);
1513 scaled[j] += 0.5;
1514 // Handle the special case when an object is exactly on the
1515 // boundary (a scaled coordinate of 1.0 is the same as
1516 // scaled coordinate of 0.0)
1517 if (scaled[j] >= 1.0) scaled[j] -= 1.0;
1518 }
1519
1520 // find xyz-indices of cell that cutoffGroup is in.
1521 whichCell.x() = int(nCells_.x() * scaled.x());
1522 whichCell.y() = int(nCells_.y() * scaled.y());
1523 whichCell.z() = int(nCells_.z() * scaled.z());
1524
1525 for (vector<Vector3i>::iterator os = cellOffsets_.begin();
1526 os != cellOffsets_.end(); ++os) {
1527 Vector3i m2v = whichCell + (*os);
1528
1529 if (m2v.x() >= nCells_.x()) {
1530 m2v.x() = 0;
1531 } else if (m2v.x() < 0) {
1532 m2v.x() = nCells_.x() - 1;
1533 }
1534
1535 if (m2v.y() >= nCells_.y()) {
1536 m2v.y() = 0;
1537 } else if (m2v.y() < 0) {
1538 m2v.y() = nCells_.y() - 1;
1539 }
1540
1541 if (m2v.z() >= nCells_.z()) {
1542 m2v.z() = 0;
1543 } else if (m2v.z() < 0) {
1544 m2v.z() = nCells_.z() - 1;
1545 }
1546 int m2 = Vlinear(m2v, nCells_);
1547#ifdef IS_MPI
1548 for (vector<int>::iterator j2 = cellListCol_[m2].begin();
1549 j2 != cellListCol_[m2].end(); ++j2) {
1550 // In parallel, we need to visit *all* pairs of row
1551 // & column indicies and will divide labor in the
1552 // force evaluation later.
1553 dr = cgColData.position[(*j2)] - rs;
1554 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(dr); }
1555 if (dr.lengthSquare() < rListSq_) {
1556 neighborList.push_back((*j2));
1557 ++len;
1558 }
1559 }
1560#else
1561 for (vector<int>::iterator j2 = cellList_[m2].begin();
1562 j2 != cellList_[m2].end(); ++j2) {
1563 // Always do this if we're in different cells or if
1564 // we're in the same cell and the global index of
1565 // the j2 cutoff group is greater than or equal to
1566 // the j1 cutoff group. Note that Rappaport's code
1567 // has a "less than" conditional here, but that
1568 // deals with atom-by-atom computation. OpenMD
1569 // allows atoms within a single cutoff group to
1570 // interact with each other.
1571
1572 if ((*j2) >= j1) {
1573 dr = snap_->cgData.position[(*j2)] - rs;
1574 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(dr); }
1575 if (dr.lengthSquare() < rListSq_) {
1576 neighborList.push_back((*j2));
1577 ++len;
1578 }
1579 }
1580 }
1581#endif
1582 }
1583 }
1584 } else {
1585 // branch to do all cutoff group pairs
1586#ifdef IS_MPI
1587 for (int j1 = 0; j1 < nGroupsInRow_; j1++) {
1588 point[j1] = len;
1589 rs = cgRowData.position[j1];
1590 for (int j2 = 0; j2 < nGroupsInCol_; j2++) {
1591 dr = cgColData.position[j2] - rs;
1592 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(dr); }
1593 if (dr.lengthSquare() < rListSq_) {
1594 neighborList.push_back(j2);
1595 ++len;
1596 }
1597 }
1598 }
1599#else
1600 // include all groups here.
1601 for (int j1 = 0; j1 < nGroups_; j1++) {
1602 point[j1] = len;
1603 rs = snap_->cgData.position[j1];
1604 // include self group interactions j2 == j1
1605 for (int j2 = j1; j2 < nGroups_; j2++) {
1606 dr = snap_->cgData.position[j2] - rs;
1607 if (usePeriodicBoundaryConditions_) { snap_->wrapVector(dr); }
1608 if (dr.lengthSquare() < rListSq_) {
1609 neighborList.push_back(j2);
1610 ++len;
1611 }
1612 }
1613 }
1614#endif
1615 }
1616
1617#ifdef IS_MPI
1618 point[nGroupsInRow_] = len;
1619#else
1620 point[nGroups_] = len;
1621#endif
1622
1623 // save the local cutoff group positions for the check that is
1624 // done on each loop:
1625 savedPositions.clear();
1626 savedPositions.reserve(nGroups_);
1627 for (int i = 0; i < nGroups_; i++)
1628 savedPositions.push_back(snap_->cgData.position[i]);
1629 }
1630
1631 int ForceMatrixDecomposition::getGlobalIDRow(int atom1) {
1632#ifdef IS_MPI
1633 return AtomRowToGlobal[atom1];
1634#else
1635 return atom1;
1636#endif
1637 }
1638
1639 int ForceMatrixDecomposition::getGlobalIDCol(int atom2) {
1640#ifdef IS_MPI
1641 return AtomColToGlobal[atom2];
1642#else
1643 return atom2;
1644#endif
1645 }
1646
1647 int ForceMatrixDecomposition::getGlobalID(int atom1) {
1648#ifdef IS_MPI
1649 return AtomLocalToGlobal[atom1];
1650#else
1651 return atom1;
1652#endif
1653 }
1654} // namespace OpenMD
vector< RealType > functional
electron density
void resize(std::size_t newSize)
Changes the size of this DataStorage.
vector< RealType > flucQFrc
fluctuating charge velocities
vector< RealType > particlePot
torque array
vector< RealType > sitePotential
fluctuating charge forces
vector< Mat3x3d > quadrupole
space-frame dipole vector
vector< RealType > functionalDerivative
density functional
vector< RealType > flucQPos
charge skipped during normal pairwise calculation
vector< Vector3d > velocity
position array
vector< RotMat3x3d > aMat
force array
vector< RealType > density
particle potential arrray
vector< RealType > skippedCharge
local electric field
void setStorageLayout(int layout)
Sets the storage layout
vector< Vector3d > torque
angular momentum array (body-fixed)
vector< Vector3d > force
velocity array
vector< Vector3d > electricField
space-frame quadrupole tensor
vector< Vector3d > dipole
derivative of functional
ForceDecomposition is an interface for passing out and collecting information from many processors at...
vector< vector< int > > toposForAtom
The topological distance between two atomic sites is handled via two vector structures for speed.
AtomType * getAtomType(const std::string &at)
getAtomType by string
bool excludeAtomPair(int atom1, int atom2)
We need to handle the interactions for atoms who are involved in the same rigid body as well as some ...
vector< int > & getAtomsInGroupRow(int cg1)
returns the list of atoms belonging to this group.
bool skipAtomPair(int atom1, int atom2, int cg1, int cg2)
We need to exclude some overcounted interactions that result from the parallel decomposition.
void distributeInitialData()
distributeInitialData is essentially a copy of the older fortran SimulationSetup
void collectSelfData()
Collects information obtained during the post-pair (and embedding functional) loops onto local data s...
InteractionManager is responsible for keeping track of the non-bonded interactions (C++)
PairList class maintains a general purpose list of atom pairs using the global indices of the atoms.
Definition PairList.hpp:63
bool hasPair(int i, int j)
Checks whether pair (i, j) is in this PairList class.
Definition PairList.cpp:140
Vector< Real, Col > getColumn(unsigned int col)
Returns a column of this matrix as a vector.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
std::vector< int > getGlobalGroupIndices()
returns a vector which maps the local cutoff group index on this processor to the global cutoff group...
Definition SimInfo.cpp:878
ForceField * getForceField()
Returns the force field.
Definition SimInfo.hpp:266
unsigned int getNAtoms()
Returns the number of local atoms.
Definition SimInfo.hpp:172
unsigned int getNCutoffGroups()
Returns the number of local cutoff groups.
Definition SimInfo.hpp:195
unsigned int getNLocalCutoffGroups()
Returns the number of effective cutoff groups on local processor.
Definition SimInfo.cpp:311
std::vector< int > getGlobalAtomIndices()
returns a vector which maps the local atom index on this processor to the global atom index.
Definition SimInfo.cpp:862
The Snapshot class is a repository storing dynamic data during a Simulation.
Definition Snapshot.hpp:147
Mat3x3d getHmat()
Returns the H-Matrix.
Definition Snapshot.cpp:214
Mat3x3d getInvHmat()
Returns the inverse H-Matrix.
Definition Snapshot.cpp:278
Mat3x3d getBoundingBox()
Returns the Bounding Box.
Definition Snapshot.cpp:281
int getNumberOfAtoms()
Returns the number of atoms.
Definition Snapshot.cpp:202
Mat3x3d getInvBoundingBox()
Returns the inverse Bounding Box.
Definition Snapshot.cpp:291
void wrapVector(Vector3d &v)
Wrapping the vector according to periodic boundary condition.
Definition Snapshot.cpp:337
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
Real & z()
Returns reference of the third element of Vector3.
Definition Vector3.hpp:120
Real & x()
Returns reference of the first element of Vector3.
Definition Vector3.hpp:96
Real & y()
Returns reference of the second element of Vector3.
Definition Vector3.hpp:108
Fix length vector class.
Definition Vector.hpp:78
void normalize()
Normalizes this vector in place.
Definition Vector.hpp:402
Real lengthSquare()
Returns the squared length of this vector.
Definition Vector.hpp:399
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
Vector3< Real > cross(const Vector3< Real > &v1, const Vector3< Real > &v2)
Returns the cross product of two Vectors.
Definition Vector3.hpp:136
Real dot(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the dot product of two DynamicVectors.
std::size_t Vlinear(const Vector2< std::size_t > &p, const Vector2< std::size_t > &s)
Returns the linear indexing for size_t vectors.
Definition Vector2.hpp:111
Vector3d conductiveHeatFlux
heat flux vector (conductive only)
Definition Snapshot.hpp:133
The InteractionData struct.
RealType sPot2
site potential on second atom
Vector3d D_1
dipole vector of first atom
RealType skippedCharge1
charge skipped on atom1 in pairwise interaction loop with atom2
RealType dVdFQ1
fluctuating charge force on atom1
RealType dVdFQ2
fluctuating charge force on atom2
RealType particlePot2
particle potential for atom2
int atid1
atomType ident for atom 1
potVec selePot
potential energies of the selected sites
RealType flucQ2
fluctuating charge on atom2
bool sameRegion
are these atoms specified to be in the same region?
bool excluded
is this excluded from direct pairwise interactions? (some indirect interactions may still apply)
RealType rho2
total electron density at second atom
RealType sw
switching function value at rij
Vector3d eField1
electric field on first atom
RotMat3x3d A2
rotation matrix of second atom
RotMat3x3d A1
rotation matrix of first atom
RealType flucQ1
fluctuating charge on atom1
Mat3x3d Q_2
quadrupole tensor of first atom
Vector3d eField2
electric field on second atom
Vector3d t1
torque on first atom
Vector3d t2
torque on second atom
RealType vpair
pair potential
RealType sPot1
site potential on first atom
int atid2
atomType ident for atom 2
Vector3d D_2
dipole vector of first atom
RealType skippedCharge2
charge skipped on atom2 in pairwise interaction loop with atom1
RealType dfrho2
derivative of functional for atom 2
potVec excludedPot
potential energy excluded from the overall calculation
RealType frho2
density functional at second atom
bool doParticlePot
should we bother with the particle pot?
RealType dfrho1
derivative of functional for atom 1
Mat3x3d Q_1
quadrupole tensor of first atom
RealType frho1
density functional at first atom
Vector3d f1
force between the two atoms
RealType particlePot1
particle potential for atom1
RealType rho1
total electron density at first atom