1 |
gezelter |
2015 |
/* |
2 |
|
|
* Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
|
|
* |
4 |
|
|
* The University of Notre Dame grants you ("Licensee") a |
5 |
|
|
* non-exclusive, royalty free, license to use, modify and |
6 |
|
|
* redistribute this software in source and binary code form, provided |
7 |
|
|
* that the following conditions are met: |
8 |
|
|
* |
9 |
|
|
* 1. Redistributions of source code must retain the above copyright |
10 |
|
|
* notice, this list of conditions and the following disclaimer. |
11 |
|
|
* |
12 |
|
|
* 2. Redistributions in binary form must reproduce the above copyright |
13 |
|
|
* notice, this list of conditions and the following disclaimer in the |
14 |
|
|
* documentation and/or other materials provided with the |
15 |
|
|
* distribution. |
16 |
|
|
* |
17 |
|
|
* This software is provided "AS IS," without a warranty of any |
18 |
|
|
* kind. All express or implied conditions, representations and |
19 |
|
|
* warranties, including any implied warranty of merchantability, |
20 |
|
|
* fitness for a particular purpose or non-infringement, are hereby |
21 |
|
|
* excluded. The University of Notre Dame and its licensors shall not |
22 |
|
|
* be liable for any damages suffered by licensee as a result of |
23 |
|
|
* using, modifying or distributing the software or its |
24 |
|
|
* derivatives. In no event will the University of Notre Dame or its |
25 |
|
|
* licensors be liable for any lost revenue, profit or data, or for |
26 |
|
|
* direct, indirect, special, consequential, incidental or punitive |
27 |
|
|
* damages, however caused and regardless of the theory of liability, |
28 |
|
|
* arising out of the use of or inability to use software, even if the |
29 |
|
|
* University of Notre Dame has been advised of the possibility of |
30 |
|
|
* such damages. |
31 |
|
|
* |
32 |
|
|
* SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
|
|
* research, please cite the appropriate papers when you publish your |
34 |
|
|
* work. Good starting points are: |
35 |
|
|
* |
36 |
|
|
* [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
|
|
* [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
|
|
* [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008). |
39 |
|
|
* [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010). |
40 |
|
|
* [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). |
41 |
|
|
* [6] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012). |
42 |
|
|
*/ |
43 |
|
|
|
44 |
|
|
#include "applications/staticProps/TetrahedralityParamXYZ.hpp" |
45 |
|
|
#include "utils/simError.h" |
46 |
|
|
#include "io/DumpReader.hpp" |
47 |
|
|
#include "primitives/Molecule.hpp" |
48 |
|
|
#include "utils/NumericConstant.hpp" |
49 |
|
|
#include <vector> |
50 |
|
|
#include <algorithm> |
51 |
|
|
#include <fstream> |
52 |
|
|
|
53 |
|
|
using namespace std; |
54 |
|
|
namespace OpenMD { |
55 |
|
|
TetrahedralityParamXYZ::TetrahedralityParamXYZ(SimInfo* info, |
56 |
|
|
const std::string& filename, |
57 |
|
|
const std::string& sele1, |
58 |
|
|
const std::string& sele2, |
59 |
|
|
RealType rCut, RealType voxelSize, |
60 |
|
|
RealType gaussWidth) |
61 |
|
|
: StaticAnalyser(info, filename), selectionScript1_(sele1), |
62 |
|
|
evaluator1_(info), seleMan1_(info), selectionScript2_(sele2), |
63 |
|
|
evaluator2_(info), seleMan2_(info), rCut_(rCut), voxelSize_(voxelSize), |
64 |
|
|
gaussWidth_(gaussWidth) { |
65 |
|
|
|
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 |
|
|
|
95 |
|
|
setOutputName(getPrefix(filename) + ".Qxyz"); |
96 |
|
|
} |
97 |
|
|
|
98 |
|
|
TetrahedralityParamXYZ::~TetrahedralityParamXYZ() { |
99 |
|
|
} |
100 |
|
|
|
101 |
|
|
void TetrahedralityParamXYZ::process() { |
102 |
|
|
Molecule* mol; |
103 |
|
|
StuntDouble* sd; |
104 |
|
|
StuntDouble* sd2; |
105 |
|
|
StuntDouble* sdi; |
106 |
|
|
StuntDouble* sdj; |
107 |
|
|
RigidBody* rb; |
108 |
|
|
int myIndex; |
109 |
|
|
SimInfo::MoleculeIterator mi; |
110 |
|
|
Molecule::RigidBodyIterator rbIter; |
111 |
|
|
Vector3d vec; |
112 |
|
|
Vector3d ri, rj, rk, rik, rkj; |
113 |
|
|
RealType r; |
114 |
|
|
RealType cospsi; |
115 |
|
|
RealType Qk; |
116 |
|
|
std::vector<std::pair<RealType,StuntDouble*> > myNeighbors; |
117 |
|
|
std::vector<std::pair<Vector3d, RealType> > qvals; |
118 |
|
|
std::vector<std::pair<Vector3d, RealType> >::iterator qiter; |
119 |
|
|
int isd1; |
120 |
|
|
int isd2; |
121 |
|
|
|
122 |
|
|
DumpReader reader(info_, dumpFilename_); |
123 |
|
|
int nFrames = reader.getNFrames(); |
124 |
|
|
|
125 |
|
|
for (int istep = 0; istep < nFrames; istep += step_) { |
126 |
|
|
reader.readFrame(istep); |
127 |
|
|
|
128 |
plouden |
2016 |
currentSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot(); |
129 |
|
|
Mat3x3d hmat = currentSnapshot_->getHmat(); |
130 |
|
|
Vector3d halfBox = Vector3d(hmat(0,0), hmat(1,1), hmat(2,2)) / 2.0; |
131 |
|
|
|
132 |
gezelter |
2015 |
if (evaluator1_.isDynamic()) { |
133 |
|
|
seleMan1_.setSelectionSet(evaluator1_.evaluate()); |
134 |
|
|
} |
135 |
|
|
|
136 |
|
|
if (evaluator2_.isDynamic()) { |
137 |
|
|
seleMan2_.setSelectionSet(evaluator2_.evaluate()); |
138 |
|
|
} |
139 |
|
|
|
140 |
|
|
// update the positions of atoms which belong to the rigidbodies |
141 |
|
|
for (mol = info_->beginMolecule(mi); mol != NULL; |
142 |
|
|
mol = info_->nextMolecule(mi)) { |
143 |
|
|
for (rb = mol->beginRigidBody(rbIter); rb != NULL; |
144 |
|
|
rb = mol->nextRigidBody(rbIter)) { |
145 |
|
|
rb->updateAtoms(); |
146 |
|
|
} |
147 |
|
|
} |
148 |
|
|
|
149 |
|
|
qvals.clear(); |
150 |
|
|
|
151 |
|
|
// outer loop is over the selected StuntDoubles: |
152 |
|
|
for (sd = seleMan1_.beginSelected(isd1); sd != NULL; |
153 |
|
|
sd = seleMan1_.nextSelected(isd1)) { |
154 |
|
|
|
155 |
|
|
myIndex = sd->getGlobalIndex(); |
156 |
|
|
|
157 |
|
|
Qk = 1.0; |
158 |
|
|
myNeighbors.clear(); |
159 |
|
|
|
160 |
|
|
for (sd2 = seleMan2_.beginSelected(isd2); sd2 != NULL; |
161 |
|
|
sd2 = seleMan2_.nextSelected(isd2)) { |
162 |
|
|
|
163 |
|
|
if (sd2->getGlobalIndex() != myIndex) { |
164 |
|
|
|
165 |
|
|
vec = sd->getPos() - sd2->getPos(); |
166 |
|
|
|
167 |
|
|
if (usePeriodicBoundaryConditions_) |
168 |
|
|
currentSnapshot_->wrapVector(vec); |
169 |
|
|
|
170 |
|
|
r = vec.length(); |
171 |
|
|
|
172 |
|
|
// Check to see if neighbor is in bond cutoff |
173 |
|
|
|
174 |
|
|
if (r < rCut_) { |
175 |
|
|
myNeighbors.push_back(std::make_pair(r,sd2)); |
176 |
|
|
} |
177 |
|
|
} |
178 |
|
|
} |
179 |
|
|
|
180 |
|
|
// Sort the vector using predicate and std::sort |
181 |
|
|
std::sort(myNeighbors.begin(), myNeighbors.end()); |
182 |
|
|
|
183 |
|
|
// Use only the 4 closest neighbors to do the rest of the work: |
184 |
|
|
|
185 |
|
|
int nbors = myNeighbors.size()> 4 ? 4 : myNeighbors.size(); |
186 |
|
|
int nang = int (0.5 * (nbors * (nbors - 1))); |
187 |
|
|
|
188 |
|
|
rk = sd->getPos(); |
189 |
|
|
|
190 |
|
|
for (int i = 0; i < nbors-1; i++) { |
191 |
|
|
|
192 |
|
|
sdi = myNeighbors[i].second; |
193 |
|
|
ri = sdi->getPos(); |
194 |
|
|
rik = rk - ri; |
195 |
|
|
if (usePeriodicBoundaryConditions_) |
196 |
|
|
currentSnapshot_->wrapVector(rik); |
197 |
|
|
|
198 |
|
|
rik.normalize(); |
199 |
|
|
|
200 |
|
|
for (int j = i+1; j < nbors; j++) { |
201 |
|
|
|
202 |
|
|
sdj = myNeighbors[j].second; |
203 |
|
|
rj = sdj->getPos(); |
204 |
|
|
rkj = rk - rj; |
205 |
|
|
if (usePeriodicBoundaryConditions_) |
206 |
|
|
currentSnapshot_->wrapVector(rkj); |
207 |
|
|
rkj.normalize(); |
208 |
|
|
|
209 |
|
|
cospsi = dot(rik,rkj); |
210 |
|
|
|
211 |
|
|
// Calculates scaled Qk for each molecule using calculated |
212 |
|
|
// angles from 4 or fewer nearest neighbors. |
213 |
|
|
Qk -= (pow(cospsi + 1.0 / 3.0, 2) * 2.25 / nang); |
214 |
|
|
} |
215 |
|
|
} |
216 |
|
|
|
217 |
|
|
if (nang > 0) { |
218 |
|
|
if (usePeriodicBoundaryConditions_) |
219 |
|
|
currentSnapshot_->wrapVector(rk); |
220 |
|
|
qvals.push_back(std::make_pair(rk, Qk)); |
221 |
|
|
} |
222 |
|
|
} |
223 |
|
|
|
224 |
|
|
for (int i = 0; i < nBins_(0); ++i) { |
225 |
|
|
for(int j = 0; j < nBins_(1); ++j) { |
226 |
|
|
for(int k = 0; k < nBins_(2); ++k) { |
227 |
plouden |
2016 |
Vector3d pos = Vector3d(i, j, k) * voxelSize_ - halfBox; |
228 |
gezelter |
2015 |
for(qiter = qvals.begin(); qiter != qvals.end(); ++qiter) { |
229 |
|
|
Vector3d d = pos - (*qiter).first; |
230 |
plouden |
2016 |
currentSnapshot_->wrapVector(d); |
231 |
gezelter |
2015 |
RealType denom = pow(2.0 * sqrt(M_PI) * gaussWidth_, 3); |
232 |
|
|
RealType exponent = -dot(d,d) / pow(2.0*gaussWidth_, 2); |
233 |
|
|
RealType weight = exp(exponent) / denom; |
234 |
|
|
count_[i][j][k] += weight; |
235 |
|
|
hist_[i][j][k] += weight * (*qiter).second; |
236 |
|
|
} |
237 |
|
|
} |
238 |
|
|
} |
239 |
|
|
} |
240 |
|
|
} |
241 |
|
|
writeQxyz(); |
242 |
|
|
} |
243 |
|
|
|
244 |
|
|
void TetrahedralityParamXYZ::writeQxyz() { |
245 |
plouden |
2016 |
|
246 |
|
|
Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat(); |
247 |
|
|
|
248 |
gezelter |
2015 |
// normalize by total weight in voxel: |
249 |
|
|
for (unsigned int i = 0; i < hist_.size(); ++i) { |
250 |
|
|
for(unsigned int j = 0; j < hist_[i].size(); ++j) { |
251 |
|
|
for(unsigned int k = 0;k < hist_[i][j].size(); ++k) { |
252 |
plouden |
2016 |
hist_[i][j][k] = hist_[i][j][k] / count_[i][j][k]; |
253 |
gezelter |
2015 |
} |
254 |
|
|
} |
255 |
|
|
} |
256 |
|
|
|
257 |
|
|
std::ofstream qXYZstream(outputFilename_.c_str()); |
258 |
|
|
if (qXYZstream.is_open()) { |
259 |
plouden |
2016 |
qXYZstream << "# AmiraMesh ASCII 1.0\n\n"; |
260 |
|
|
qXYZstream << "# Dimensions in x-, y-, and z-direction\n"; |
261 |
|
|
qXYZstream << " define Lattice " << hist_.size() << " " << hist_[0].size() << " " << hist_[0][0].size() << "\n"; |
262 |
|
|
|
263 |
|
|
qXYZstream << "Parameters {\n"; |
264 |
|
|
qXYZstream << " CoordType \"uniform\",\n"; |
265 |
|
|
qXYZstream << " # BoundingBox is xmin xmax ymin ymax zmin zmax\n"; |
266 |
|
|
qXYZstream << " BoundingBox 0.0 " << hmat(0,0) << |
267 |
|
|
" 0.0 " << hmat(1,1) << |
268 |
|
|
" 0.0 " << hmat(2,2) << "\n"; |
269 |
|
|
qXYZstream << "}\n"; |
270 |
|
|
|
271 |
|
|
qXYZstream << "Lattice { double ScalarField } = @1\n"; |
272 |
|
|
|
273 |
|
|
qXYZstream << "@1\n"; |
274 |
gezelter |
2015 |
|
275 |
plouden |
2016 |
int xsize = hist_.size(); |
276 |
|
|
int ysize = hist_[0].size(); |
277 |
|
|
int zsize = hist_[0][0].size(); |
278 |
|
|
|
279 |
|
|
for (unsigned int k = 0; k < zsize; ++k) { |
280 |
|
|
for(unsigned int j = 0; j < ysize; ++j) { |
281 |
|
|
for(unsigned int i = 0; i < xsize; ++i) { |
282 |
|
|
qXYZstream << hist_[i][j][k] << " "; |
283 |
|
|
|
284 |
|
|
//qXYZstream.write(reinterpret_cast<char *>( &hist_[i][j][k] ), |
285 |
|
|
// sizeof( hist_[i][j][k] )); |
286 |
gezelter |
2015 |
} |
287 |
|
|
} |
288 |
|
|
} |
289 |
|
|
|
290 |
|
|
} else { |
291 |
|
|
sprintf(painCave.errMsg, "TetrahedralityParamXYZ: unable to open %s\n", |
292 |
|
|
outputFilename_.c_str()); |
293 |
|
|
painCave.isFatal = 1; |
294 |
|
|
simError(); |
295 |
|
|
} |
296 |
|
|
qXYZstream.close(); |
297 |
|
|
} |
298 |
|
|
} |
299 |
|
|
|
300 |
|
|
|
301 |
|
|
|