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