| 1 |
/* |
| 2 |
* Copyright (c) 2014 The University of Notre Dame. All Rights Reserved. |
| 3 |
* |
| 4 |
* The University of Notre Dame grants you ("Licensee") a |
| 5 |
* non-exclusive, royalty free, license to use, modify and |
| 6 |
* redistribute this software in source and binary code form, provided |
| 7 |
* that the following conditions are met: |
| 8 |
* |
| 9 |
* 1. Redistributions of source code must retain the above copyright |
| 10 |
* notice, this list of conditions and the following disclaimer. |
| 11 |
* |
| 12 |
* 2. Redistributions in binary form must reproduce the above copyright |
| 13 |
* notice, this list of conditions and the following disclaimer in the |
| 14 |
* documentation and/or other materials provided with the |
| 15 |
* distribution. |
| 16 |
* |
| 17 |
* This software is provided "AS IS," without a warranty of any |
| 18 |
* kind. All express or implied conditions, representations and |
| 19 |
* warranties, including any implied warranty of merchantability, |
| 20 |
* fitness for a particular purpose or non-infringement, are hereby |
| 21 |
* excluded. The University of Notre Dame and its licensors shall not |
| 22 |
* be liable for any damages suffered by licensee as a result of |
| 23 |
* using, modifying or distributing the software or its |
| 24 |
* derivatives. In no event will the University of Notre Dame or its |
| 25 |
* licensors be liable for any lost revenue, profit or data, or for |
| 26 |
* direct, indirect, special, consequential, incidental or punitive |
| 27 |
* damages, however caused and regardless of the theory of liability, |
| 28 |
* arising out of the use of or inability to use software, even if the |
| 29 |
* University of Notre Dame has been advised of the possibility of |
| 30 |
* such damages. |
| 31 |
* |
| 32 |
* SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
| 33 |
* research, please cite the appropriate papers when you publish your |
| 34 |
* work. Good starting points are: |
| 35 |
* |
| 36 |
* [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
| 37 |
* [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
| 38 |
* [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008). |
| 39 |
* [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010). |
| 40 |
* [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). |
| 41 |
*/ |
| 42 |
|
| 43 |
#include "perturbations/UniformGradient.hpp" |
| 44 |
#include "types/FixedChargeAdapter.hpp" |
| 45 |
#include "types/FluctuatingChargeAdapter.hpp" |
| 46 |
#include "types/MultipoleAdapter.hpp" |
| 47 |
#include "primitives/Molecule.hpp" |
| 48 |
#include "nonbonded/NonBondedInteraction.hpp" |
| 49 |
#include "utils/PhysicalConstants.hpp" |
| 50 |
|
| 51 |
namespace OpenMD { |
| 52 |
|
| 53 |
UniformGradient::UniformGradient(SimInfo* info) : initialized(false), |
| 54 |
doUniformGradient(false), |
| 55 |
doParticlePot(false), |
| 56 |
info_(info) { |
| 57 |
simParams = info_->getSimParams(); |
| 58 |
} |
| 59 |
|
| 60 |
void UniformGradient::initialize() { |
| 61 |
|
| 62 |
bool haveA = false; |
| 63 |
bool haveB = false; |
| 64 |
bool haveG = false; |
| 65 |
|
| 66 |
if (simParams->haveUniformGradientDirection1()) { |
| 67 |
std::vector<RealType> d1 = simParams->getUniformGradientDirection1(); |
| 68 |
if (d1.size() != 3) { |
| 69 |
sprintf(painCave.errMsg, |
| 70 |
"uniformGradientDirection1: Incorrect number of parameters\n" |
| 71 |
"\tspecified. There should be 3 parameters, but %lu were\n" |
| 72 |
"\tspecified.\n", d1.size()); |
| 73 |
painCave.isFatal = 1; |
| 74 |
simError(); |
| 75 |
} |
| 76 |
a_.x() = d1[0]; |
| 77 |
a_.y() = d1[1]; |
| 78 |
a_.z() = d1[2]; |
| 79 |
|
| 80 |
a_.normalize(); |
| 81 |
haveA = true; |
| 82 |
} |
| 83 |
|
| 84 |
if (simParams->haveUniformGradientDirection2()) { |
| 85 |
std::vector<RealType> d2 = simParams->getUniformGradientDirection2(); |
| 86 |
if (d2.size() != 3) { |
| 87 |
sprintf(painCave.errMsg, |
| 88 |
"uniformGradientDirection2: Incorrect number of parameters\n" |
| 89 |
"\tspecified. There should be 3 parameters, but %lu were\n" |
| 90 |
"\tspecified.\n", d2.size()); |
| 91 |
painCave.isFatal = 1; |
| 92 |
simError(); |
| 93 |
} |
| 94 |
b_.x() = d2[0]; |
| 95 |
b_.y() = d2[1]; |
| 96 |
b_.z() = d2[2]; |
| 97 |
|
| 98 |
b_.normalize(); |
| 99 |
haveB = true; |
| 100 |
} |
| 101 |
|
| 102 |
if (simParams->haveUniformGradientStrength()) { |
| 103 |
g_ = simParams->getUniformGradientStrength(); |
| 104 |
haveG = true; |
| 105 |
} |
| 106 |
|
| 107 |
if (haveA && haveB && haveG) { |
| 108 |
doUniformGradient = true; |
| 109 |
cpsi_ = dot(a_, b_); |
| 110 |
|
| 111 |
Grad_(0,0) = 2.0 * (a_.x()*b_.x() - cpsi_ / 3.0); |
| 112 |
Grad_(0,1) = a_.x()*b_.y() + a_.y()*b_.x(); |
| 113 |
Grad_(0,2) = a_.x()*b_.z() + a_.z()*b_.x(); |
| 114 |
Grad_(1,0) = Grad_(0,1); |
| 115 |
Grad_(1,1) = 2.0 * (a_.y()*b_.y() - cpsi_ / 3.0); |
| 116 |
Grad_(1,2) = a_.y()*b_.z() + a_.z()*b_.y(); |
| 117 |
Grad_(2,0) = Grad_(0,2); |
| 118 |
Grad_(2,1) = Grad_(1,2); |
| 119 |
Grad_(2,2) = 2.0 * (a_.z()*b_.z() - cpsi_ / 3.0); |
| 120 |
|
| 121 |
Grad_ *= g_ / 2.0; |
| 122 |
|
| 123 |
} else { |
| 124 |
if (!haveA) { |
| 125 |
sprintf(painCave.errMsg, |
| 126 |
"UniformGradient: uniformGradientDirection1 not specified.\n"); |
| 127 |
painCave.isFatal = 1; |
| 128 |
simError(); |
| 129 |
} |
| 130 |
if (!haveB) { |
| 131 |
sprintf(painCave.errMsg, |
| 132 |
"UniformGradient: uniformGradientDirection2 not specified.\n"); |
| 133 |
painCave.isFatal = 1; |
| 134 |
simError(); |
| 135 |
} |
| 136 |
if (!haveG) { |
| 137 |
sprintf(painCave.errMsg, |
| 138 |
"UniformGradient: uniformGradientStrength not specified.\n"); |
| 139 |
painCave.isFatal = 1; |
| 140 |
simError(); |
| 141 |
} |
| 142 |
} |
| 143 |
|
| 144 |
int storageLayout_ = info_->getSnapshotManager()->getStorageLayout(); |
| 145 |
if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot = true; |
| 146 |
initialized = true; |
| 147 |
} |
| 148 |
|
| 149 |
void UniformGradient::applyPerturbation() { |
| 150 |
|
| 151 |
if (!initialized) initialize(); |
| 152 |
|
| 153 |
SimInfo::MoleculeIterator i; |
| 154 |
Molecule::AtomIterator j; |
| 155 |
Molecule* mol; |
| 156 |
Atom* atom; |
| 157 |
AtomType* atype; |
| 158 |
potVec longRangePotential(0.0); |
| 159 |
|
| 160 |
RealType C; |
| 161 |
Vector3d D; |
| 162 |
Mat3x3d Q; |
| 163 |
|
| 164 |
RealType U; |
| 165 |
RealType fPot; |
| 166 |
Vector3d t; |
| 167 |
Vector3d f; |
| 168 |
|
| 169 |
Vector3d r; |
| 170 |
Vector3d EF; |
| 171 |
|
| 172 |
bool isCharge; |
| 173 |
|
| 174 |
if (doUniformGradient) { |
| 175 |
|
| 176 |
U = 0.0; |
| 177 |
fPot = 0.0; |
| 178 |
|
| 179 |
for (mol = info_->beginMolecule(i); mol != NULL; |
| 180 |
mol = info_->nextMolecule(i)) { |
| 181 |
|
| 182 |
for (atom = mol->beginAtom(j); atom != NULL; |
| 183 |
atom = mol->nextAtom(j)) { |
| 184 |
|
| 185 |
isCharge = false; |
| 186 |
C = 0.0; |
| 187 |
|
| 188 |
atype = atom->getAtomType(); |
| 189 |
|
| 190 |
r = atom->getPos(); |
| 191 |
EF = Grad_ * r; |
| 192 |
|
| 193 |
if (atype->isElectrostatic()) { |
| 194 |
atom->addElectricField(EF * PhysicalConstants::chargeFieldConvert); |
| 195 |
} |
| 196 |
|
| 197 |
FixedChargeAdapter fca = FixedChargeAdapter(atype); |
| 198 |
if ( fca.isFixedCharge() ) { |
| 199 |
isCharge = true; |
| 200 |
C = fca.getCharge(); |
| 201 |
} |
| 202 |
|
| 203 |
FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter(atype); |
| 204 |
if ( fqa.isFluctuatingCharge() ) { |
| 205 |
isCharge = true; |
| 206 |
C += atom->getFlucQPos(); |
| 207 |
atom->addFlucQFrc( dot(r, EF) |
| 208 |
* PhysicalConstants::chargeFieldConvert ); |
| 209 |
} |
| 210 |
|
| 211 |
if (isCharge) { |
| 212 |
f = EF * C * PhysicalConstants::chargeFieldConvert; |
| 213 |
atom->addFrc(f); |
| 214 |
|
| 215 |
U = -dot(r, f); |
| 216 |
if (doParticlePot) { |
| 217 |
atom->addParticlePot(U); |
| 218 |
} |
| 219 |
fPot += U; |
| 220 |
} |
| 221 |
|
| 222 |
MultipoleAdapter ma = MultipoleAdapter(atype); |
| 223 |
if (ma.isDipole() ) { |
| 224 |
D = atom->getDipole() * PhysicalConstants::dipoleFieldConvert; |
| 225 |
|
| 226 |
f = D * Grad_; |
| 227 |
atom->addFrc(f); |
| 228 |
|
| 229 |
t = cross(D, EF); |
| 230 |
atom->addTrq(t); |
| 231 |
|
| 232 |
U = -dot(D, EF); |
| 233 |
if (doParticlePot) { |
| 234 |
atom->addParticlePot(U); |
| 235 |
} |
| 236 |
fPot += U; |
| 237 |
} |
| 238 |
|
| 239 |
if (ma.isQuadrupole() ) { |
| 240 |
Q = atom->getQuadrupole() * PhysicalConstants::dipoleFieldConvert; |
| 241 |
|
| 242 |
t = 2.0 * mCross(Q, Grad_); |
| 243 |
atom->addTrq(t); |
| 244 |
|
| 245 |
U = -doubleDot(Q, Grad_); |
| 246 |
if (doParticlePot) { |
| 247 |
atom->addParticlePot(U); |
| 248 |
} |
| 249 |
fPot += U; |
| 250 |
} |
| 251 |
} |
| 252 |
} |
| 253 |
|
| 254 |
#ifdef IS_MPI |
| 255 |
MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE, |
| 256 |
MPI_SUM, MPI_COMM_WORLD); |
| 257 |
#endif |
| 258 |
|
| 259 |
Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot(); |
| 260 |
longRangePotential = snap->getLongRangePotentials(); |
| 261 |
longRangePotential[ELECTROSTATIC_FAMILY] += fPot; |
| 262 |
snap->setLongRangePotential(longRangePotential); |
| 263 |
} |
| 264 |
} |
| 265 |
} |