OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
UniformGradient.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
46
47#ifdef IS_MPI
48#include <mpi.h>
49#endif
50
51#include "brains/ForceModifier.hpp"
52#include "nonbonded/NonBondedInteraction.hpp"
54#include "types/FixedChargeAdapter.hpp"
55#include "types/FluctuatingChargeAdapter.hpp"
56#include "types/MultipoleAdapter.hpp"
57#include "utils/Constants.hpp"
58
59namespace OpenMD {
60 UniformGradient::UniformGradient(SimInfo* info) :
61 ForceModifier {info}, initialized {false}, doUniformGradient {false},
62 doParticlePot {false} {
63 simParams = info_->getSimParams();
64 }
65
66 void UniformGradient::initialize() {
67 bool haveA = false;
68 bool haveB = false;
69 bool haveG = false;
70
71 if (simParams->haveUniformGradientDirection1()) {
72 std::vector<RealType> d1 = simParams->getUniformGradientDirection1();
73 if (d1.size() != 3) {
74 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
75 "uniformGradientDirection1: Incorrect number of parameters\n"
76 "\tspecified. There should be 3 parameters, but %lu were\n"
77 "\tspecified.\n",
78 d1.size());
79 painCave.isFatal = 1;
80 simError();
81 }
82 a_.x() = d1[0];
83 a_.y() = d1[1];
84 a_.z() = d1[2];
85
86 a_.normalize();
87 haveA = true;
88 }
89
90 if (simParams->haveUniformGradientDirection2()) {
91 std::vector<RealType> d2 = simParams->getUniformGradientDirection2();
92 if (d2.size() != 3) {
93 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
94 "uniformGradientDirection2: Incorrect number of parameters\n"
95 "\tspecified. There should be 3 parameters, but %lu were\n"
96 "\tspecified.\n",
97 d2.size());
98 painCave.isFatal = 1;
99 simError();
100 }
101 b_.x() = d2[0];
102 b_.y() = d2[1];
103 b_.z() = d2[2];
104
105 b_.normalize();
106 haveB = true;
107 }
108
109 if (simParams->haveUniformGradientStrength()) {
110 g_ = simParams->getUniformGradientStrength();
111 haveG = true;
112 }
113
114 if (haveA && haveB && haveG) {
115 doUniformGradient = true;
116 cpsi_ = dot(a_, b_);
117
118 Grad_(0, 0) = 2.0 * (a_.x() * b_.x() - cpsi_ / 3.0);
119 Grad_(0, 1) = a_.x() * b_.y() + a_.y() * b_.x();
120 Grad_(0, 2) = a_.x() * b_.z() + a_.z() * b_.x();
121 Grad_(1, 0) = Grad_(0, 1);
122 Grad_(1, 1) = 2.0 * (a_.y() * b_.y() - cpsi_ / 3.0);
123 Grad_(1, 2) = a_.y() * b_.z() + a_.z() * b_.y();
124 Grad_(2, 0) = Grad_(0, 2);
125 Grad_(2, 1) = Grad_(1, 2);
126 Grad_(2, 2) = 2.0 * (a_.z() * b_.z() - cpsi_ / 3.0);
127
128 Grad_ *= g_ / 2.0;
129
130 } else {
131 if (!haveA) {
132 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
133 "UniformGradient: uniformGradientDirection1 not specified.\n");
134 painCave.isFatal = 1;
135 simError();
136 }
137 if (!haveB) {
138 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
139 "UniformGradient: uniformGradientDirection2 not specified.\n");
140 painCave.isFatal = 1;
141 simError();
142 }
143 if (!haveG) {
144 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
145 "UniformGradient: uniformGradientStrength not specified.\n");
146 painCave.isFatal = 1;
147 simError();
148 }
149 }
150
151 int storageLayout_ = info_->getSnapshotManager()->getAtomStorageLayout();
152 if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot = true;
153 initialized = true;
154 }
155
156 void UniformGradient::modifyForces() {
157 if (!initialized) initialize();
158
159 SimInfo::MoleculeIterator i;
160 Molecule::AtomIterator j;
161 Molecule* mol;
162 Atom* atom;
163 AtomType* atype;
164 potVec longRangePotential(0.0);
165
166 RealType C;
167 Vector3d D;
168 Mat3x3d Q;
169
170 RealType U;
171 RealType fPot;
172 Vector3d t;
173 Vector3d f;
174
175 Vector3d r;
176 Vector3d EF;
177
178 bool isCharge;
179
180 if (doUniformGradient) {
181 U = 0.0;
182 fPot = 0.0;
183
184 for (mol = info_->beginMolecule(i); mol != NULL;
185 mol = info_->nextMolecule(i)) {
186 for (atom = mol->beginAtom(j); atom != NULL; atom = mol->nextAtom(j)) {
187 isCharge = false;
188 C = 0.0;
189
190 atype = atom->getAtomType();
191
192 r = atom->getPos();
194
195 EF = Grad_ * r;
196
197 atom->addElectricField(EF * Constants::chargeFieldConvert);
198
200 if (fca.isFixedCharge()) {
201 isCharge = true;
202 C = fca.getCharge();
203 }
204
206 if (fqa.isFluctuatingCharge()) {
207 isCharge = true;
208 C += atom->getFlucQPos();
209 atom->addFlucQFrc(dot(r, EF) * Constants::chargeFieldConvert);
210 }
211
212 if (isCharge) {
213 f = EF * C * Constants::chargeFieldConvert;
214 atom->addFrc(f);
215
216 U = -dot(r, f);
217 if (doParticlePot) { atom->addParticlePot(U); }
218 fPot += U;
219 }
220
222 if (ma.isDipole()) {
223 D = atom->getDipole() * Constants::dipoleFieldConvert;
224
225 f = D * Grad_;
226 atom->addFrc(f);
227
228 t = cross(D, EF);
229 atom->addTrq(t);
230
231 U = -dot(D, EF);
232 if (doParticlePot) { atom->addParticlePot(U); }
233 fPot += U;
234 }
235
236 if (ma.isQuadrupole()) {
237 Q = atom->getQuadrupole() * Constants::dipoleFieldConvert;
238
239 t = 2.0 * mCross(Q, Grad_);
240 atom->addTrq(t);
241
242 U = -doubleDot(Q, Grad_);
243 if (doParticlePot) { atom->addParticlePot(U); }
244 fPot += U;
245 }
246 }
247 }
248
249#ifdef IS_MPI
250 MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE, MPI_SUM,
251 MPI_COMM_WORLD);
252#endif
253
255 longRangePotential = snap->getLongRangePotentials();
256 longRangePotential[ELECTROSTATIC_FAMILY] += fPot;
257 snap->setLongRangePotentials(longRangePotential);
258 }
259 }
260} // namespace OpenMD
Uniform Electric Field Gradient perturbation.
AtomType * getAtomType()
Returns the AtomType of this Atom.
Definition A.hpp:93
AtomType is what OpenMD looks to for unchanging data about an atom.
Definition AtomType.hpp:66
Abstract class for external ForceModifier classes.
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
The Snapshot class is a repository storing dynamic data during a Simulation.
Definition Snapshot.hpp:147
void wrapVector(Vector3d &v)
Wrapping the vector according to periodic boundary condition.
Definition Snapshot.cpp:337
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
void addFlucQFrc(RealType cfrc)
Adds charge force into the current charge force of this stuntDouble.
void addElectricField(const Vector3d &eField)
Adds electric field into the current electric field of this stuntDouble.
Mat3x3d getQuadrupole()
Returns the current quadrupole tensor of this stuntDouble.
Vector3d getPos()
Returns the current position of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
Vector3d getDipole()
Returns the current dipole vector of this stuntDouble.
void addTrq(const Vector3d &trq)
Adds torque into the current torque of this stuntDouble.
void addParticlePot(const RealType &particlePot)
Adds particlePot into the current particlePot of this stuntDouble.
void addFrc(const Vector3d &frc)
Adds force into the current force of this stuntDouble.
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
void normalize()
Normalizes this vector in place.
Definition Vector.hpp:402
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
Real doubleDot(const RectMatrix< Real, Row, Col > &t1, const RectMatrix< Real, Row, Col > &t2)
Returns the tensor contraction (double dot product) of two rank 2 tensors (or Matrices)
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.
Vector< Real, Row > mCross(const RectMatrix< Real, Row, Col > &t1, const RectMatrix< Real, Row, Col > &t2)
Returns the vector (cross) product of two matrices.
@ ELECTROSTATIC_FAMILY
Coulombic and point-multipole interactions.