OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
UniformField.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
61 UniformField::UniformField(SimInfo* info) :
62 ForceModifier {info}, initialized {false}, doUniformField {false},
63 doParticlePot {false} {
64 simParams = info_->getSimParams();
65 }
66
67 void UniformField::initialize() {
68 std::vector<RealType> ef;
69
70 if (simParams->haveElectricField()) {
71 doUniformField = true;
72 ef = simParams->getElectricField();
73 }
74 if (simParams->haveUniformField()) {
75 doUniformField = true;
76 ef = simParams->getUniformField();
77 }
78 if (ef.size() != 3) {
79 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
80 "UniformField: Incorrect number of parameters specified.\n"
81 "\tthere should be 3 parameters, but %lu were specified.\n",
82 ef.size());
83 painCave.isFatal = 1;
84 simError();
85 }
86 EF.x() = ef[0];
87 EF.y() = ef[1];
88 EF.z() = ef[2];
89
90 int storageLayout_ = info_->getSnapshotManager()->getAtomStorageLayout();
91 if (storageLayout_ & DataStorage::dslParticlePot) doParticlePot = true;
92 initialized = true;
93 }
94
95 void UniformField::modifyForces() {
96 if (!initialized) initialize();
97
98 SimInfo::MoleculeIterator i;
99 Molecule::AtomIterator j;
100 Molecule* mol;
101 Atom* atom;
102 AtomType* atype;
103 potVec longRangePotential(0.0);
104
105 RealType C;
106 Vector3d D;
107 RealType U;
108 RealType fPot;
109 Vector3d t;
110 Vector3d f;
111 Vector3d r;
112
113 bool isCharge;
114
115 if (doUniformField) {
116 U = 0.0;
117 fPot = 0.0;
118
119 for (mol = info_->beginMolecule(i); mol != NULL;
120 mol = info_->nextMolecule(i)) {
121 for (atom = mol->beginAtom(j); atom != NULL; atom = mol->nextAtom(j)) {
122 isCharge = false;
123 C = 0.0;
124
125 atype = atom->getAtomType();
126
127 // ad-hoc choice of the origin for potential calculation and
128 // fluctuating charge force:
129
130 r = atom->getPos();
132
133 atom->addElectricField(EF * Constants::chargeFieldConvert);
134
136 if (fca.isFixedCharge()) {
137 isCharge = true;
138 C = fca.getCharge();
139 }
140
142 if (fqa.isFluctuatingCharge()) {
143 isCharge = true;
144 C += atom->getFlucQPos();
145 atom->addFlucQFrc(dot(r, EF) * Constants::chargeFieldConvert);
146 }
147
148 if (isCharge) {
149 f = EF * C * Constants::chargeFieldConvert;
150 atom->addFrc(f);
151 U = -dot(r, f);
152
153 if (doParticlePot) { atom->addParticlePot(U); }
154 fPot += U;
155 }
156
158 if (ma.isDipole()) {
159 D = atom->getDipole() * Constants::dipoleFieldConvert;
160
161 t = cross(D, EF);
162 atom->addTrq(t);
163
164 U = -dot(D, EF);
165
166 if (doParticlePot) { atom->addParticlePot(U); }
167 fPot += U;
168 }
169 }
170 }
171
172#ifdef IS_MPI
173 MPI_Allreduce(MPI_IN_PLACE, &fPot, 1, MPI_REALTYPE, MPI_SUM,
174 MPI_COMM_WORLD);
175#endif
176
178 longRangePotential = snap->getLongRangePotentials();
179 longRangePotential[ELECTROSTATIC_FAMILY] += fPot;
180 snap->setLongRangePotentials(longRangePotential);
181 }
182 }
183} // namespace OpenMD
Uniform Electric Field perturbation.
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.
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
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
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.