ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/applications/staticProps/RNEMDStats.cpp
Revision: 1953
Committed: Thu Dec 5 18:19:26 2013 UTC (11 years, 5 months ago) by gezelter
File size: 16338 byte(s)
Log Message:
Rewrote much of selection module, added a bond correlation function

File Contents

# User Rev Content
1 gezelter 1865 /*
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     */
41    
42    
43     #include <algorithm>
44     #include <fstream>
45     #include "applications/staticProps/RNEMDStats.hpp"
46 gezelter 1881 #include "primitives/Molecule.hpp"
47 gezelter 1865 #include "utils/PhysicalConstants.hpp"
48    
49     namespace OpenMD {
50    
51     RNEMDZ::RNEMDZ(SimInfo* info, const std::string& filename,
52     const std::string& sele, int nzbins)
53     : SlabStatistics(info, filename, sele, nzbins) {
54    
55     setOutputName(getPrefix(filename) + ".rnemdZ");
56    
57     temperature = new OutputData;
58     temperature->units = "K";
59     temperature->title = "Temperature";
60     temperature->dataType = odtReal;
61     temperature->dataHandling = odhAverage;
62     temperature->accumulator.reserve(nBins_);
63     for (int i = 0; i < nBins_; i++)
64     temperature->accumulator.push_back( new Accumulator() );
65     data_.push_back(temperature);
66    
67     velocity = new OutputData;
68     velocity->units = "angstroms/fs";
69     velocity->title = "Velocity";
70     velocity->dataType = odtVector3;
71     velocity->dataHandling = odhAverage;
72     velocity->accumulator.reserve(nBins_);
73     for (int i = 0; i < nBins_; i++)
74     velocity->accumulator.push_back( new VectorAccumulator() );
75     data_.push_back(velocity);
76    
77     density = new OutputData;
78     density->units = "g cm^-3";
79     density->title = "Density";
80     density->dataType = odtReal;
81     density->dataHandling = odhAverage;
82     density->accumulator.reserve(nBins_);
83     for (int i = 0; i < nBins_; i++)
84     density->accumulator.push_back( new Accumulator() );
85     data_.push_back(density);
86     }
87    
88 gezelter 1881 void RNEMDZ::processFrame(int istep) {
89 gezelter 1884 RealType z;
90    
91     hmat_ = currentSnapshot_->getHmat();
92     for (int i = 0; i < nBins_; i++) {
93     z = (((RealType)i + 0.5) / (RealType)nBins_) * hmat_(2,2);
94     dynamic_cast<Accumulator*>(z_->accumulator[i])->add(z);
95     }
96     volume_ = currentSnapshot_->getVolume();
97    
98    
99 gezelter 1881 Molecule* mol;
100     RigidBody* rb;
101     StuntDouble* sd;
102     SimInfo::MoleculeIterator mi;
103     Molecule::RigidBodyIterator rbIter;
104     int i;
105 gezelter 1865
106 gezelter 1881 vector<RealType> binMass(nBins_, 0.0);
107 gezelter 1944 vector<Vector3d> binP(nBins_, V3Zero);
108 gezelter 1881 vector<RealType> binKE(nBins_, 0.0);
109 gezelter 1883 vector<unsigned int> binDof(nBins_, 0);
110 gezelter 1881
111     for (mol = info_->beginMolecule(mi); mol != NULL;
112     mol = info_->nextMolecule(mi)) {
113    
114     // change the positions of atoms which belong to the rigidbodies
115    
116     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
117     rb = mol->nextRigidBody(rbIter)) {
118 gezelter 1892 rb->updateAtomVel();
119 gezelter 1865 }
120     }
121 gezelter 1884
122 gezelter 1881 if (evaluator_.isDynamic()) {
123     seleMan_.setSelectionSet(evaluator_.evaluate());
124     }
125 gezelter 1865
126 gezelter 1881 // loop over the selected atoms:
127    
128     for (sd = seleMan_.beginSelected(i); sd != NULL;
129     sd = seleMan_.nextSelected(i)) {
130    
131     // figure out where that object is:
132     Vector3d pos = sd->getPos();
133 gezelter 1883 Vector3d vel = sd->getVel();
134     RealType m = sd->getMass();
135    
136 gezelter 1881 int bin = getBin(pos);
137 gezelter 1884
138 gezelter 1881 binMass[bin] += m;
139 gezelter 1944 binP[bin] += m * vel;
140 gezelter 1881 binKE[bin] += 0.5 * (m * vel.lengthSquare());
141     binDof[bin] += 3;
142    
143     if (sd->isDirectional()) {
144     Vector3d angMom = sd->getJ();
145     Mat3x3d I = sd->getI();
146     if (sd->isLinear()) {
147     int i = sd->linearAxis();
148     int j = (i + 1) % 3;
149     int k = (i + 2) % 3;
150     binKE[bin] += 0.5 * (angMom[j] * angMom[j] / I(j, j) +
151     angMom[k] * angMom[k] / I(k, k));
152     binDof[bin] += 2;
153     } else {
154     binKE[bin] += 0.5 * (angMom[0] * angMom[0] / I(0, 0) +
155     angMom[1] * angMom[1] / I(1, 1) +
156     angMom[2] * angMom[2] / I(2, 2));
157     binDof[bin] += 3;
158     }
159     }
160     }
161    
162 gezelter 1953 for (int i = 0; i < nBins_; i++) {
163 gezelter 1885
164 gezelter 1882 if (binDof[i] > 0) {
165     RealType temp = 2.0 * binKE[i] / (binDof[i] * PhysicalConstants::kb *
166     PhysicalConstants::energyConvert);
167     RealType den = binMass[i] * nBins_ * PhysicalConstants::densityConvert
168     / volume_;
169 gezelter 1944 Vector3d vel = binP[i] / binMass[i];
170 gezelter 1892
171 gezelter 1882 dynamic_cast<Accumulator *>(temperature->accumulator[i])->add(temp);
172     dynamic_cast<VectorAccumulator *>(velocity->accumulator[i])->add(vel);
173     dynamic_cast<Accumulator *>(density->accumulator[i])->add(den);
174     dynamic_cast<Accumulator *>(counts_->accumulator[i])->add(1);
175     }
176 gezelter 1881 }
177 gezelter 1865 }
178 gezelter 1881
179     void RNEMDZ::processStuntDouble(StuntDouble* sd, int bin) {
180     }
181 gezelter 1865
182     RNEMDR::RNEMDR(SimInfo* info, const std::string& filename,
183     const std::string& sele, int nrbins)
184     : ShellStatistics(info, filename, sele, nrbins) {
185    
186    
187     setOutputName(getPrefix(filename) + ".rnemdR");
188    
189     temperature = new OutputData;
190     temperature->units = "K";
191     temperature->title = "Temperature";
192     temperature->dataType = odtReal;
193     temperature->dataHandling = odhAverage;
194     temperature->accumulator.reserve(nBins_);
195     for (int i = 0; i < nBins_; i++)
196     temperature->accumulator.push_back( new Accumulator() );
197     data_.push_back(temperature);
198    
199     angularVelocity = new OutputData;
200     angularVelocity->units = "angstroms^2/fs";
201 gezelter 1953 angularVelocity->title = "Angular Velocity";
202 gezelter 1865 angularVelocity->dataType = odtVector3;
203     angularVelocity->dataHandling = odhAverage;
204     angularVelocity->accumulator.reserve(nBins_);
205     for (int i = 0; i < nBins_; i++)
206     angularVelocity->accumulator.push_back( new VectorAccumulator() );
207     data_.push_back(angularVelocity);
208    
209     density = new OutputData;
210     density->units = "g cm^-3";
211     density->title = "Density";
212     density->dataType = odtReal;
213     density->dataHandling = odhAverage;
214     density->accumulator.reserve(nBins_);
215     for (int i = 0; i < nBins_; i++)
216     density->accumulator.push_back( new Accumulator() );
217     data_.push_back(density);
218     }
219    
220    
221 gezelter 1887 void RNEMDR::processFrame(int istep) {
222 gezelter 1865
223 gezelter 1887 Molecule* mol;
224     RigidBody* rb;
225     StuntDouble* sd;
226     SimInfo::MoleculeIterator mi;
227     Molecule::RigidBodyIterator rbIter;
228     int i;
229    
230     vector<RealType> binMass(nBins_, 0.0);
231 gezelter 1944 vector<Mat3x3d> binI(nBins_);
232     vector<Vector3d> binL(nBins_, V3Zero);
233 gezelter 1887 vector<RealType> binKE(nBins_, 0.0);
234 gezelter 1944 vector<int> binDof(nBins_, 0);
235 gezelter 1887
236     for (mol = info_->beginMolecule(mi); mol != NULL;
237     mol = info_->nextMolecule(mi)) {
238    
239     // change the positions of atoms which belong to the rigidbodies
240    
241     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
242     rb = mol->nextRigidBody(rbIter)) {
243 gezelter 1892 rb->updateAtomVel();
244 gezelter 1865 }
245     }
246 gezelter 1887
247     if (evaluator_.isDynamic()) {
248     seleMan_.setSelectionSet(evaluator_.evaluate());
249     }
250 gezelter 1865
251 gezelter 1887 // loop over the selected atoms:
252    
253     for (sd = seleMan_.beginSelected(i); sd != NULL;
254     sd = seleMan_.nextSelected(i)) {
255 gezelter 1944
256 gezelter 1887 // figure out where that object is:
257 gezelter 1945 int bin = getBin( sd->getPos() );
258 gezelter 1865
259 gezelter 1944 if (bin >= 0 && bin < nBins_) {
260 gezelter 1887
261 gezelter 1944 Vector3d rPos = sd->getPos() - coordinateOrigin_;
262     Vector3d vel = sd->getVel();
263     RealType m = sd->getMass();
264     Vector3d L = m * cross(rPos, vel);
265     Mat3x3d I(0.0);
266     I = outProduct(rPos, rPos) * m;
267     RealType r2 = rPos.lengthSquare();
268     I(0, 0) += m * r2;
269     I(1, 1) += m * r2;
270     I(2, 2) += m * r2;
271 gezelter 1887
272 gezelter 1944 binMass[bin] += m;
273     binI[bin] += I;
274     binL[bin] += L;
275     binKE[bin] += 0.5 * (m * vel.lengthSquare());
276     binDof[bin] += 3;
277    
278     if (sd->isDirectional()) {
279     Vector3d angMom = sd->getJ();
280     Mat3x3d Ia = sd->getI();
281     if (sd->isLinear()) {
282     int i = sd->linearAxis();
283     int j = (i + 1) % 3;
284     int k = (i + 2) % 3;
285     binKE[bin] += 0.5 * (angMom[j] * angMom[j] / Ia(j, j) +
286     angMom[k] * angMom[k] / Ia(k, k));
287     binDof[bin] += 2;
288     } else {
289     binKE[bin] += 0.5 * (angMom[0] * angMom[0] / Ia(0, 0) +
290     angMom[1] * angMom[1] / Ia(1, 1) +
291     angMom[2] * angMom[2] / Ia(2, 2));
292     binDof[bin] += 3;
293     }
294 gezelter 1887 }
295     }
296     }
297 gezelter 1865
298 gezelter 1953 for (int i = 0; i < nBins_; i++) {
299 gezelter 1887 RealType rinner = (RealType)i * binWidth_;
300     RealType router = (RealType)(i+1) * binWidth_;
301     if (binDof[i] > 0) {
302     RealType temp = 2.0 * binKE[i] / (binDof[i] * PhysicalConstants::kb *
303     PhysicalConstants::energyConvert);
304     RealType den = binMass[i] * 3.0 * PhysicalConstants::densityConvert
305 gezelter 1944 / (4.0 * M_PI * (pow(router,3) - pow(rinner,3)));
306    
307     Vector3d omega = binI[i].inverse() * binL[i];
308    
309 gezelter 1887 dynamic_cast<Accumulator *>(temperature->accumulator[i])->add(temp);
310 gezelter 1944 dynamic_cast<VectorAccumulator *>(angularVelocity->accumulator[i])->add(omega);
311 gezelter 1887 dynamic_cast<Accumulator *>(density->accumulator[i])->add(den);
312     dynamic_cast<Accumulator *>(counts_->accumulator[i])->add(1);
313     }
314     }
315     }
316 gezelter 1865
317 gezelter 1887
318     void RNEMDR::processStuntDouble(StuntDouble* sd, int bin) {
319 gezelter 1865 }
320 gezelter 1953
321     RNEMDRTheta::RNEMDRTheta(SimInfo* info, const std::string& filename,
322     const std::string& sele, int nrbins, int nangleBins)
323     : ShellStatistics(info, filename, sele, nrbins), nAngleBins_(nangleBins) {
324    
325     Globals* simParams = info->getSimParams();
326     RNEMDParameters* rnemdParams = simParams->getRNEMDParameters();
327     bool hasAngularMomentumFluxVector = rnemdParams->haveAngularMomentumFluxVector();
328    
329     if (hasAngularMomentumFluxVector) {
330     fluxVector_ = rnemdParams->getAngularMomentumFluxVector();
331     } else {
332    
333     std::string fluxStr = rnemdParams->getFluxType();
334     if (fluxStr.find("Lx") != std::string::npos) {
335     fluxVector_ = V3X;
336     } else if (fluxStr.find("Ly") != std::string::npos) {
337     fluxVector_ = V3Y;
338     } else {
339     fluxVector_ = V3Z;
340     }
341     }
342    
343     fluxVector_.normalize();
344    
345     setOutputName(getPrefix(filename) + ".rnemdRTheta");
346    
347     angularVelocity = new OutputData;
348     angularVelocity->units = "angstroms^2/fs";
349     angularVelocity->title = "Angular Velocity";
350     angularVelocity->dataType = odtArray2d;
351     angularVelocity->dataHandling = odhAverage;
352     angularVelocity->accumulatorArray2d.reserve(nBins_);
353     for (int i = 0; i < nBins_; i++) {
354     angularVelocity->accumulatorArray2d[i].reserve(nAngleBins_);
355     for (int j = 0 ; j < nAngleBins_; j++) {
356     angularVelocity->accumulatorArray2d[i][j] = new Accumulator();
357     }
358     }
359     data_.push_back(angularVelocity);
360    
361     }
362    
363    
364     std::pair<int,int> RNEMDRTheta::getBins(Vector3d pos) {
365     std::pair<int,int> result;
366    
367     Vector3d rPos = pos - coordinateOrigin_;
368     RealType cosAngle= dot(rPos, fluxVector_) / rPos.length();
369    
370     result.first = int(rPos.length() / binWidth_);
371     result.second = int( (nAngleBins_ - 1) * 0.5 * (cosAngle + 1.0) );
372     return result;
373     }
374    
375     void RNEMDRTheta::processStuntDouble(StuntDouble* sd, int bin) {
376     }
377    
378     void RNEMDRTheta::processFrame(int istep) {
379    
380     Molecule* mol;
381     RigidBody* rb;
382     StuntDouble* sd;
383     SimInfo::MoleculeIterator mi;
384     Molecule::RigidBodyIterator rbIter;
385     int i;
386    
387     vector<vector<Mat3x3d> > binI;
388     vector<vector<Vector3d> > binL;
389     vector<vector<int> > binCount;
390    
391     for (mol = info_->beginMolecule(mi); mol != NULL;
392     mol = info_->nextMolecule(mi)) {
393    
394     // change the positions of atoms which belong to the rigidbodies
395    
396     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
397     rb = mol->nextRigidBody(rbIter)) {
398     rb->updateAtomVel();
399     }
400     }
401    
402     if (evaluator_.isDynamic()) {
403     seleMan_.setSelectionSet(evaluator_.evaluate());
404     }
405    
406     // loop over the selected atoms:
407    
408     for (sd = seleMan_.beginSelected(i); sd != NULL;
409     sd = seleMan_.nextSelected(i)) {
410    
411     // figure out where that object is:
412     std::pair<int,int> bins = getBins( sd->getPos() );
413    
414     if (bins.first >= 0 && bins.first < nBins_) {
415     if (bins.second >= 0 && bins.second < nAngleBins_) {
416    
417     Vector3d rPos = sd->getPos() - coordinateOrigin_;
418     Vector3d vel = sd->getVel();
419     RealType m = sd->getMass();
420     Vector3d L = m * cross(rPos, vel);
421     Mat3x3d I(0.0);
422     I = outProduct(rPos, rPos) * m;
423     RealType r2 = rPos.lengthSquare();
424     I(0, 0) += m * r2;
425     I(1, 1) += m * r2;
426     I(2, 2) += m * r2;
427    
428     binI[bins.first][bins.second] += I;
429     binL[bins.first][bins.second] += L;
430     binCount[bins.first][bins.second]++;
431     }
432     }
433     }
434    
435    
436     for (int i = 0; i < nBins_; i++) {
437     for (int j = 0; j < nAngleBins_; j++) {
438    
439     if (binCount[i][j] > 0) {
440     Vector3d omega = binI[i][j].inverse() * binL[i][j];
441     RealType omegaProj = dot(omega, fluxVector_);
442    
443     dynamic_cast<Accumulator *>(angularVelocity->accumulatorArray2d[i][j])->add(omegaProj);
444     }
445     }
446     }
447     }
448    
449     void RNEMDRTheta::writeOutput() {
450    
451     vector<OutputData*>::iterator i;
452     OutputData* outputData;
453    
454     ofstream outStream(outputFilename_.c_str());
455     if (outStream.is_open()) {
456    
457     //write title
458     outStream << "# SPATIAL STATISTICS\n";
459     outStream << "#";
460    
461     for(outputData = beginOutputData(i); outputData;
462     outputData = nextOutputData(i)) {
463     outStream << "\t" << outputData->title <<
464     "(" << outputData->units << ")";
465     // add some extra tabs for column alignment
466     if (outputData->dataType == odtVector3) outStream << "\t\t";
467     }
468    
469     outStream << std::endl;
470    
471     outStream.precision(8);
472    
473     for (int j = 0; j < nBins_; j++) {
474    
475     int counts = counts_->accumulator[j]->count();
476    
477     if (counts > 0) {
478     for(outputData = beginOutputData(i); outputData;
479     outputData = nextOutputData(i)) {
480    
481     int n = outputData->accumulator[j]->count();
482     if (n != 0) {
483     writeData( outStream, outputData, j );
484     }
485     }
486     outStream << std::endl;
487     }
488     }
489     }
490     }
491 gezelter 1865 }
492    

Properties

Name Value
svn:eol-style native