OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
ChargeOrientationCorrFunc.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/ChargeOrientationCorrFunc.hpp"
46
49#include "types/FixedChargeAdapter.hpp"
50#include "types/FluctuatingChargeAdapter.hpp"
51
52namespace OpenMD {
53 ChargeOrientationCorrFunc::ChargeOrientationCorrFunc(
54 SimInfo* info, const std::string& filename, const std::string& sele1,
55 const std::string& sele2, const RealType dipoleX, const RealType dipoleY,
56 const RealType dipoleZ, const RealType cutOff, const int axis) :
57 ObjectCCF<RealType>(info, filename, sele1, sele2),
58 axis_(axis) {
59 setCorrFuncType(
60 "Charge - Orientation Order Parameter Cross Correlation Function");
61 setOutputName(getPrefix(dumpFilename_) + ".QScorr");
62
63 charges_.resize(nFrames_);
64 CosTheta_.resize(nFrames_);
65
66 sumCharge_ = 0;
67 sumCosTheta_ = 0;
68 chargeCount_ = 0;
69 CosThetaCount_ = 0;
70
71 dipoleVector_ = Vector3d(dipoleX, dipoleY, dipoleZ);
72 dipoleVector_.normalize();
73
74 switch (axis_) {
75 case 0:
76 axisLabel_ = "x";
77 refAxis_ = Vector3d(1, 0, 0);
78 break;
79 case 1:
80 axisLabel_ = "y";
81 refAxis_ = Vector3d(0, 1, 0);
82 break;
83 case 2:
84 default:
85 axisLabel_ = "z";
86 refAxis_ = Vector3d(0, 0, 1);
87 break;
88 }
89 }
90
91 void ChargeOrientationCorrFunc::validateSelection(SelectionManager& seleMan) {
92 StuntDouble* sd;
93 int i;
94
95 for (sd = seleMan.beginSelected(i); sd != NULL;
96 sd = seleMan.nextSelected(i)) {
97 Atom* atom = static_cast<Atom*>(sd);
98 AtomType* atomType = atom->getAtomType();
100
101 if (!sd->isDirectional() && !fqa.isFluctuatingCharge()) {
102 snprintf(
103 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
104 "ChargeOrientationCorrFunc::validateSelection Error: selection "
105 "%d (%s)\n"
106 "\t is not a Directional object\n",
107 sd->getGlobalIndex(), sd->getType().c_str());
108 painCave.isFatal = 1;
109 simError();
110 }
111 }
112 }
113
114 int ChargeOrientationCorrFunc::computeProperty1(int frame, StuntDouble* sd) {
115 RealType q = 0.0;
116 Atom* atom = static_cast<Atom*>(sd);
117
118 AtomType* atomType = atom->getAtomType();
119
121 if (fca.isFixedCharge()) { q += fca.getCharge(); }
122
124 if (fqa.isFluctuatingCharge()) { q += atom->getFlucQPos(); }
125
126 propertyTemp = q;
127 charges_[frame].push_back(propertyTemp);
128 sumCharge_ += propertyTemp;
129 chargeCount_++;
130 return charges_[frame].size() - 1;
131 }
132
133 int ChargeOrientationCorrFunc::computeProperty2(int frame, StuntDouble* sd) {
135 Vector3d rotatedDipoleVector;
136 RealType ctheta(0.0);
137
138 rotMat = sd->getA();
139 rotatedDipoleVector = rotMat * dipoleVector_;
140 rotatedDipoleVector.normalize();
141 ctheta = dot(rotatedDipoleVector, refAxis_);
142
143 propertyTemp = ctheta;
144 CosTheta_[frame].push_back(propertyTemp);
145 sumCosTheta_ += propertyTemp;
146 CosThetaCount_++;
147 return CosTheta_[frame].size() - 1;
148 }
149
150 RealType ChargeOrientationCorrFunc::calcCorrVal(int frame1, int frame2,
151 int id1, int id2) {
152 return charges_[frame1][id1] * CosTheta_[frame2][id2];
153 }
154
155 void ChargeOrientationCorrFunc::postCorrelate() {
156 // gets the average of the charges
157 sumCharge_ /= RealType(chargeCount_);
158
159 // gets the average of the CosTheta
160 sumCosTheta_ /= RealType(CosThetaCount_);
161
162 RealType correlationOfAverages_ = sumCharge_ * sumCosTheta_;
163 for (unsigned int i = 0; i < nTimeBins_; ++i) {
164 if (count_[i] > 0) {
165 histogram_[i] /= RealType(count_[i]);
166
167 // The correlation of the averages is subtracted
168 // from the correlation value:
169 histogram_[i] -= correlationOfAverages_;
170 } else {
171 histogram_[i] = 0;
172 }
173 }
174 }
175} // namespace OpenMD
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
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!"
RotMat3x3d getA()
Returns the current rotation matrix of this stuntDouble.
virtual std::string getType()=0
Returns the name of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
int getGlobalIndex()
Returns the global index of this stuntDouble.
bool isDirectional()
Tests if this stuntDouble is a directional one.
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 dot(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the dot product of two DynamicVectors.
std::string getPrefix(const std::string &str)