OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
FluctuatingChargeObjectiveFunction.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 "flucq/FluctuatingChargeObjectiveFunction.hpp"
46
47#ifdef IS_MPI
48#include <mpi.h>
49#endif
50
51namespace OpenMD {
52
53 FluctuatingChargeObjectiveFunction::FluctuatingChargeObjectiveFunction(
54 SimInfo* info, ForceManager* forceMan,
55 FluctuatingChargeConstraints* fqConstraints) :
56 info_(info),
57 forceMan_(forceMan), fqConstraints_(fqConstraints), thermo(info) {}
58
60 const DynamicVector<RealType>& x) {
61 setCoor(x);
62 forceMan_->calcForces();
63 fqConstraints_->applyConstraints();
64 return thermo.getPotential();
65 }
66
69 setCoor(x);
70 forceMan_->calcForces();
71 fqConstraints_->applyConstraints();
72 getGrad(grad);
73 }
74
77 setCoor(x);
78 forceMan_->calcForces();
79 fqConstraints_->applyConstraints();
80 getGrad(grad);
81 return thermo.getPotential();
82 }
83
84 void FluctuatingChargeObjectiveFunction::setCoor(
85 const DynamicVector<RealType>& x) const {
86 SimInfo::MoleculeIterator i;
87 Molecule::FluctuatingChargeIterator j;
88 Molecule* mol;
89 Atom* atom;
90
91 info_->getSnapshotManager()->advance();
92
93 int index;
94#ifdef IS_MPI
95 index = displacements_[myrank_];
96#else
97 index = 0;
98#endif
99
100 for (mol = info_->beginMolecule(i); mol != NULL;
101 mol = info_->nextMolecule(i)) {
102 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
103 atom = mol->nextFluctuatingCharge(j)) {
104 atom->setFlucQPos(x[index++]);
105 }
106 }
107 }
108
109 void FluctuatingChargeObjectiveFunction::getGrad(
111 SimInfo::MoleculeIterator i;
112 Molecule::FluctuatingChargeIterator j;
113 Molecule* mol;
114 Atom* atom;
115 grad.setZero();
116
117 int index;
118#ifdef IS_MPI
119 index = displacements_[myrank_];
120#else
121 index = 0;
122#endif
123
124 for (mol = info_->beginMolecule(i); mol != NULL;
125 mol = info_->nextMolecule(i)) {
126 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
127 atom = mol->nextFluctuatingCharge(j)) {
128 grad[index++] = -atom->getFlucQFrc();
129 }
130 }
131
132#ifdef IS_MPI
133 MPI_Allreduce(MPI_IN_PLACE, &grad[0], nFlucQ_, MPI_REALTYPE, MPI_SUM,
134 MPI_COMM_WORLD);
135#endif
136 }
137
139 FluctuatingChargeObjectiveFunction::setInitialCoords() {
140#ifdef IS_MPI
141 MPI_Comm_size(MPI_COMM_WORLD, &nproc_);
142 MPI_Comm_rank(MPI_COMM_WORLD, &myrank_);
143 std::vector<int> flucqOnProc_(nproc_, 0);
144
145 displacements_.clear();
146 displacements_.resize(nproc_, 0);
147#endif
148
149 SimInfo::MoleculeIterator i;
150 Molecule::FluctuatingChargeIterator j;
151 Molecule* mol;
152 Atom* atom;
153
154 nFlucQ_ = 0;
155
156 for (mol = info_->beginMolecule(i); mol != NULL;
157 mol = info_->nextMolecule(i)) {
158 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
159 atom = mol->nextFluctuatingCharge(j)) {
160 nFlucQ_++;
161 }
162 }
163
164#ifdef IS_MPI
165 MPI_Allgather(&nFlucQ_, 1, MPI_INT, &flucqOnProc_[0], 1, MPI_INT,
166 MPI_COMM_WORLD);
167
168 nFlucQ_ = 0;
169 for (int iproc = 0; iproc < nproc_; iproc++) {
170 nFlucQ_ += flucqOnProc_[iproc];
171 }
172
173 displacements_[0] = 0;
174 for (int iproc = 1; iproc < nproc_; iproc++) {
175 displacements_[iproc] =
176 displacements_[iproc - 1] + flucqOnProc_[iproc - 1];
177 }
178#endif
179
180 DynamicVector<RealType> initCoords(nFlucQ_);
181
182 int index;
183#ifdef IS_MPI
184 index = displacements_[myrank_];
185#else
186 index = 0;
187#endif
188
189 for (mol = info_->beginMolecule(i); mol != NULL;
190 mol = info_->nextMolecule(i)) {
191 for (atom = mol->beginFluctuatingCharge(j); atom != NULL;
192 atom = mol->nextFluctuatingCharge(j)) {
193 initCoords[index++] = atom->getFlucQPos();
194 }
195 }
196
197#ifdef IS_MPI
198 MPI_Allgatherv(MPI_IN_PLACE, 0, MPI_DATATYPE_NULL, &initCoords[0],
199 &flucqOnProc_[0], &displacements_[0], MPI_REALTYPE,
200 MPI_COMM_WORLD);
201#endif
202
203 return initCoords;
204 }
205} // namespace OpenMD
Dynamically-sized vector class.
void setZero()
zero out the vector
virtual RealType value(const DynamicVector< RealType > &x)
method to overload to compute the objective function value in x
virtual void gradient(DynamicVector< RealType > &grad, const DynamicVector< RealType > &x)
method to overload to compute grad_f, the first derivative of
virtual RealType valueAndGradient(DynamicVector< RealType > &grad, const DynamicVector< RealType > &x)
method to overload to compute grad_f, the first derivative
ForceManager is responsible for calculating both the short range (bonded) interactions and long range...
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
void setFlucQPos(RealType charge)
Sets the current fluctuating charge of this stuntDouble.
RealType getFlucQFrc()
Returns the current charge force of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.