OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
Light.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::Perturbations {
60 Light::Light(SimInfo* info) :
61 ForceModifier {info}, initialized {false}, doLight {false},
62 doParticlePot {false}, info_(info) {
63 lightParams = info_->getSimParams()->getLightParameters();
64 }
65
66 void Light::initialize() {
67 bool haveE0 = false;
68 bool haveDirection = false;
69 bool haveFrequency = false;
70 bool havePolarization = false;
71
72 if (lightParams->haveWaveVector()) {
73 std::vector<RealType> k = lightParams->getWaveVector();
74 // wave vectors are input in inverse angstroms, so no unit conversion:
75 k_.x() = k[0];
76 k_.y() = k[1];
77 k_.z() = k[2];
78 kmag_ = k_.length();
79 lambda_ = 2.0 * Constants::PI / kmag_;
80 omega_ = 2.0 * Constants::PI * Constants::c / lambda_;
81 haveFrequency = true;
82 khat_ = k_;
83 khat_.normalize();
84 haveDirection = true;
85 }
86
87 if (lightParams->havePropagationDirection()) {
88 if (haveDirection) {
89 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
90 "light: please specify either waveVector or "
91 "propagationDirection, but not both.\n");
92 painCave.isFatal = 1;
93 simError();
94 }
95 std::vector<RealType> pd = lightParams->getPropagationDirection();
96 khat_.x() = pd[0];
97 khat_.y() = pd[1];
98 khat_.z() = pd[2];
99 khat_.normalize();
100 haveDirection = true;
101 }
102
103 if (lightParams->haveWavelength()) {
104 if (haveFrequency) {
105 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
106 "light: please specify one of: waveVector, wavelength, or"
107 "frequency (but only one of these).\n");
108 painCave.isFatal = 1;
109 simError();
110 }
111 // wavelengths are entered in nm to work with experimentalists.
112 // Convert to angstroms:
113 lambda_ = lightParams->getWavelength() * 10.0;
114 omega_ = 2.0 * Constants::PI * Constants::c / lambda_;
115 kmag_ = 2.0 * Constants::PI / lambda_;
116 haveFrequency = true;
117 }
118
119 if (lightParams->haveFrequency()) {
120 if (haveFrequency) {
121 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
122 "light: please specify one of: waveVector, wavelength, or"
123 "frequency (but only one of these).\n");
124 painCave.isFatal = 1;
125 simError();
126 }
127 // frequencies are entered in Hz to work with experimentalists.
128 // Convert to fs^-1
129 omega_ = lightParams->getFrequency() * 1.0e-15;
130 lambda_ = 2.0 * Constants::PI * Constants::c / omega_;
131 kmag_ = 2.0 * Constants::PI / lambda_;
132 haveFrequency = true;
133 }
134
135 if (haveFrequency && haveDirection) { k_ = khat_ * kmag_; }
136
137 if (lightParams->haveIntensity()) {
138 RealType intense = lightParams->getIntensity();
139 // intensities are input in W/cm^2
140 intense *= 1.439326e-11;
141 E0_ = std::sqrt(2.0 * intense / (Constants::epsilon0 * Constants::c));
142 // E0 now has units of kcal/mol e^-1 Angstroms^-1
143 haveE0 = true;
144 }
145
146 // Determine Polarization Type
147 jones_.clear();
148 jones_.reserve(2);
149 std::map<std::string, LightPolarization> stringToPolarization;
150
151 stringToPolarization["X"] = lightX;
152 stringToPolarization["Y"] = lightY;
153 stringToPolarization["+"] = lightPlus;
154 stringToPolarization["-"] = lightMinus;
155
156 if (lightParams->havePolarization()) {
157 std::string lpl = lightParams->getPolarization();
158 LightPolarization lp = stringToPolarization.find(lpl)->second;
159 switch (lp) {
160 case lightX:
161 jones_[0] = {1.0, 0.0};
162 jones_[1] = {0.0, 0.0};
163 havePolarization = true;
164 break;
165 case lightY:
166 jones_[0] = {0.0, 0.0};
167 jones_[1] = {1.0, 0.0};
168 havePolarization = true;
169 break;
170 case lightPlus:
171 jones_[0] = {1.0, 0.0};
172 jones_[1] = {0.0, 1.0};
173 havePolarization = true;
174 break;
175 case lightMinus:
176 jones_[0] = {1.0, 0.0};
177 jones_[1] = {0.0, -1.0};
178 havePolarization = true;
179 break;
180 default:
181 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
182 "Light: Unknown polarization type\n");
183 painCave.isFatal = 1;
184 painCave.severity = OPENMD_ERROR;
185 simError();
186 break;
187 }
188
189 } else {
190 std::string allowedPolarizations;
191 int currentLineLength = 0;
192
193 for (std::map<std::string, LightPolarization>::iterator polStrIter =
194 stringToPolarization.begin();
195 polStrIter != stringToPolarization.end(); ++polStrIter) {
196 allowedPolarizations += polStrIter->first + ", ";
197 currentLineLength += polStrIter->first.length() + 2;
198
199 if (currentLineLength >= 50) {
200 allowedPolarizations += "\n\t\t";
201 currentLineLength = 0;
202 }
203 }
204
205 allowedPolarizations.erase(allowedPolarizations.length() - 2, 2);
206
207 snprintf(
208 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
209 "Light: No polarization was set in the omd file. This parameter\n"
210 "\tmust be set to use Light, and can take any of these values:\n"
211 "\t\t%s.\n",
212 allowedPolarizations.c_str());
213 painCave.isFatal = 1;
214 painCave.severity = OPENMD_ERROR;
215 simError();
216 }
217
218 if (haveE0 && haveDirection && haveFrequency && havePolarization) {
219 doLight = true;
220 } else {
221 if (!haveDirection) {
222 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
223 "Light: could not determine direction of propagation.\n");
224 painCave.isFatal = 1;
225 simError();
226 }
227 if (!haveE0) {
228 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
229 "Light: intensity not specified.\n");
230 painCave.isFatal = 1;
231 simError();
232 }
233 if (!haveFrequency) {
234 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
235 "Light: could not determine frequency or wavelength.\n");
236 painCave.isFatal = 1;
237 simError();
238 }
239 if (!havePolarization) {
240 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
241 "Light: polarization not specifieid.\n");
242 painCave.isFatal = 1;
243 simError();
244 }
245 }
246
247 int storageLayout_ = info_->getSnapshotManager()->getAtomStorageLayout();
248 if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot = true;
249
250 // Relatively simple Euler angles between khat_ and lab frame:
251
252 RealType psi = 0.0;
253 RealType theta =
254 acos(std::min((RealType)1.0, std::max((RealType)-1.0, khat_[2])));
255 RealType phi = std::atan2(-khat_[1], khat_[0]);
256
257 if (phi < 0) phi += 2.0 * Constants::PI;
258
259 A_.setupRotMat(phi, theta, psi);
260 Ainv_ = A_.inverse();
261
262 initialized = true;
263 }
264
265 void Light::modifyForces() {
266 if (!initialized) initialize();
267
268 SimInfo::MoleculeIterator i;
269 Molecule::AtomIterator j;
270 Molecule* mol;
271 Atom* atom;
272 AtomType* atype;
273 potVec longRangePotential(0.0);
274 int l, m, n;
275 RealType C {}, U {}, fPot {};
276 Vector3d r {}, rp {}, v {}, f {}, trq {}, D {}, J {}, av {};
277 Vector3d EFk {}, EF {}, BF {};
278 Mat3x3d I {};
279
280 bool isCharge;
281
282 if (doLight) {
283 RealType t = info_->getSnapshotManager()->getCurrentSnapshot()->getTime();
284 U = 0.0;
285 fPot = 0.0;
286
287 for (mol = info_->beginMolecule(i); mol != NULL;
288 mol = info_->nextMolecule(i)) {
289 for (atom = mol->beginAtom(j); atom != NULL; atom = mol->nextAtom(j)) {
290 isCharge = false;
291 C = 0.0;
292
293 atype = atom->getAtomType();
294
295 // We are not wrapping coordinates for light interactions:
296 r = atom->getPos();
297 v = atom->getVel();
298
299 rp = A_ * r; // atom's position in frame of light propagation
300
301 // e^{ i (k*z - omega * t) } is the main oscillatory component:
302 std::complex<RealType> argument(0.0, kmag_ * rp.z() - omega_ * t);
303 std::complex<RealType> Ex = E0_ * jones_[0] * std::exp(argument);
304 std::complex<RealType> Ey = E0_ * jones_[1] * std::exp(argument);
305
306 EFk.x() = Ex.real();
307 EFk.y() = Ey.real();
308 EFk.z() = 0.0;
309
310 EF = Ainv_ * EFk; // electric field rotated back into lab coordinates
311
312 // The magnetic field (BF) is perpendicular to both electric field
313 // and light propagation direction:
314
315 BF = cross(EF, khat_) / Constants::c;
316
317 atom->addElectricField(EF);
318
320 if (fca.isFixedCharge()) {
321 isCharge = true;
322 C = fca.getCharge();
323 }
324
326 if (fqa.isFluctuatingCharge()) {
327 isCharge = true;
328 C += atom->getFlucQPos();
329 atom->addFlucQFrc(dot(r, EF));
330 }
331
332 if (isCharge) {
333 f = C * (EF + cross(v, BF));
334 atom->addFrc(f);
335 U = -dot(r, f);
336 if (doParticlePot) { atom->addParticlePot(U); }
337 fPot += U;
338 }
339
341 if (ma.isDipole()) {
342 D = atom->getDipole() * Constants::dipoleConvert;
343
344 trq += cross(D, EF) + cross(D, cross(v, BF));
345
346 atom->addTrq(trq);
347
348 J = atom->getJ();
349 I = atom->getI();
350 if (atom->isLinear()) {
351 l = atom->linearAxis();
352 m = (l + 1) % 3;
353 n = (l + 2) % 3;
354 av[l] = 0;
355 av[m] = J[m] / I(m, m);
356 av[n] = J[n] / I(n, n);
357 } else {
358 av = I.inverse() * J;
359 }
360
361 f = cross(cross(av, D), BF);
362 atom->addFrc(f);
363
364 U = -dot(D, EF);
365 if (doParticlePot) { atom->addParticlePot(U); }
366 fPot += U;
367 }
368 }
369 }
370
371#ifdef IS_MPI
372 MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE, MPI_SUM,
373 MPI_COMM_WORLD);
374#endif
375
376 Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot();
377 longRangePotential = snap->getLongRangePotentials();
378 longRangePotential[ELECTROSTATIC_FAMILY] += fPot;
379 snap->setLongRangePotentials(longRangePotential);
380 }
381 }
382} // namespace OpenMD::Perturbations
Rotating Electric Field perturbation.
virtual Mat3x3d getI()
Returns the inertia tensor of this stuntdouble.
Definition Atom.cpp:62
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
The Snapshot class is a repository storing dynamic data during a Simulation.
Definition Snapshot.hpp:147
void addFlucQFrc(RealType cfrc)
Adds charge force into the current charge force of this stuntDouble.
Vector3d getVel()
Returns the current velocity of this stuntDouble.
void addElectricField(const Vector3d &eField)
Adds electric field into the current electric field of this stuntDouble.
int linearAxis()
Returns the linear axis of the rigidbody, atom and directional atom will always return -1.
bool isLinear()
Tests the if this stuntDouble is a linear rigidbody.
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.
Vector3d getJ()
Returns the current angular momentum of this stuntDouble (body -fixed).
void addFrc(const Vector3d &frc)
Adds force into the current force of this stuntDouble.
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.
@ ELECTROSTATIC_FAMILY
Coulombic and point-multipole interactions.