OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
PotentialEnergyObjectiveFunction.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 "optimization/PotentialEnergyObjectiveFunction.hpp"
46
47#ifdef IS_MPI
48#include <mpi.h>
49#endif
50
51namespace OpenMD {
52
53 PotentialEnergyObjectiveFunction::PotentialEnergyObjectiveFunction(
54 SimInfo* info, ForceManager* forceMan) :
55 info_(info),
56 forceMan_(forceMan), thermo(info), hasFlucQ_(false) {
57 shake_ = new Shake(info_);
58
59 if (info_->usesFluctuatingCharges()) {
60 if (info_->getNFluctuatingCharges() > 0) {
61 hasFlucQ_ = true;
62 fqConstraints_ = new FluctuatingChargeConstraints(info_);
63 bool cr = info_->getSimParams()
64 ->getFluctuatingChargeParameters()
65 ->getConstrainRegions();
66 fqConstraints_->setConstrainRegions(cr);
67 }
68 }
69 }
70
72 const DynamicVector<RealType>& x) {
73 setCoor(x);
74 shake_->constraintR();
75 forceMan_->calcForces();
76 if (hasFlucQ_) fqConstraints_->applyConstraints();
77 shake_->constraintF();
78 return thermo.getPotential();
79 }
80
83 setCoor(x);
84 shake_->constraintR();
85 forceMan_->calcForces();
86 if (hasFlucQ_) fqConstraints_->applyConstraints();
87 shake_->constraintF();
88 getGrad(grad);
89 }
90
93 setCoor(x);
94 shake_->constraintR();
95 forceMan_->calcForces();
96 if (hasFlucQ_) fqConstraints_->applyConstraints();
97 shake_->constraintF();
98 getGrad(grad);
99 return thermo.getPotential();
100 }
101
102 void PotentialEnergyObjectiveFunction::setCoor(
103 const DynamicVector<RealType>& x) const {
104 Vector3d position;
105 Vector3d eulerAngle;
106 SimInfo::MoleculeIterator i;
107 Molecule::IntegrableObjectIterator j;
108 Molecule::AtomIterator ai;
109 Molecule* mol;
110 StuntDouble* sd;
111 Atom* atom;
112
113 info_->getSnapshotManager()->advance();
114
115 int index;
116#ifdef IS_MPI
117 index = displacements_[myrank_];
118#else
119 index = 0;
120#endif
121
122 info_->getSnapshotManager()->advance();
123
124 for (mol = info_->beginMolecule(i); mol != NULL;
125 mol = info_->nextMolecule(i)) {
126 for (sd = mol->beginIntegrableObject(j); sd != NULL;
127 sd = mol->nextIntegrableObject(j)) {
128 position[0] = x[index++];
129 position[1] = x[index++];
130 position[2] = x[index++];
131
132 sd->setPos(position);
133
134 if (sd->isDirectional()) {
135 eulerAngle[0] = x[index++];
136 eulerAngle[1] = x[index++];
137 eulerAngle[2] = x[index++];
138
139 sd->setEuler(eulerAngle);
140
141 if (sd->isRigidBody()) {
142 RigidBody* rb = static_cast<RigidBody*>(sd);
143 rb->updateAtoms();
144 }
145 }
146 }
147 }
148
149 if (hasFlucQ_) {
150 for (mol = info_->beginMolecule(i); mol != NULL;
151 mol = info_->nextMolecule(i)) {
152 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
153 atom = mol->nextFluctuatingCharge(ai)) {
154 atom->setFlucQPos(x[index++]);
155 }
156 }
157 }
158 }
159
160 void PotentialEnergyObjectiveFunction::getGrad(
162 SimInfo::MoleculeIterator i;
163 Molecule::IntegrableObjectIterator j;
164 Molecule::AtomIterator ai;
165 Molecule* mol;
166 StuntDouble* sd;
167 Atom* atom;
168 std::vector<RealType> myGrad;
169
170 int index;
171#ifdef IS_MPI
172 index = displacements_[myrank_];
173 grad.setZero();
174#else
175 index = 0;
176#endif
177
178 for (mol = info_->beginMolecule(i); mol != NULL;
179 mol = info_->nextMolecule(i)) {
180 for (sd = mol->beginIntegrableObject(j); sd != NULL;
181 sd = mol->nextIntegrableObject(j)) {
182 myGrad = sd->getGrad();
183
184 for (size_t k = 0; k < myGrad.size(); ++k) {
185 grad[index++] = myGrad[k];
186 }
187 }
188 }
189
190 if (hasFlucQ_) {
191 for (mol = info_->beginMolecule(i); mol != NULL;
192 mol = info_->nextMolecule(i)) {
193 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
194 atom = mol->nextFluctuatingCharge(ai)) {
195 grad[index++] = -atom->getFlucQFrc();
196 }
197 }
198 }
199#ifdef IS_MPI
200 MPI_Allreduce(MPI_IN_PLACE, &grad[0], ndf_, MPI_REALTYPE, MPI_SUM,
201 MPI_COMM_WORLD);
202#endif
203 }
204
205 DynamicVector<RealType> PotentialEnergyObjectiveFunction::setInitialCoords() {
206#ifdef IS_MPI
207 MPI_Comm_size(MPI_COMM_WORLD, &nproc_);
208 MPI_Comm_rank(MPI_COMM_WORLD, &myrank_);
209 std::vector<int> onProc(nproc_, 0);
210
211 displacements_.clear();
212 displacements_.resize(nproc_, 0);
213#endif
214
215 SimInfo::MoleculeIterator i;
216 Molecule::IntegrableObjectIterator j;
217 Molecule::AtomIterator ai;
218 Molecule* mol;
219 StuntDouble* sd;
220 Atom* atom;
221
222 Vector3d pos;
223 Vector3d eulerAngle;
224
225 ndf_ = 0;
226
227 for (mol = info_->beginMolecule(i); mol != NULL;
228 mol = info_->nextMolecule(i)) {
229 for (sd = mol->beginIntegrableObject(j); sd != NULL;
230 sd = mol->nextIntegrableObject(j)) {
231 ndf_ += 3;
232
233 if (sd->isDirectional()) { ndf_ += 3; }
234 }
235 }
236
237 if (hasFlucQ_) {
238 for (mol = info_->beginMolecule(i); mol != NULL;
239 mol = info_->nextMolecule(i)) {
240 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
241 atom = mol->nextFluctuatingCharge(ai)) {
242 ndf_++;
243 }
244 }
245 }
246
247#ifdef IS_MPI
248 MPI_Allgather(&ndf_, 1, MPI_INT, &onProc[0], 1, MPI_INT, MPI_COMM_WORLD);
249
250 ndf_ = 0;
251 for (int iproc = 0; iproc < nproc_; iproc++) {
252 ndf_ += onProc[iproc];
253 }
254
255 displacements_[0] = 0;
256 for (int iproc = 1; iproc < nproc_; iproc++) {
257 displacements_[iproc] = displacements_[iproc - 1] + onProc[iproc - 1];
258 }
259#endif
260
261 DynamicVector<RealType> xinit(ndf_, 0.0);
262
263 int index;
264#ifdef IS_MPI
265 index = displacements_[myrank_];
266#else
267 index = 0;
268#endif
269
270 for (mol = info_->beginMolecule(i); mol != NULL;
271 mol = info_->nextMolecule(i)) {
272 for (sd = mol->beginIntegrableObject(j); sd != NULL;
273 sd = mol->nextIntegrableObject(j)) {
274 pos = sd->getPos();
275 xinit[index++] = pos[0];
276 xinit[index++] = pos[1];
277 xinit[index++] = pos[2];
278
279 if (sd->isDirectional()) {
280 eulerAngle = sd->getEuler();
281 xinit[index++] = eulerAngle[0];
282 xinit[index++] = eulerAngle[1];
283 xinit[index++] = eulerAngle[2];
284 }
285 }
286 }
287
288 if (hasFlucQ_) {
289 for (mol = info_->beginMolecule(i); mol != NULL;
290 mol = info_->nextMolecule(i)) {
291 for (atom = mol->beginFluctuatingCharge(ai); atom != NULL;
292 atom = mol->nextFluctuatingCharge(ai)) {
293 xinit[index++] = atom->getFlucQPos();
294 }
295 }
296 }
297
298 return xinit;
299 }
300} // namespace OpenMD
Dynamically-sized vector class.
void setZero()
zero out the vector
ForceManager is responsible for calculating both the short range (bonded) interactions and long range...
RealType value(const DynamicVector< RealType > &x)
method to overload to compute the objective function value in x
void gradient(DynamicVector< RealType > &grad, const DynamicVector< RealType > &x)
method to overload to compute grad_f, the first derivative of
RealType valueAndGradient(DynamicVector< RealType > &grad, const DynamicVector< RealType > &x)
method to overload to compute grad_f, the first derivative
void updateAtoms()
update the positions of atoms belong to this rigidbody
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Definition SimInfo.cpp:240
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
Definition SimInfo.cpp:245
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
Definition SimInfo.hpp:248
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getEuler()
Returns the current euler angles of this stuntDouble.
void setFlucQPos(RealType charge)
Sets the current fluctuating charge of this stuntDouble.
Vector3d getPos()
Returns the current position of this stuntDouble.
RealType getFlucQFrc()
Returns the current charge force of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
void setPos(const Vector3d &pos)
Sets the current position of this stuntDouble.
void setEuler(const Vector3d &euler)
Sets the current euler angles of this stuntDouble.
bool isRigidBody()
Tests if this stuntDouble is a rigid body.
bool isDirectional()
Tests if this stuntDouble is a directional one.
virtual std::vector< RealType > getGrad()=0
Returns the gradient of this stuntDouble.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.