OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ForceDecomposition.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
45#include "parallel/ForceDecomposition.hpp"
46
47#ifdef IS_MPI
48#include <mpi.h>
49#endif
50
51#include "math/Vector3.hpp"
52
53using namespace std;
54namespace OpenMD {
55
56 ForceDecomposition::ForceDecomposition(SimInfo* info,
57 InteractionManager* iMan) :
58 info_(info),
59 interactionMan_(iMan), needVelocities_(false) {
60 sman_ = info_->getSnapshotManager();
61 atomStorageLayout_ = sman_->getAtomStorageLayout();
62 rigidBodyStorageLayout_ = sman_->getRigidBodyStorageLayout();
63 cutoffGroupStorageLayout_ = sman_->getCutoffGroupStorageLayout();
64 ff_ = info_->getForceField();
65
66 usePeriodicBoundaryConditions_ =
67 info->getSimParams()->getUsePeriodicBoundaryConditions();
68
69 Globals* simParams_ = info_->getSimParams();
70 if (simParams_->havePrintHeatFlux()) {
71 if (simParams_->getPrintHeatFlux()) { needVelocities_ = true; }
72 }
73
74 if (simParams_->haveSkinThickness()) {
75 skinThickness_ = simParams_->getSkinThickness();
76 } else {
77 skinThickness_ = 1.0;
78 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
79 "ForceDecomposition: No value was set for the skinThickness.\n"
80 "\tOpenMD will use a default value of %f Angstroms\n"
81 "\tfor this simulation\n",
82 skinThickness_);
83 painCave.severity = OPENMD_INFO;
84 painCave.isFatal = 0;
85 simError();
86 }
87
88 // cellOffsets are the partial space for the cell lists used in
89 // constructing the neighbor lists
90 cellOffsets_.clear();
91 cellOffsets_.push_back(Vector3i(0, 0, 0));
92 cellOffsets_.push_back(Vector3i(1, 0, 0));
93 cellOffsets_.push_back(Vector3i(1, 1, 0));
94 cellOffsets_.push_back(Vector3i(0, 1, 0));
95 cellOffsets_.push_back(Vector3i(-1, 1, 0));
96 cellOffsets_.push_back(Vector3i(0, 0, 1));
97 cellOffsets_.push_back(Vector3i(1, 0, 1));
98 cellOffsets_.push_back(Vector3i(1, 1, 1));
99 cellOffsets_.push_back(Vector3i(0, 1, 1));
100 cellOffsets_.push_back(Vector3i(-1, 1, 1));
101 cellOffsets_.push_back(Vector3i(-1, 0, 1));
102 cellOffsets_.push_back(Vector3i(-1, -1, 1));
103 cellOffsets_.push_back(Vector3i(0, -1, 1));
104 cellOffsets_.push_back(Vector3i(1, -1, 1));
105 }
106
107 void ForceDecomposition::setCutoffRadius(RealType rcut) {
108 rCut_ = rcut;
109 rList_ = rCut_ + skinThickness_;
110 rListSq_ = rList_ * rList_;
111 }
112
113 void ForceDecomposition::fillPreForceData(SelfData& sdat, int atom) {
114 sdat.atid = idents[atom];
115 sdat.selfPot = selfPot;
116 sdat.selePot = selectedSelfPot;
117
118 if (atomStorageLayout_ & DataStorage::dslDensity) {
119 sdat.rho = snap_->atomData.density[atom];
120 }
121
122 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
123 sdat.particlePot = snap_->atomData.particlePot[atom];
124 }
125 }
126
127 void ForceDecomposition::fillSelfData(SelfData& sdat, int atom) {
128 sdat.atid = idents[atom];
129 sdat.selfPot = selfPot;
130 sdat.selePot = selectedSelfPot;
131
132 if (atomStorageLayout_ & DataStorage::dslSkippedCharge) {
133 sdat.skippedCharge = snap_->atomData.skippedCharge[atom];
134 }
135
136 if (atomStorageLayout_ & DataStorage::dslParticlePot) {
137 sdat.particlePot = snap_->atomData.particlePot[atom];
138 }
139
140 if (atomStorageLayout_ & DataStorage::dslFlucQPosition) {
141 sdat.flucQ = snap_->atomData.flucQPos[atom];
142 }
143
144 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
145 sdat.flucQfrc = snap_->atomData.flucQFrc[atom];
146 }
147 }
148
149 void ForceDecomposition::unpackPreForceData(SelfData& sdat, int atom) {
150 selfPot = sdat.selfPot;
151 selectedSelfPot = sdat.selePot;
152
153 if (atomStorageLayout_ & DataStorage::dslFunctional) {
154 snap_->atomData.functional[atom] += sdat.frho;
155 }
156
157 if (atomStorageLayout_ & DataStorage::dslFunctionalDerivative) {
158 snap_->atomData.functionalDerivative[atom] += sdat.dfrhodrho;
159 }
160
161 if (sdat.doParticlePot &&
162 (atomStorageLayout_ & DataStorage::dslParticlePot)) {
163 snap_->atomData.particlePot[atom] = sdat.particlePot;
164 }
165 }
166
167 void ForceDecomposition::unpackSelfData(SelfData& sdat, int atom) {
168 selfPot = sdat.selfPot;
169 selectedSelfPot = sdat.selePot;
170
171 if (atomStorageLayout_ & DataStorage::dslFlucQForce) {
172 snap_->atomData.flucQFrc[atom] = sdat.flucQfrc;
173 }
174 }
175
176 bool ForceDecomposition::checkNeighborList(vector<Vector3d> savedPositions) {
177 RealType st2 = pow(skinThickness_ / 2.0, 2);
178 std::size_t nGroups = snap_->cgData.position.size();
179 if (needVelocities_)
180 snap_->cgData.setStorageLayout(DataStorage::dslPosition |
181 DataStorage::dslVelocity);
182
183 // if we have changed the group identities or haven't set up the
184 // saved positions we automatically will need a neighbor list update:
185
186 if (savedPositions.size() != nGroups) return true;
187
188 RealType dispmax = 0.0;
189 Vector3d disp;
190
191 for (std::size_t i = 0; i < nGroups; i++) {
192 disp = snap_->cgData.position[i] - savedPositions[i];
193 dispmax = max(dispmax, disp.lengthSquare());
194 }
195
196#ifdef IS_MPI
197 MPI_Allreduce(MPI_IN_PLACE, &dispmax, 1, MPI_REALTYPE, MPI_MAX,
198 MPI_COMM_WORLD);
199#endif
200
201 return (dispmax > st2) ? true : false;
202 }
203
204 void ForceDecomposition::addToHeatFlux(Vector3d hf) {
205 Vector3d chf = snap_->getConductiveHeatFlux();
206 chf += hf;
207 snap_->setConductiveHeatFlux(chf);
208 }
209 void ForceDecomposition::setHeatFlux(Vector3d hf) {
210 snap_->setConductiveHeatFlux(hf);
211 }
212} // namespace OpenMD
vector< RealType > functional
electron density
vector< RealType > flucQFrc
fluctuating charge velocities
vector< RealType > particlePot
torque array
vector< RealType > functionalDerivative
density functional
vector< RealType > flucQPos
charge skipped during normal pairwise calculation
vector< RealType > density
particle potential arrray
vector< RealType > skippedCharge
local electric field
void setStorageLayout(int layout)
Sets the storage layout
RealType skinThickness_
Verlet neighbor list skin thickness.
InteractionManager is responsible for keeping track of the non-bonded interactions (C++)
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
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.
The SelfData struct.
potVec selePot
potential energy of the selected site
RealType flucQfrc
fluctuating charge derivative
RealType frho
value of density functional for atom
RealType dfrhodrho
derivative of density functional for atom
potVec selfPot
total potential (including embedding energy)
RealType particlePot
contribution to potential from this particle
RealType rho
electron density
RealType flucQ
current value of atom's fluctuating charge
RealType skippedCharge
charge skipped in pairwise interaction loop
bool doParticlePot
should we bother with the particle pot?
int atid
atomType ident for the atom