ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/Snapshot.cpp
Revision: 1879
Committed: Sun Jun 16 15:15:42 2013 UTC (11 years, 10 months ago) by gezelter
File size: 18847 byte(s)
Log Message:
MERGE OpenMD development 1783:1878 into trunk

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * 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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * 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 gezelter 1390 *
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 gezelter 1879 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1782 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 246 */
42    
43 gezelter 507 /**
44     * @file Snapshot.cpp
45     * @author tlin
46     * @date 11/11/2004
47     * @version 1.0
48     */
49 gezelter 246
50     #include "brains/Snapshot.hpp"
51     #include "utils/NumericConstant.hpp"
52     #include "utils/simError.h"
53     #include "utils/Utility.hpp"
54 jmarr 1401 #include <cstdio>
55    
56 gezelter 1390 namespace OpenMD {
57 gezelter 246
58 gezelter 1782 Snapshot::Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups) :
59     atomData(nAtoms), rigidbodyData(nRigidbodies),
60     cgData(nCutoffGroups, DataStorage::dslPosition),
61     orthoTolerance_(1e-6) {
62 gezelter 246
63 gezelter 1782 frameData.id = -1;
64     frameData.currentTime = 0;
65     frameData.hmat = Mat3x3d(0.0);
66     frameData.invHmat = Mat3x3d(0.0);
67     frameData.orthoRhombic = false;
68     frameData.bondPotential = 0.0;
69     frameData.bendPotential = 0.0;
70     frameData.torsionPotential = 0.0;
71     frameData.inversionPotential = 0.0;
72     frameData.lrPotentials = potVec(0.0);
73     frameData.excludedPotentials = potVec(0.0);
74     frameData.restraintPotential = 0.0;
75     frameData.rawPotential = 0.0;
76     frameData.xyArea = 0.0;
77     frameData.volume = 0.0;
78     frameData.thermostat = make_pair(0.0, 0.0);
79     frameData.electronicThermostat = make_pair(0.0, 0.0);
80     frameData.barostat = Mat3x3d(0.0);
81     frameData.stressTensor = Mat3x3d(0.0);
82     frameData.conductiveHeatFlux = Vector3d(0.0, 0.0, 0.0);
83 chuckv 1112
84 gezelter 1782 clearDerivedProperties();
85     }
86    
87     Snapshot::Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups,
88     int storageLayout) :
89     atomData(nAtoms, storageLayout),
90     rigidbodyData(nRigidbodies, storageLayout),
91     cgData(nCutoffGroups, DataStorage::dslPosition),
92     orthoTolerance_(1e-6) {
93    
94     frameData.id = -1;
95     frameData.currentTime = 0;
96     frameData.hmat = Mat3x3d(0.0);
97 gezelter 1879 frameData.invHmat = Mat3x3d(0.0);
98     frameData.bBox = Mat3x3d(0.0);
99     frameData.invBbox = Mat3x3d(0.0);
100 gezelter 1782 frameData.orthoRhombic = false;
101     frameData.bondPotential = 0.0;
102     frameData.bendPotential = 0.0;
103     frameData.torsionPotential = 0.0;
104     frameData.inversionPotential = 0.0;
105     frameData.lrPotentials = potVec(0.0);
106     frameData.excludedPotentials = potVec(0.0);
107     frameData.restraintPotential = 0.0;
108     frameData.rawPotential = 0.0;
109     frameData.xyArea = 0.0;
110     frameData.volume = 0.0;
111     frameData.thermostat = make_pair(0.0, 0.0);
112     frameData.electronicThermostat = make_pair(0.0, 0.0);
113     frameData.barostat = Mat3x3d(0.0);
114     frameData.stressTensor = Mat3x3d(0.0);
115     frameData.conductiveHeatFlux = Vector3d(0.0, 0.0, 0.0);
116 gezelter 246
117 gezelter 1782 clearDerivedProperties();
118     }
119    
120     void Snapshot::clearDerivedProperties() {
121     frameData.totalEnergy = 0.0;
122     frameData.translationalKinetic = 0.0;
123     frameData.rotationalKinetic = 0.0;
124     frameData.kineticEnergy = 0.0;
125     frameData.potentialEnergy = 0.0;
126     frameData.shortRangePotential = 0.0;
127     frameData.longRangePotential = 0.0;
128     frameData.pressure = 0.0;
129     frameData.temperature = 0.0;
130     frameData.pressureTensor = Mat3x3d(0.0);
131     frameData.systemDipole = Vector3d(0.0);
132     frameData.convectiveHeatFlux = Vector3d(0.0, 0.0, 0.0);
133     frameData.electronicTemperature = 0.0;
134     frameData.COM = V3Zero;
135     frameData.COMvel = V3Zero;
136 gezelter 1879 frameData.COMw = V3Zero;
137 gezelter 1782
138     hasTotalEnergy = false;
139     hasTranslationalKineticEnergy = false;
140     hasRotationalKineticEnergy = false;
141     hasKineticEnergy = false;
142     hasShortRangePotential = false;
143     hasLongRangePotential = false;
144     hasPotentialEnergy = false;
145     hasXYarea = false;
146     hasVolume = false;
147     hasPressure = false;
148     hasTemperature = false;
149     hasElectronicTemperature = false;
150     hasCOM = false;
151     hasCOMvel = false;
152     hasCOMw = false;
153     hasPressureTensor = false;
154     hasSystemDipole = false;
155     hasConvectiveHeatFlux = false;
156     hasInertiaTensor = false;
157     hasGyrationalVolume = false;
158     hasHullVolume = false;
159     hasConservedQuantity = false;
160 gezelter 1879 hasBoundingBox = false;
161 gezelter 1782 }
162    
163     /** Returns the id of this Snapshot */
164     int Snapshot::getID() {
165     return frameData.id;
166     }
167    
168     /** Sets the id of this Snapshot */
169     void Snapshot::setID(int id) {
170     frameData.id = id;
171     }
172    
173     int Snapshot::getSize() {
174     return atomData.getSize() + rigidbodyData.getSize();
175     }
176    
177     /** Returns the number of atoms */
178     int Snapshot::getNumberOfAtoms() {
179     return atomData.getSize();
180     }
181    
182     /** Returns the number of rigid bodies */
183     int Snapshot::getNumberOfRigidBodies() {
184     return rigidbodyData.getSize();
185     }
186    
187     /** Returns the number of rigid bodies */
188     int Snapshot::getNumberOfCutoffGroups() {
189     return cgData.getSize();
190     }
191    
192     /** Returns the H-Matrix */
193     Mat3x3d Snapshot::getHmat() {
194     return frameData.hmat;
195     }
196    
197     /** Sets the H-Matrix */
198     void Snapshot::setHmat(const Mat3x3d& m) {
199     hasVolume = false;
200     frameData.hmat = m;
201     frameData.invHmat = frameData.hmat.inverse();
202    
203 gezelter 246 //determine whether the box is orthoTolerance or not
204 gezelter 1782 bool oldOrthoRhombic = frameData.orthoRhombic;
205 gezelter 246
206 gezelter 1782 RealType smallDiag = fabs(frameData.hmat(0, 0));
207     if(smallDiag > fabs(frameData.hmat(1, 1))) smallDiag = fabs(frameData.hmat(1, 1));
208     if(smallDiag > fabs(frameData.hmat(2, 2))) smallDiag = fabs(frameData.hmat(2, 2));
209 gezelter 1021 RealType tol = smallDiag * orthoTolerance_;
210 gezelter 246
211 gezelter 1782 frameData.orthoRhombic = true;
212 gezelter 246
213     for (int i = 0; i < 3; i++ ) {
214 gezelter 507 for (int j = 0 ; j < 3; j++) {
215     if (i != j) {
216 gezelter 1782 if (frameData.orthoRhombic) {
217     if ( fabs(frameData.hmat(i, j)) >= tol)
218     frameData.orthoRhombic = false;
219 gezelter 507 }
220     }
221     }
222 gezelter 246 }
223 gezelter 1782
224     if( oldOrthoRhombic != frameData.orthoRhombic){
225    
226 gezelter 1879 // It is finally time to suppress these warnings once and for
227     // all. They were annoying and not very informative.
228    
229     // if( frameData.orthoRhombic ) {
230     // sprintf( painCave.errMsg,
231     // "OpenMD is switching from the default Non-Orthorhombic\n"
232     // "\tto the faster Orthorhombic periodic boundary computations.\n"
233     // "\tThis is usually a good thing, but if you want the\n"
234     // "\tNon-Orthorhombic computations, make the orthoBoxTolerance\n"
235     // "\tvariable ( currently set to %G ) smaller.\n",
236     // orthoTolerance_);
237     // painCave.severity = OPENMD_INFO;
238     // simError();
239     // }
240     // else {
241     // sprintf( painCave.errMsg,
242     // "OpenMD is switching from the faster Orthorhombic to the more\n"
243     // "\tflexible Non-Orthorhombic periodic boundary computations.\n"
244     // "\tThis is usually because the box has deformed under\n"
245     // "\tNPTf integration. If you want to live on the edge with\n"
246     // "\tthe Orthorhombic computations, make the orthoBoxTolerance\n"
247     // "\tvariable ( currently set to %G ) larger.\n",
248     // orthoTolerance_);
249     // painCave.severity = OPENMD_WARNING;
250     // simError();
251     // }
252 gezelter 246 }
253 gezelter 1782 }
254    
255     /** Returns the inverse H-Matrix */
256     Mat3x3d Snapshot::getInvHmat() {
257     return frameData.invHmat;
258     }
259 gezelter 246
260 gezelter 1879 /** Returns the Bounding Box */
261     Mat3x3d Snapshot::getBoundingBox() {
262     return frameData.bBox;
263     }
264    
265     /** Sets the Bounding Box */
266     void Snapshot::setBoundingBox(const Mat3x3d& m) {
267     frameData.bBox = m;
268     frameData.invBbox = frameData.bBox.inverse();
269     hasBoundingBox = true;
270     }
271    
272     /** Returns the inverse Bounding Box */
273     Mat3x3d Snapshot::getInvBoundingBox() {
274     return frameData.invBbox;
275     }
276    
277 gezelter 1782 RealType Snapshot::getXYarea() {
278     if (!hasXYarea) {
279     Vector3d x = frameData.hmat.getColumn(0);
280     Vector3d y = frameData.hmat.getColumn(1);
281     frameData.xyArea = cross(x,y).length();
282     hasXYarea = true;
283     }
284     return frameData.xyArea;
285 gezelter 507 }
286 gezelter 246
287 gezelter 1782 RealType Snapshot::getVolume() {
288     if (!hasVolume) {
289     frameData.volume = frameData.hmat.determinant();
290     hasVolume = true;
291     }
292     return frameData.volume;
293     }
294 gezelter 246
295 gezelter 1782 void Snapshot::setVolume(RealType vol) {
296     hasVolume = true;
297     frameData.volume = vol;
298     }
299    
300     /** Wrap a vector according to periodic boundary conditions */
301 gezelter 507 void Snapshot::wrapVector(Vector3d& pos) {
302 gezelter 1782
303     Vector3d scaled = scaleVector(pos);
304    
305     for (int i = 0; i < 3; i++)
306     scaled[i] -= roundMe(scaled[i]);
307    
308     if( !frameData.orthoRhombic )
309     pos = frameData.hmat * scaled;
310     else {
311    
312     // calc the wrapped real coordinates from the wrapped scaled coordinates
313     for (int i=0; i<3; i++) {
314     pos[i] = scaled[i] * frameData.hmat(i, i);
315     }
316     }
317     }
318 gezelter 246
319 gezelter 1782 /** Scaling a vector to multiples of the periodic box */
320     inline Vector3d Snapshot::scaleVector(Vector3d& pos) {
321    
322 gezelter 246 Vector3d scaled;
323    
324 gezelter 1782 if( !frameData.orthoRhombic )
325     scaled = frameData.invHmat * pos;
326     else {
327 gezelter 507 // calc the scaled coordinates.
328 gezelter 1782 for (int i=0; i<3; i++)
329     scaled[i] = pos[i] * frameData.invHmat(i, i);
330     }
331 gezelter 246
332 gezelter 1782 return scaled;
333     }
334 gezelter 246
335 gezelter 1782 void Snapshot::setCOM(const Vector3d& com) {
336     frameData.COM = com;
337     hasCOM = true;
338     }
339    
340     void Snapshot::setCOMvel(const Vector3d& comVel) {
341     frameData.COMvel = comVel;
342     hasCOMvel = true;
343     }
344    
345     void Snapshot::setCOMw(const Vector3d& comw) {
346     frameData.COMw = comw;
347     hasCOMw = true;
348     }
349    
350     Vector3d Snapshot::getCOM() {
351     return frameData.COM;
352     }
353    
354     Vector3d Snapshot::getCOMvel() {
355     return frameData.COMvel;
356     }
357    
358     Vector3d Snapshot::getCOMw() {
359     return frameData.COMw;
360     }
361    
362     RealType Snapshot::getTime() {
363     return frameData.currentTime;
364     }
365    
366     void Snapshot::increaseTime(RealType dt) {
367     setTime(getTime() + dt);
368     }
369    
370     void Snapshot::setTime(RealType time) {
371     frameData.currentTime = time;
372     }
373    
374     void Snapshot::setBondPotential(RealType bp) {
375     frameData.bondPotential = bp;
376     }
377    
378     void Snapshot::setBendPotential(RealType bp) {
379     frameData.bendPotential = bp;
380     }
381    
382     void Snapshot::setTorsionPotential(RealType tp) {
383     frameData.torsionPotential = tp;
384     }
385    
386     void Snapshot::setInversionPotential(RealType ip) {
387     frameData.inversionPotential = ip;
388     }
389 gezelter 246
390    
391 gezelter 1782 RealType Snapshot::getBondPotential() {
392     return frameData.bondPotential;
393     }
394     RealType Snapshot::getBendPotential() {
395     return frameData.bendPotential;
396     }
397     RealType Snapshot::getTorsionPotential() {
398     return frameData.torsionPotential;
399     }
400     RealType Snapshot::getInversionPotential() {
401     return frameData.inversionPotential;
402     }
403    
404     RealType Snapshot::getShortRangePotential() {
405     if (!hasShortRangePotential) {
406     frameData.shortRangePotential = frameData.bondPotential;
407     frameData.shortRangePotential += frameData.bendPotential;
408     frameData.shortRangePotential += frameData.torsionPotential;
409     frameData.shortRangePotential += frameData.inversionPotential;
410     hasShortRangePotential = true;
411 gezelter 246 }
412 gezelter 1782 return frameData.shortRangePotential;
413     }
414 gezelter 246
415 gezelter 1782 void Snapshot::setLongRangePotential(potVec lrPot) {
416     frameData.lrPotentials = lrPot;
417 gezelter 507 }
418 gezelter 1782
419     RealType Snapshot::getLongRangePotential() {
420     if (!hasLongRangePotential) {
421     for (int i = 0; i < N_INTERACTION_FAMILIES; i++) {
422     frameData.longRangePotential += frameData.lrPotentials[i];
423     }
424     hasLongRangePotential = true;
425     }
426     return frameData.longRangePotential;
427     }
428 gezelter 246
429 gezelter 1782 potVec Snapshot::getLongRangePotentials() {
430     return frameData.lrPotentials;
431     }
432    
433     RealType Snapshot::getPotentialEnergy() {
434     if (!hasPotentialEnergy) {
435     frameData.potentialEnergy = this->getLongRangePotential();
436     frameData.potentialEnergy += this->getShortRangePotential();
437     hasPotentialEnergy = true;
438 gezelter 1104 }
439 gezelter 1782 return frameData.potentialEnergy;
440 gezelter 1104 }
441 gezelter 1782
442     void Snapshot::setExcludedPotentials(potVec exPot) {
443     frameData.excludedPotentials = exPot;
444     }
445    
446     potVec Snapshot::getExcludedPotentials() {
447     return frameData.excludedPotentials;
448     }
449    
450     void Snapshot::setRestraintPotential(RealType rp) {
451     frameData.restraintPotential = rp;
452     }
453 gezelter 1104
454 gezelter 1782 RealType Snapshot::getRestraintPotential() {
455     return frameData.restraintPotential;
456 gezelter 1104 }
457    
458 gezelter 1782 void Snapshot::setRawPotential(RealType rp) {
459     frameData.rawPotential = rp;
460 gezelter 1104 }
461 gezelter 1782
462     RealType Snapshot::getRawPotential() {
463     return frameData.rawPotential;
464     }
465    
466     RealType Snapshot::getTranslationalKineticEnergy() {
467     return frameData.translationalKinetic;
468     }
469    
470     RealType Snapshot::getRotationalKineticEnergy() {
471     return frameData.rotationalKinetic;
472     }
473    
474     RealType Snapshot::getKineticEnergy() {
475     return frameData.kineticEnergy;
476     }
477    
478     void Snapshot::setTranslationalKineticEnergy(RealType tke) {
479     hasTranslationalKineticEnergy = true;
480     frameData.translationalKinetic = tke;
481     }
482    
483     void Snapshot::setRotationalKineticEnergy(RealType rke) {
484     hasRotationalKineticEnergy = true;
485     frameData.rotationalKinetic = rke;
486     }
487    
488     void Snapshot::setKineticEnergy(RealType ke) {
489     hasKineticEnergy = true;
490     frameData.kineticEnergy = ke;
491     }
492    
493     RealType Snapshot::getTotalEnergy() {
494     return frameData.totalEnergy;
495     }
496    
497     void Snapshot::setTotalEnergy(RealType te) {
498     hasTotalEnergy = true;
499     frameData.totalEnergy = te;
500     }
501    
502     RealType Snapshot::getConservedQuantity() {
503     return frameData.conservedQuantity;
504     }
505    
506     void Snapshot::setConservedQuantity(RealType cq) {
507     hasConservedQuantity = true;
508     frameData.conservedQuantity = cq;
509     }
510    
511     RealType Snapshot::getTemperature() {
512     return frameData.temperature;
513     }
514    
515     void Snapshot::setTemperature(RealType temp) {
516     hasTemperature = true;
517     frameData.temperature = temp;
518     }
519    
520     RealType Snapshot::getElectronicTemperature() {
521     return frameData.electronicTemperature;
522     }
523    
524     void Snapshot::setElectronicTemperature(RealType eTemp) {
525     hasElectronicTemperature = true;
526     frameData.electronicTemperature = eTemp;
527     }
528    
529     RealType Snapshot::getPressure() {
530     return frameData.pressure;
531     }
532    
533     void Snapshot::setPressure(RealType pressure) {
534     hasPressure = true;
535     frameData.pressure = pressure;
536     }
537    
538     Mat3x3d Snapshot::getPressureTensor() {
539     return frameData.pressureTensor;
540     }
541    
542    
543     void Snapshot::setPressureTensor(const Mat3x3d& pressureTensor) {
544     hasPressureTensor = true;
545     frameData.pressureTensor = pressureTensor;
546     }
547    
548     void Snapshot::setStressTensor(const Mat3x3d& stressTensor) {
549     frameData.stressTensor = stressTensor;
550     }
551    
552     Mat3x3d Snapshot::getStressTensor() {
553     return frameData.stressTensor;
554     }
555    
556     void Snapshot::setConductiveHeatFlux(const Vector3d& chf) {
557     frameData.conductiveHeatFlux = chf;
558     }
559    
560     Vector3d Snapshot::getConductiveHeatFlux() {
561     return frameData.conductiveHeatFlux;
562     }
563    
564     Vector3d Snapshot::getConvectiveHeatFlux() {
565     return frameData.convectiveHeatFlux;
566     }
567    
568     void Snapshot::setConvectiveHeatFlux(const Vector3d& chf) {
569     hasConvectiveHeatFlux = true;
570     frameData.convectiveHeatFlux = chf;
571     }
572    
573     Vector3d Snapshot::getHeatFlux() {
574     // BE CAREFUL WITH UNITS
575     return getConductiveHeatFlux() + getConvectiveHeatFlux();
576     }
577    
578     Vector3d Snapshot::getSystemDipole() {
579     return frameData.systemDipole;
580     }
581    
582     void Snapshot::setSystemDipole(const Vector3d& bd) {
583     hasSystemDipole = true;
584     frameData.systemDipole = bd;
585     }
586    
587     void Snapshot::setThermostat(const pair<RealType, RealType>& thermostat) {
588     frameData.thermostat = thermostat;
589     }
590    
591     pair<RealType, RealType> Snapshot::getThermostat() {
592     return frameData.thermostat;
593     }
594    
595     void Snapshot::setElectronicThermostat(const pair<RealType, RealType>& eTherm) {
596     frameData.electronicThermostat = eTherm;
597     }
598    
599     pair<RealType, RealType> Snapshot::getElectronicThermostat() {
600     return frameData.electronicThermostat;
601     }
602 gezelter 1104
603 gezelter 1782 void Snapshot::setBarostat(const Mat3x3d& barostat) {
604     frameData.barostat = barostat;
605     }
606    
607     Mat3x3d Snapshot::getBarostat() {
608     return frameData.barostat;
609     }
610    
611     void Snapshot::setInertiaTensor(const Mat3x3d& inertiaTensor) {
612     frameData.inertiaTensor = inertiaTensor;
613     hasInertiaTensor = true;
614     }
615    
616     Mat3x3d Snapshot::getInertiaTensor() {
617     return frameData.inertiaTensor;
618     }
619    
620     void Snapshot::setGyrationalVolume(const RealType gyrationalVolume) {
621     frameData.gyrationalVolume = gyrationalVolume;
622     hasGyrationalVolume = true;
623     }
624    
625     RealType Snapshot::getGyrationalVolume() {
626     return frameData.gyrationalVolume;
627     }
628    
629     void Snapshot::setHullVolume(const RealType hullVolume) {
630     frameData.hullVolume = hullVolume;
631     hasHullVolume = true;
632     }
633    
634     RealType Snapshot::getHullVolume() {
635     return frameData.hullVolume;
636     }
637    
638     void Snapshot::setOrthoTolerance(RealType ot) {
639     orthoTolerance_ = ot;
640     }
641 gezelter 246 }

Properties

Name Value
svn:keywords Author Id Revision Date