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