ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/Snapshot.cpp
Revision: 1858
Committed: Wed Apr 3 21:32:13 2013 UTC (12 years, 1 month ago) by gezelter
File size: 18792 byte(s)
Log Message:
Some bugfixes for cell-linked-list-style neighbor lists when the simulation
doesn't use periodic boundary conditions.   Cleaning up some of the Hull
stuff while we're in there.

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

Properties

Name Value
svn:keywords Author Id Revision Date