OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ChargeKineticCorrFunc.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 "applications/dynamicProps/ChargeKineticCorrFunc.hpp"
46
48#include "types/FixedChargeAdapter.hpp"
49#include "types/FluctuatingChargeAdapter.hpp"
50
51namespace OpenMD {
52 ChargeKineticCorrFunc::ChargeKineticCorrFunc(SimInfo* info,
53 const std::string& filename,
54 const std::string& sele1,
55 const std::string& sele2,
56 const RealType cutoff) :
57 ObjectCCF<RealType>(info, filename, sele1, sele2),
58 rcut_(cutoff) {
59 setCorrFuncType("Charge - Kinetic Cross Correlation Function");
60 setOutputName(getPrefix(dumpFilename_) + ".QKcorr");
61
62 charges_.resize(nFrames_);
63 kinetic_.resize(nFrames_);
64
65 sumCharge_ = 0;
66 sumKinetic_ = 0;
67 chargeCount_ = 0;
68 kineticCount_ = 0;
69 }
70
71 void ChargeKineticCorrFunc::validateSelection(SelectionManager& seleMan) {
72 StuntDouble* sd;
73 int i;
74
75 for (sd = seleMan.beginSelected(i); sd != NULL;
76 sd = seleMan.nextSelected(i)) {
77 Atom* atom = static_cast<Atom*>(sd);
78 AtomType* atomType = atom->getAtomType();
80
81 if (!sd->isDirectional() && !fqa.isFluctuatingCharge()) {
82 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
83 "ChargeKineticCorrFunc::validateSelection Error: selection "
84 "%d (%s)\n"
85 "\t is not a Directional object\n",
86 sd->getGlobalIndex(), sd->getType().c_str());
87 painCave.isFatal = 1;
88 simError();
89 }
90 }
91 }
92
93 int ChargeKineticCorrFunc::computeProperty1(int frame, StuntDouble* sd) {
94 RealType q = 0.0;
95 pos1_ = sd->getPos();
96 Atom* atom = static_cast<Atom*>(sd);
97
98 AtomType* atomType = atom->getAtomType();
99
101 if (fca.isFixedCharge()) { q += fca.getCharge(); }
102
104 if (fqa.isFluctuatingCharge()) { q += atom->getFlucQPos(); }
105
106 propertyTemp = q;
107 charges_[frame].push_back(propertyTemp);
108 sumCharge_ += propertyTemp;
109 chargeCount_++;
110 return charges_[frame].size() - 1;
111 }
112
113 int ChargeKineticCorrFunc::computeProperty2(int frame, StuntDouble* sd) {
114 RealType kinetic(0.0);
115 pos2_ = sd->getPos();
116
117 RealType mass = sd->getMass();
118 Vector3d vel = sd->getVel();
119
120 kinetic += mass * (vel[0] * vel[0] + vel[1] * vel[1] + vel[2] * vel[2]);
121
122 Vector3d angMom;
123 Mat3x3d I;
124 int i, j, k;
125
126 if (sd->isDirectional()) {
127 angMom = sd->getJ();
128 I = sd->getI();
129
130 if (sd->isLinear()) {
131 i = sd->linearAxis();
132 j = (i + 1) % 3;
133 k = (i + 2) % 3;
134 kinetic +=
135 angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k);
136 } else {
137 kinetic += angMom[0] * angMom[0] / I(0, 0) +
138 angMom[1] * angMom[1] / I(1, 1) +
139 angMom[2] * angMom[2] / I(2, 2);
140 }
141 }
142 propertyTemp = 0.5 * kinetic;
143 kinetic_[frame].push_back(propertyTemp);
144 sumKinetic_ += propertyTemp;
145 kineticCount_++;
146 return kinetic_[frame].size() - 1;
147 }
148
149 RealType ChargeKineticCorrFunc::calcCorrVal(int frame1, int frame2, int id1,
150 int id2) {
151 Vector3d sep = pos1_ - pos2_;
152 RealType distance = sep.length();
153 if (distance <= rcut_)
154 return 0;
155 else
156 return charges_[frame1][id1] * kinetic_[frame2][id2];
157 }
158
159 void ChargeKineticCorrFunc::postCorrelate() {
160 // gets the average of the charges
161 sumCharge_ /= RealType(chargeCount_);
162
163 // gets the average of the kinetic energy
164 sumKinetic_ /= RealType(kineticCount_);
165
166 RealType correlationOfAverages_ = sumCharge_ * sumKinetic_;
167 for (unsigned int i = 0; i < nTimeBins_; ++i) {
168 if (count_[i] > 0) {
169 histogram_[i] /= RealType(count_[i]);
170
171 // The correlation of the averages is subtracted
172 // from the correlation value:
173 histogram_[i] -= correlationOfAverages_;
174 } else {
175 histogram_[i] = 0;
176 }
177 }
178 }
179} // namespace OpenMD
AtomType is what OpenMD looks to for unchanging data about an atom.
Definition AtomType.hpp:66
StuntDouble * nextSelected(int &i)
Finds the next selected StuntDouble in the selection.
StuntDouble * beginSelected(int &i)
Finds the first selected StuntDouble in the selection.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getVel()
Returns the current velocity of this stuntDouble.
int linearAxis()
Returns the linear axis of the rigidbody, atom and directional atom will always return -1.
RealType getMass()
Returns the mass of this stuntDouble.
virtual Mat3x3d getI()=0
Returns the inertia tensor of this stuntDouble.
virtual std::string getType()=0
Returns the name of this stuntDouble.
bool isLinear()
Tests the if this stuntDouble is a linear rigidbody.
Vector3d getPos()
Returns the current position of this stuntDouble.
Vector3d getJ()
Returns the current angular momentum of this stuntDouble (body -fixed).
int getGlobalIndex()
Returns the global index of this stuntDouble.
bool isDirectional()
Tests if this stuntDouble is a directional one.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
std::string getPrefix(const std::string &str)
Real distance(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the distance between two DynamicVectors.