ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/Snapshot.hpp
Revision: 1723
Committed: Thu May 24 20:59:54 2012 UTC (12 years, 11 months ago) by gezelter
File size: 11043 byte(s)
Log Message:
Bug fixes for heat flux import

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 tim 122 *
4 gezelter 246 * 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     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (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 tim 122 */
42 gezelter 246
43 gezelter 507 /**
44     * @file Snapshot.hpp
45     * @author tlin
46     * @date 10/20/2004
47     * @time 23:56am
48     * @version 1.0
49     */
50 tim 146
51 tim 122 #ifndef BRAINS_SNAPSHOT_HPP
52     #define BRAINS_SNAPSHOT_HPP
53 tim 137
54     #include <vector>
55    
56 tim 155 #include "brains/DataStorage.hpp"
57 gezelter 246 #include "brains/Stats.hpp"
58 tim 137
59 gezelter 1390 namespace OpenMD{
60 tim 122
61 gezelter 1715 struct FrameData {
62     int id; /**< identification number of the snapshot */
63     RealType currentTime; /**< current time */
64     Mat3x3d hmat; /**< axes of the periodic box in matrix form */
65     Mat3x3d invHmat; /**< the inverse of the Hmat matrix */
66     bool orthoRhombic; /**< is this an orthorhombic periodic box? */
67     RealType volume; /**< total volume of this frame */
68     RealType pressure; /**< pressure of this frame */
69     RealType totalEnergy; /**< total energy of this frame */
70     RealType kineticEnergy; /**< kinetic energy of this frame */
71     RealType potentialEnergy; /**< potential energy of this frame */
72     RealType temperature; /**< temperature of this frame */
73     RealType chi; /**< thermostat velocity */
74     RealType integralOfChiDt; /**< the actual thermostat */
75     RealType electronicTemperature; /**< temperature of the electronic degrees of freedom */
76     RealType chiQ; /**< fluctuating charge thermostat velocity */
77     RealType integralOfChiQDt; /**< the actual fluctuating charge thermostat */
78     Mat3x3d eta; /**< barostat matrix */
79     Vector3d COM; /**< location of center of mass */
80     Vector3d COMvel; /**< system center of mass velocity */
81     Vector3d COMw; /**< system center of mass angular velocity */
82 gezelter 1723 Mat3x3d stressTensor; /**< stress tensor */
83 gezelter 1715 Mat3x3d pressureTensor; /**< pressure tensor */
84     Vector3d systemDipole; /**< total system dipole moment */
85 gezelter 1723 Vector3d conductiveHeatFlux; /**< heat flux vector (conductive only) */
86 gezelter 1715 };
87    
88    
89 gezelter 507 /**
90     * @class Snapshot Snapshot.hpp "brains/Snapshot.hpp"
91     * @brief Snapshot class is a repository class for storing dynamic data during
92     * Simulation
93 gezelter 1540 * Every snapshot class will contain one DataStorage for atoms and one DataStorage
94 gezelter 507 * for rigid bodies.
95     */
96     class Snapshot {
97     public:
98 tim 146
99 gezelter 1540 Snapshot(int nAtoms, int nRigidbodies,
100     int nCutoffGroups) : atomData(nAtoms),
101     rigidbodyData(nRigidbodies),
102 gezelter 1715 cgData(nCutoffGroups, DataStorage::dslPosition),
103     orthoTolerance_(1e-6), hasCOM_(false), hasVolume_(false){
104    
105     frameData.id = -1;
106     frameData.currentTime = 0;
107     frameData.hmat = Mat3x3d(0.0);
108     frameData.invHmat = Mat3x3d(0.0);
109     frameData.orthoRhombic = false;
110     frameData.volume = 0.0;
111     frameData.pressure = 0.0;
112     frameData.totalEnergy = 0.0;
113     frameData.kineticEnergy = 0.0;
114     frameData.potentialEnergy = 0.0;
115     frameData.temperature = 0.0;
116     frameData.chi = 0.0;
117     frameData.integralOfChiDt = 0.0;
118     frameData.electronicTemperature = 0.0;
119     frameData.chiQ = 0.0;
120     frameData.integralOfChiQDt = 0.0;
121     frameData.eta = Mat3x3d(0.0);
122     frameData.COM = V3Zero;
123     frameData.COMvel = V3Zero;
124     frameData.COMw = V3Zero;
125 gezelter 1723 frameData.stressTensor = Mat3x3d(0.0);
126 gezelter 1715 frameData.pressureTensor = Mat3x3d(0.0);
127     frameData.systemDipole = V3Zero;
128 gezelter 1723 frameData.conductiveHeatFlux = V3Zero;
129 gezelter 507 }
130 tim 137
131 gezelter 1540 Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups,
132     int storageLayout) : atomData(nAtoms, storageLayout),
133     rigidbodyData(nRigidbodies, storageLayout),
134     cgData(nCutoffGroups, DataStorage::dslPosition),
135 gezelter 1715 orthoTolerance_(1e-6),
136     hasCOM_(false),
137 gezelter 1708 hasVolume_(false) {
138 gezelter 1715 frameData.id = -1;
139     frameData.currentTime = 0;
140     frameData.hmat = Mat3x3d(0.0);
141     frameData.invHmat = Mat3x3d(0.0);
142     frameData.orthoRhombic = false;
143     frameData.volume = 0.0;
144     frameData.pressure = 0.0;
145     frameData.totalEnergy = 0.0;
146     frameData.kineticEnergy = 0.0;
147     frameData.potentialEnergy = 0.0;
148     frameData.temperature = 0.0;
149     frameData.chi = 0.0;
150     frameData.integralOfChiDt = 0.0;
151     frameData.electronicTemperature = 0.0;
152     frameData.chiQ = 0.0;
153     frameData.integralOfChiQDt = 0.0;
154     frameData.eta = Mat3x3d(0.0);
155     frameData.COM = V3Zero;
156     frameData.COMvel = V3Zero;
157     frameData.COMw = V3Zero;
158 gezelter 1723 frameData.stressTensor = Mat3x3d(0.0);
159 gezelter 1715 frameData.pressureTensor = Mat3x3d(0.0);
160     frameData.systemDipole = V3Zero;
161 gezelter 1723 frameData.conductiveHeatFlux = V3Zero;
162 gezelter 1540 }
163    
164 gezelter 507 /** Returns the id of this Snapshot */
165     int getID() {
166 gezelter 1715 return frameData.id;
167 gezelter 507 }
168 tim 137
169 gezelter 507 /** Sets the id of this Snapshot */
170     void setID(int id) {
171 gezelter 1715 frameData.id = id;
172 gezelter 507 }
173 tim 137
174 gezelter 507 int getSize() {
175     return atomData.getSize() + rigidbodyData.getSize();
176     }
177 tim 152
178 gezelter 507 /** Returns the number of atoms */
179     int getNumberOfAtoms() {
180     return atomData.getSize();
181     }
182 tim 152
183 gezelter 507 /** Returns the number of rigid bodies */
184     int getNumberOfRigidBodies() {
185     return rigidbodyData.getSize();
186     }
187 tim 152
188 gezelter 1544 /** Returns the number of rigid bodies */
189     int getNumberOfCutoffGroups() {
190     return cgData.getSize();
191     }
192    
193 gezelter 507 /** Returns the H-Matrix */
194     Mat3x3d getHmat() {
195 gezelter 1715 return frameData.hmat;
196 gezelter 507 }
197 tim 152
198 gezelter 507 /** Sets the H-Matrix */
199     void setHmat(const Mat3x3d& m);
200 gezelter 246
201 tim 963 RealType getVolume() {
202 chuckv 1302 if (hasVolume_){
203 gezelter 1715 return frameData.volume;
204 chuckv 1302 }else{
205 gezelter 1715 return frameData.hmat.determinant();
206 chuckv 1302 }
207 gezelter 507 }
208 tim 152
209 chuckv 1302 void setVolume(RealType volume){
210     hasVolume_=true;
211 gezelter 1715 frameData.volume = volume;
212 chuckv 1302 }
213    
214 gezelter 507 /** Returns the inverse H-Matrix */
215     Mat3x3d getInvHmat() {
216 gezelter 1715 return frameData.invHmat;
217 gezelter 507 }
218 tim 152
219 gezelter 507 /** Wrapping the vector according to periodic boundary condition*/
220     void wrapVector(Vector3d& v);
221 gezelter 1562 /** Scaling a vector to multiples of the periodic box */
222     Vector3d scaleVector(Vector3d &v);
223    
224    
225 gezelter 1104 Vector3d getCOM();
226     Vector3d getCOMvel();
227     Vector3d getCOMw();
228 gezelter 246
229 tim 963 RealType getTime() {
230 gezelter 1715 return frameData.currentTime;
231 gezelter 507 }
232 tim 152
233 tim 963 void increaseTime(RealType dt) {
234 gezelter 507 setTime(getTime() + dt);
235     }
236 tim 152
237 tim 963 void setTime(RealType time) {
238 gezelter 1715 frameData.currentTime =time;
239 gezelter 507 //time at statData is redundant
240 gezelter 1715 statData[Stats::TIME] = frameData.currentTime;
241 gezelter 507 }
242 tim 152
243 tim 963 RealType getChi() {
244 gezelter 1715 return frameData.chi;
245 gezelter 507 }
246 gezelter 246
247 tim 963 void setChi(RealType chi) {
248 gezelter 1715 frameData.chi = chi;
249 gezelter 507 }
250 gezelter 246
251 tim 963 RealType getIntegralOfChiDt() {
252 gezelter 1715 return frameData.integralOfChiDt;
253 gezelter 507 }
254 gezelter 246
255 tim 963 void setIntegralOfChiDt(RealType integralOfChiDt) {
256 gezelter 1715 frameData.integralOfChiDt = integralOfChiDt;
257 gezelter 507 }
258 gezelter 246
259 gezelter 1715 RealType getChiElectronic() {
260     return frameData.chiQ;
261     }
262 gezelter 1021
263 gezelter 1715 void setChiElectronic(RealType chiQ) {
264     frameData.chiQ = chiQ;
265     }
266    
267     RealType getIntegralOfChiElectronicDt() {
268     return frameData.integralOfChiQDt;
269     }
270    
271     void setIntegralOfChiElectronicDt(RealType integralOfChiQDt) {
272     frameData.integralOfChiQDt = integralOfChiQDt;
273     }
274    
275    
276 gezelter 1021 void setOrthoTolerance(RealType orthoTolerance) {
277     orthoTolerance_ = orthoTolerance;
278     }
279    
280 gezelter 507 Mat3x3d getEta() {
281 gezelter 1715 return frameData.eta;
282 gezelter 507 }
283 gezelter 246
284 gezelter 507 void setEta(const Mat3x3d& eta) {
285 gezelter 1715 frameData.eta = eta;
286 gezelter 507 }
287 gezelter 1104
288 gezelter 1723 Mat3x3d getStressTensor() {
289     return frameData.stressTensor;
290 gezelter 1709 }
291    
292 gezelter 1723 void setStressTensor(const Mat3x3d& stressTensor) {
293     frameData.stressTensor = stressTensor;
294 gezelter 1709 }
295    
296 gezelter 1723 Vector3d getConductiveHeatFlux() {
297     return frameData.conductiveHeatFlux;
298     }
299    
300     void setConductiveHeatFlux(const Vector3d& heatFlux) {
301     frameData.conductiveHeatFlux = heatFlux;
302     }
303    
304 gezelter 1104 bool hasCOM() {
305     return hasCOM_;
306     }
307    
308     void setCOMprops(const Vector3d& COM, const Vector3d& COMvel, const Vector3d& COMw) {
309 gezelter 1715 frameData.COM = COM;
310     frameData.COMvel = COMvel;
311     frameData.COMw = COMw;
312 gezelter 1104 hasCOM_ = true;
313     }
314 gezelter 1540
315 gezelter 507 DataStorage atomData;
316     DataStorage rigidbodyData;
317 gezelter 1540 DataStorage cgData;
318 gezelter 1715 FrameData frameData;
319 gezelter 1541 Stats statData;
320 gezelter 1540
321 gezelter 507 private:
322 gezelter 1021 RealType orthoTolerance_;
323 gezelter 1104 bool hasCOM_;
324 gezelter 1715 bool hasVolume_;
325 gezelter 507 };
326 tim 122
327 gezelter 507 typedef DataStorage (Snapshot::*DataStoragePointer);
328 tim 122 }
329     #endif //BRAINS_SNAPSHOT_HPP

Properties

Name Value
svn:keywords Author Id Revision Date