OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
DipoleOrientation.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/staticProps/DipoleOrientation.hpp"
46
47#include <string>
48#include <vector>
49
50#include "applications/staticProps/SpatialStatistics.hpp"
51#include "brains/SimInfo.hpp"
53#include "math/Vector3.hpp"
55#include "utils/Accumulator.hpp"
56#include "utils/AccumulatorView.hpp"
57#include "utils/BaseAccumulator.hpp"
58#include "utils/StringUtils.hpp"
59
60using namespace OpenMD::Utils;
61
62namespace OpenMD {
63
64 DipoleOrientation::DipoleOrientation(
65 SimInfo* info, const std::string& filename, const std::string& sele,
66 const RealType dipoleX, const RealType dipoleY, const RealType dipoleZ,
67 int nzbins, int axis) :
68 SlabStatistics(info, filename, sele, nzbins, axis),
69 axis_(axis) {
70 switch (axis_) {
71 case 0:
72 axisLabel_ = "x";
73 refAxis_ = Vector3d(1, 0, 0);
74 break;
75 case 1:
76 axisLabel_ = "y";
77 refAxis_ = Vector3d(0, 1, 0);
78 break;
79 case 2:
80 default:
81 axisLabel_ = "z";
82 refAxis_ = Vector3d(0, 0, 1);
83 break;
84 }
85 setOutputName(getPrefix(filename) + ".Sz");
86
87 dipoleVector_ = Vector3d(dipoleX, dipoleY, dipoleZ);
88 dipoleVector_.normalize();
89
90 // Pre-load the OutputData
91 data_.resize(DipoleOrientation::ENDINDEX);
92
93 OutputData z;
94 z.units = "Angstroms";
95 z.title = axisLabel_;
96 z.dataHandling = DataHandling::Average;
97 for (unsigned int i = 0; i < nBins_; i++)
98 z.accumulator.push_back(
99 std::make_unique<AccumulatorView<RealAccumulator>>());
100 data_[Z] = std::move(z);
101
102 OutputData orderS;
103 orderS.units = "";
104 orderS.title = "Orientational Order parameter";
105 orderS.dataHandling = DataHandling::Average;
106 for (unsigned int i = 0; i < nBins_; i++)
107 orderS.accumulator.push_back(
108 std::make_unique<AccumulatorView<RealAccumulator>>());
109 data_[ORDERS] = std::move(orderS);
110
111 OutputData orderSCos;
112 orderSCos.units = "";
113 orderSCos.title = "Orientational Order parameter cosine Theta";
114 orderSCos.dataHandling = DataHandling::Average;
115 for (unsigned int i = 0; i < nBins_; i++)
116 orderSCos.accumulator.push_back(
117 std::make_unique<AccumulatorView<RealAccumulator>>());
118 data_[ORDERSCOS] = std::move(orderSCos);
119 }
120
121 void DipoleOrientation::processFrame(int) {
122 RealType z;
123
124 hmat_ = currentSnapshot_->getHmat();
125
126 for (unsigned int i = 0; i < nBins_; i++) {
127 z = (((RealType)i + 0.5) / (RealType)nBins_) * hmat_(axis_, axis_);
128 data_[Z].accumulator[i]->add(z);
129 }
130
131 volume_ = currentSnapshot_->getVolume();
132
133 StuntDouble* sd;
134 int i;
135
136 std::vector<RealType> binS(nBins_, 0.0);
137 std::vector<RealType> binSCos(nBins_, 0.0);
138
139 std::vector<int> count(nBins_, 0.0);
140
141 if (evaluator_.isDynamic()) {
142 seleMan_.setSelectionSet(evaluator_.evaluate());
143 }
144
145 // loop over the selected atoms:
147 Vector3d rotatedDipoleVector;
148 RealType ctheta;
149 RealType orderParameter;
150
151 for (sd = seleMan_.beginSelected(i); sd != NULL;
152 sd = seleMan_.nextSelected(i)) {
153 // figure out where that object is:
154 Vector3d pos = sd->getPos();
155
156 int bin = getBin(pos);
157
158 if (sd->isDirectional() || sd->isRigidBody()) {
159 rotMat = sd->getA();
160 rotatedDipoleVector = rotMat * dipoleVector_;
161 rotatedDipoleVector.normalize();
162 ctheta = dot(rotatedDipoleVector, refAxis_);
163
164 orderParameter = (3 * (ctheta * ctheta) - 1) / 2;
165
166 binS[bin] += orderParameter;
167 binSCos[bin] += ctheta;
168
169 count[bin] += 1;
170 }
171 }
172
173 for (unsigned int i = 0; i < nBins_; i++) {
174 count[i] != 0 ? data_[ORDERS].accumulator[i]->add(binS[i] / count[i]) :
175 data_[ORDERS].accumulator[i]->add(binS[i]);
176 count[i] != 0 ?
177 data_[ORDERSCOS].accumulator[i]->add(binSCos[i] / count[i]) :
178 data_[ORDERSCOS].accumulator[i]->add(binSCos[i]);
179 }
180 }
181} // namespace OpenMD
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.
Vector3d getPos()
Returns the current position of this stuntDouble.
bool isRigidBody()
Tests if this stuntDouble is a rigid body.
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)