OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
TetrahedralityParamXYZ.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/TetrahedralityParamXYZ.hpp"
46
47#include <algorithm>
48#include <fstream>
49#include <vector>
50
51#include "io/DumpReader.hpp"
53#include "utils/Constants.hpp"
54#include "utils/simError.h"
55
56using namespace std;
57namespace OpenMD {
58 TetrahedralityParamXYZ::TetrahedralityParamXYZ(
59 SimInfo* info, const std::string& filename, const std::string& sele1,
60 const std::string& sele2, RealType rCut, RealType voxelSize,
61 RealType gaussWidth) :
62 StaticAnalyser(info, filename, 1),
63 selectionScript1_(sele1), selectionScript2_(sele2), seleMan1_(info),
64 seleMan2_(info), evaluator1_(info), evaluator2_(info), rCut_(rCut),
65 voxelSize_(voxelSize), gaussWidth_(gaussWidth) {
66 evaluator1_.loadScriptString(sele1);
67 if (!evaluator1_.isDynamic()) {
68 seleMan1_.setSelectionSet(evaluator1_.evaluate());
69 }
70 evaluator2_.loadScriptString(sele2);
71 if (!evaluator2_.isDynamic()) {
72 seleMan2_.setSelectionSet(evaluator2_.evaluate());
73 }
74
75 Mat3x3d hmat = info->getSnapshotManager()->getCurrentSnapshot()->getHmat();
76
77 nBins_(0) = int(hmat(0, 0) / voxelSize);
78 nBins_(1) = int(hmat(1, 1) / voxelSize);
79 nBins_(2) = int(hmat(2, 2) / voxelSize);
80
81 hist_.resize(nBins_(0));
82 count_.resize(nBins_(0));
83 for (int i = 0; i < nBins_(0); ++i) {
84 hist_[i].resize(nBins_(1));
85 count_[i].resize(nBins_(1));
86 for (int j = 0; j < nBins_(1); ++j) {
87 hist_[i][j].resize(nBins_(2));
88 count_[i][j].resize(nBins_(2));
89 std::fill(hist_[i][j].begin(), hist_[i][j].end(), 0.0);
90 std::fill(count_[i][j].begin(), count_[i][j].end(), 0.0);
91 }
92 }
93
94 setOutputName(getPrefix(filename) + ".Qxyz");
95 }
96
97 void TetrahedralityParamXYZ::process() {
98 StuntDouble* sd;
99 StuntDouble* sd2;
100 StuntDouble* sdi;
101 StuntDouble* sdj;
102 int myIndex;
103 Vector3d vec;
104 Vector3d ri, rj, rk, rik, rkj;
105 RealType r;
106 RealType cospsi;
107 RealType Qk;
108 std::vector<std::pair<RealType, StuntDouble*>> myNeighbors;
109 // std::vector<std::pair<Vector3d, RealType> > qvals;
110 // std::vector<std::pair<Vector3d, RealType> >::iterator qiter;
111 int isd1;
112 int isd2;
113 bool usePeriodicBoundaryConditions_ =
114 info_->getSimParams()->getUsePeriodicBoundaryConditions();
115
116 int kMax = int(5.0 * gaussWidth_ / voxelSize_);
117 int kSqLim = kMax * kMax;
118 cerr << "gw = " << gaussWidth_ << " vS = " << voxelSize_
119 << " kMax = " << kMax << " kSqLim = " << kSqLim << "\n";
120
121 DumpReader reader(info_, dumpFilename_);
122 int nFrames = reader.getNFrames();
123
124 for (int istep = 0; istep < nFrames; istep += step_) {
125 reader.readFrame(istep);
126
127 currentSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot();
128 Mat3x3d hmat = currentSnapshot_->getHmat();
129 Vector3d halfBox = Vector3d(hmat(0, 0), hmat(1, 1), hmat(2, 2)) / 2.0;
130
131 if (evaluator1_.isDynamic()) {
132 seleMan1_.setSelectionSet(evaluator1_.evaluate());
133 }
134
135 if (evaluator2_.isDynamic()) {
136 seleMan2_.setSelectionSet(evaluator2_.evaluate());
137 }
138
139 // qvals.clear();
140
141 // outer loop is over the selected StuntDoubles:
142 for (sd = seleMan1_.beginSelected(isd1); sd != NULL;
143 sd = seleMan1_.nextSelected(isd1)) {
144 myIndex = sd->getGlobalIndex();
145
146 Qk = 1.0;
147 myNeighbors.clear();
148
149 for (sd2 = seleMan2_.beginSelected(isd2); sd2 != NULL;
150 sd2 = seleMan2_.nextSelected(isd2)) {
151 if (sd2->getGlobalIndex() != myIndex) {
152 vec = sd->getPos() - sd2->getPos();
153
154 if (usePeriodicBoundaryConditions_)
155 currentSnapshot_->wrapVector(vec);
156
157 r = vec.length();
158
159 // Check to see if neighbor is in bond cutoff
160
161 if (r < rCut_) { myNeighbors.push_back(std::make_pair(r, sd2)); }
162 }
163 }
164
165 // Sort the vector using predicate and std::sort
166 std::sort(myNeighbors.begin(), myNeighbors.end());
167
168 // Use only the 4 closest neighbors to do the rest of the work:
169
170 int nbors = myNeighbors.size() > 4 ? 4 : myNeighbors.size();
171 int nang = int(0.5 * (nbors * (nbors - 1)));
172
173 rk = sd->getPos();
174
175 for (int i = 0; i < nbors - 1; i++) {
176 sdi = myNeighbors[i].second;
177 ri = sdi->getPos();
178 rik = rk - ri;
179 if (usePeriodicBoundaryConditions_) currentSnapshot_->wrapVector(rik);
180
181 rik.normalize();
182
183 for (int j = i + 1; j < nbors; j++) {
184 sdj = myNeighbors[j].second;
185 rj = sdj->getPos();
186 rkj = rk - rj;
187 if (usePeriodicBoundaryConditions_)
188 currentSnapshot_->wrapVector(rkj);
189 rkj.normalize();
190
191 cospsi = dot(rik, rkj);
192
193 // Calculates scaled Qk for each molecule using calculated
194 // angles from 4 or fewer nearest neighbors.
195 Qk -= (pow(cospsi + 1.0 / 3.0, 2) * 2.25 / nang);
196 }
197 }
198
199 if (nang > 0) {
200 if (usePeriodicBoundaryConditions_) currentSnapshot_->wrapVector(rk);
201 // qvals.push_back(std::make_pair(rk, Qk));
202
203 Vector3d pos = rk + halfBox;
204
205 Vector3i whichVoxel(int(pos[0] / voxelSize_),
206 int(pos[1] / voxelSize_),
207 int(pos[2] / voxelSize_));
208
209 for (int l = -kMax; l <= kMax; l++) {
210 for (int m = -kMax; m <= kMax; m++) {
211 for (int n = -kMax; n <= kMax; n++) {
212 int kk = l * l + m * m + n * n;
213 if (kk <= kSqLim) {
214 int ll = (whichVoxel[0] + l) % nBins_(0);
215 ll = ll < 0 ? nBins_(0) + ll : ll;
216 int mm = (whichVoxel[1] + m) % nBins_(1);
217 mm = mm < 0 ? nBins_(1) + mm : mm;
218 int nn = (whichVoxel[2] + n) % nBins_(2);
219 nn = nn < 0 ? nBins_(2) + nn : nn;
220
221 Vector3d bPos = Vector3d(ll, mm, nn) * voxelSize_ - halfBox;
222 Vector3d d = bPos - rk;
223 currentSnapshot_->wrapVector(d);
224 RealType denom =
225 pow(2.0 * sqrt(Constants::PI) * gaussWidth_, 3);
226 RealType exponent = -dot(d, d) / pow(2.0 * gaussWidth_, 2);
227 RealType weight = exp(exponent) / denom;
228 count_[ll][mm][nn] += weight;
229 hist_[ll][mm][nn] += weight * Qk;
230 }
231 }
232 }
233 }
234 }
235 }
236
237 // for (int i = 0; i < nBins_(0); ++i) {
238 // for(int j = 0; j < nBins_(1); ++j) {
239 // for(int k = 0; k < nBins_(2); ++k) {
240 // Vector3d pos = Vector3d(i, j, k) * voxelSize_ - halfBox;
241 // for(qiter = qvals.begin(); qiter != qvals.end(); ++qiter) {
242 // Vector3d d = pos - (*qiter).first;
243 // currentSnapshot_->wrapVector(d);
244 // RealType denom = pow(2.0 * sqrt(Constants::PI) * gaussWidth_,
245 // 3); RealType exponent = -dot(d,d) / pow(2.0*gaussWidth_, 2);
246 // RealType weight = exp(exponent) / denom;
247 // count_[i][j][k] += weight;
248 // hist_[i][j][k] += weight * (*qiter).second;
249 // }
250 // }
251 // }
252 // }
253 }
254 // writeQxyz();
255 writeVTKGrid();
256 }
257
258 void TetrahedralityParamXYZ::writeVTKGrid() {
259 Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat();
260
261 // normalize by total weight in voxel:
262 for (unsigned int i = 0; i < hist_.size(); ++i) {
263 for (unsigned int j = 0; j < hist_[i].size(); ++j) {
264 for (unsigned int k = 0; k < hist_[i][j].size(); ++k) {
265 hist_[i][j][k] = hist_[i][j][k] / count_[i][j][k];
266 }
267 }
268 }
269
270 std::ofstream qXYZstream(outputFilename_.c_str());
271 if (qXYZstream.is_open()) {
272 qXYZstream << "# vtk DataFile Version 2.0\n";
273 qXYZstream << "Tetrahedrality Parameter volume rendering\n";
274 qXYZstream << "ASCII\n";
275 qXYZstream << "DATASET RECTILINEAR_GRID\n";
276 qXYZstream << "DIMENSIONS " << hist_.size() << " " << hist_[0].size()
277 << " " << hist_[0][0].size() << "\n";
278 qXYZstream << "X_COORDINATES " << hist_.size() << " float\n";
279 for (std::size_t i = 0; i < hist_.size(); ++i) {
280 qXYZstream << (RealType(i) / nBins_.x()) * hmat(0, 0) << " ";
281 }
282 qXYZstream << "\n";
283 qXYZstream << "Y_COORDINATES " << hist_[0].size() << " float\n";
284 for (std::size_t j = 0; j < hist_[0].size(); ++j) {
285 qXYZstream << (RealType(j) / nBins_.y()) * hmat(1, 1) << " ";
286 }
287 qXYZstream << "\n";
288 qXYZstream << "Z_COORDINATES " << hist_[0][0].size() << " float\n";
289 for (std::size_t k = 0; k < hist_[0][0].size(); ++k) {
290 qXYZstream << (RealType(k) / nBins_.z()) * hmat(2, 2) << " ";
291 }
292 qXYZstream << "\n";
293 qXYZstream << "POINT_DATA "
294 << hist_.size() * hist_[0].size() * hist_[0][0].size() << "\n";
295 qXYZstream << "SCALARS scalars float\n";
296 qXYZstream << "LOOKUP_TABLE default\n";
297
298 for (std::size_t k = 0; k < hist_[0][0].size(); ++k) {
299 for (std::size_t j = 0; j < hist_[0].size(); ++j) {
300 for (std::size_t i = 0; i < hist_.size(); ++i) {
301 qXYZstream << hist_[i][j][k] << " ";
302
303 // qXYZstream.write(reinterpret_cast<char *>( &hist_[i][j][k] ),
304 // sizeof( hist_[i][j][k] ));
305 }
306 }
307 }
308 qXYZstream << "\n";
309 } else {
310 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
311 "TetrahedralityParamXYZ: unable to open %s\n",
312 outputFilename_.c_str());
313 painCave.isFatal = 1;
314 simError();
315 }
316 qXYZstream.close();
317 }
318
319 void TetrahedralityParamXYZ::writeQxyz() {
320 Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat();
321
322 // normalize by total weight in voxel:
323 for (unsigned int i = 0; i < hist_.size(); ++i) {
324 for (unsigned int j = 0; j < hist_[i].size(); ++j) {
325 for (unsigned int k = 0; k < hist_[i][j].size(); ++k) {
326 hist_[i][j][k] = hist_[i][j][k] / count_[i][j][k];
327 }
328 }
329 }
330
331 std::ofstream qXYZstream(outputFilename_.c_str());
332 if (qXYZstream.is_open()) {
333 qXYZstream << "# AmiraMesh ASCII 1.0\n\n";
334 qXYZstream << "# Dimensions in x-, y-, and z-direction\n";
335 qXYZstream << " define Lattice " << hist_.size() << " "
336 << hist_[0].size() << " " << hist_[0][0].size() << "\n";
337
338 qXYZstream << "Parameters {\n";
339 qXYZstream << " CoordType \"uniform\",\n";
340 qXYZstream << " # BoundingBox is xmin xmax ymin ymax zmin zmax\n";
341 qXYZstream << " BoundingBox 0.0 " << hmat(0, 0) << " 0.0 "
342 << hmat(1, 1) << " 0.0 " << hmat(2, 2) << "\n";
343 qXYZstream << "}\n";
344
345 qXYZstream << "Lattice { double ScalarField } = @1\n";
346
347 qXYZstream << "@1\n";
348
349 for (std::size_t k = 0; k < hist_[0][0].size(); ++k) {
350 for (std::size_t j = 0; j < hist_[0].size(); ++j) {
351 for (std::size_t i = 0; i < hist_.size(); ++i) {
352 qXYZstream << hist_[i][j][k] << " ";
353
354 // qXYZstream.write(reinterpret_cast<char *>( &hist_[i][j][k] ),
355 // sizeof( hist_[i][j][k] ));
356 }
357 }
358 }
359
360 } else {
361 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
362 "TetrahedralityParamXYZ: unable to open %s\n",
363 outputFilename_.c_str());
364 painCave.isFatal = 1;
365 simError();
366 }
367 qXYZstream.close();
368 }
369} // namespace OpenMD
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)