ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/Snapshot.hpp
(Generate patch)

Comparing:
trunk/src/brains/Snapshot.hpp (file contents), Revision 963 by tim, Wed May 17 21:51:42 2006 UTC vs.
branches/development/src/brains/Snapshot.hpp (file contents), Revision 1774 by gezelter, Wed Aug 8 16:03:02 2012 UTC

# Line 6 | Line 6
6   * redistribute this software in source and binary code form, provided
7   * that the following conditions are met:
8   *
9 < * 1. Acknowledgement of the program authors must be made in any
10 < *    publication of scientific results based in part on use of the
11 < *    program.  An acceptable form of acknowledgement is citation of
12 < *    the article in which the program was described (Matthew
13 < *    A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 < *    J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 < *    Parallel Simulation Engine for Molecular Dynamics,"
16 < *    J. Comput. Chem. 26, pp. 252-271 (2005))
17 < *
18 < * 2. Redistributions of source code must retain the above copyright
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
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.
# Line 37 | Line 28
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, 24107 (2008).          
39 + * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 + * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42    
42 /**
43 * @file Snapshot.hpp
44 * @author tlin
45 * @date 10/20/2004
46 * @time 23:56am
47 * @version 1.0
48 */
49  
43   #ifndef BRAINS_SNAPSHOT_HPP
44   #define BRAINS_SNAPSHOT_HPP
45  
46   #include <vector>
47  
48   #include "brains/DataStorage.hpp"
49 + #include "nonbonded/NonBondedInteraction.hpp"
50   #include "brains/Stats.hpp"
57 #include "UseTheForce/DarkSide/simulation_interface.h"
51  
52 < namespace oopse{
52 > using namespace std;
53 > namespace OpenMD{
54  
55    /**
56 <   * @class Snapshot Snapshot.hpp "brains/Snapshot.hpp"
57 <   * @brief Snapshot class is a repository class for storing dynamic data during
64 <   *  Simulation
65 <   * Every snapshot class will contain one DataStorage  for atoms and one DataStorage
66 <   *  for rigid bodies.
56 >   * FrameData is a structure for holding system-wide dynamic data
57 >   * about the simulation.
58     */
59 <  class Snapshot {
60 <  public:
61 <            
62 <    Snapshot(int nAtoms, int nRigidbodies) : atomData(nAtoms), rigidbodyData(nRigidbodies),
63 <                                             currentTime_(0), orthoRhombic_(0), chi_(0.0), integralOfChiDt_(0.0), eta_(0.0), id_(-1) {
59 >  
60 >  struct FrameData {
61 >    int id;                       /**< identification number of the snapshot */
62 >    RealType currentTime;         /**< current time */
63 >    Mat3x3d  hmat;                /**< axes of the periodic box in matrix form */
64 >    Mat3x3d  invHmat;             /**< the inverse of the Hmat matrix */
65 >    bool     orthoRhombic;        /**< is this an orthorhombic periodic box? */
66 >    RealType totalEnergy;         /**< total energy of this frame */
67 >    RealType translationalKinetic; /**< translational kinetic energy of this frame */
68 >    RealType rotationalKinetic;   /**< rotational kinetic energy of this frame */
69 >    RealType kineticEnergy;       /**< kinetic energy of this frame */
70 >    RealType potentialEnergy;     /**< potential energy of this frame */
71 >    RealType shortRangePotential; /**< short-range contributions to the potential*/
72 >    RealType longRangePotential;  /**< long-range contributions to the potential */
73 >    RealType bondPotential;       /**< bonded contribution to the potential */
74 >    RealType bendPotential;       /**< angle-bending contribution to the potential */
75 >    RealType torsionPotential;    /**< dihedral (torsion angle) contribution to the potential */
76 >    RealType inversionPotential;  /**< inversion (planarity) contribution to the potential */
77 >    potVec   lrPotentials;        /**< breakdown of long-range potentials by family */
78 >    potVec   excludedPotentials;  /**< breakdown of excluded potentials by family */
79 >    RealType restraintPotential;  /**< potential energy of restraints */
80 >    RealType rawPotential;        /**< unrestrained potential energy (when restraints are applied) */
81 >    RealType xyArea;              /**< XY area of this frame */
82 >    RealType volume;              /**< total volume of this frame */
83 >    RealType pressure;            /**< pressure of this frame */
84 >    RealType temperature;         /**< temperature of this frame */
85 >    pair<RealType, RealType> thermostat;    /**< thermostat variables */
86 >    RealType electronicTemperature; /**< temperature of the electronic degrees of freedom */
87 >    pair<RealType, RealType> electronicThermostat; /**< thermostat variables for electronic degrees of freedom */
88 >    Mat3x3d  barostat;            /**< barostat matrix */
89 >    Vector3d COM;                 /**< location of system center of mass */
90 >    Vector3d COMvel;              /**< system center of mass velocity */
91 >    Vector3d COMw;                /**< system center of mass angular velocity */
92 >    Mat3x3d  inertiaTensor;       /**< inertia tensor for entire system */
93 >    RealType gyrationalVolume;    /**< gyrational volume for entire system */
94 >    RealType hullVolume;          /**< hull volume for entire system */
95 >    Mat3x3d  stressTensor;        /**< stress tensor */
96 >    Mat3x3d  pressureTensor;      /**< pressure tensor */
97 >    Vector3d systemDipole;        /**< total system dipole moment */
98 >    Vector3d conductiveHeatFlux;  /**< heat flux vector (conductive only) */
99 >    Vector3d convectiveHeatFlux;  /**< heat flux vector (convective only) */
100 >    RealType conservedQuantity;   /**< anything conserved by the integrator */
101 >  };
102  
74    }
103  
104 <    Snapshot(int nAtoms, int nRigidbodies, int storageLayout)
105 <      : atomData(nAtoms, storageLayout), rigidbodyData(nRigidbodies, storageLayout),
106 <        currentTime_(0), orthoRhombic_(0), chi_(0.0), integralOfChiDt_(0.0), eta_(0.0), id_(-1) {
104 >  /**
105 >   * @class Snapshot
106 >   * @brief The Snapshot class is a repository storing dynamic data during a
107 >   * Simulation.  Every Snapshot contains FrameData (for global information)
108 >   * as well as DataStorage (one for Atoms, one for RigidBodies, and one for
109 >   * CutoffGroups).
110 >   */
111 >  class Snapshot {
112  
113 <      }
114 <            
113 >  public:            
114 >    Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups);
115 >    Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups, int storageLayout);    
116      /** Returns the id of this Snapshot */
117 <    int getID() {
84 <      return id_;
85 <    }
86 <
117 >    int      getID();
118      /** Sets the id of this Snapshot */
119 <    void setID(int id) {
89 <      id_ = id;
90 <    }
119 >    void     setID(int id);
120  
121 <    int getSize() {
122 <      return atomData.getSize() + rigidbodyData.getSize();
94 <    }
121 >    /** sets the state of the computed properties to false */
122 >    void     clearDerivedProperties();
123  
124 +    int      getSize();
125      /** Returns the number of atoms */
126 <    int getNumberOfAtoms() {
98 <      return atomData.getSize();
99 <    }
100 <
126 >    int      getNumberOfAtoms();
127      /** Returns the number of rigid bodies */
128 <    int getNumberOfRigidBodies() {
129 <      return rigidbodyData.getSize();
130 <    }
128 >    int      getNumberOfRigidBodies();
129 >    /** Returns the number of rigid bodies */
130 >    int      getNumberOfCutoffGroups();
131  
132      /** Returns the H-Matrix */
133 <    Mat3x3d getHmat() {
108 <      return hmat_;
109 <    }
110 <
133 >    Mat3x3d  getHmat();
134      /** Sets the H-Matrix */
135 <    void setHmat(const Mat3x3d& m);
135 >    void     setHmat(const Mat3x3d& m);
136 >    /** Returns the inverse H-Matrix */
137 >    Mat3x3d  getInvHmat();
138              
139 <    RealType getVolume() {
140 <      return hmat_.determinant();
141 <    }
139 >    RealType getVolume();
140 >    RealType getXYarea();
141 >    void     setVolume(const RealType vol);
142  
118    /** Returns the inverse H-Matrix */
119    Mat3x3d getInvHmat() {
120      return invHmat_;
121    }
122
143      /** Wrapping the vector according to periodic boundary condition*/
144 <    void wrapVector(Vector3d& v);
144 >    void     wrapVector(Vector3d& v);
145  
146 +    /** Scaling a vector to multiples of the periodic box */
147 +    Vector3d scaleVector(Vector3d &v);
148 +
149 +    void     setCOM(const Vector3d &com);
150 +    void     setCOMvel(const Vector3d &comVel);
151 +    void     setCOMw(const Vector3d &comw);
152 +
153 +    Vector3d getCOM();
154 +    Vector3d getCOMvel();
155 +    Vector3d getCOMw();
156              
157 <    RealType getTime() {
158 <      return currentTime_;
159 <    }
157 >    RealType getTime();
158 >    void     increaseTime(const RealType dt);
159 >    void     setTime(const RealType time);
160  
161 <    void increaseTime(RealType dt) {
162 <      setTime(getTime() + dt);
163 <    }
161 >    void     setBondPotential(const RealType bp);
162 >    void     setBendPotential(const RealType bp);
163 >    void     setTorsionPotential(const RealType tp);
164 >    void     setInversionPotential(const RealType ip);
165 >    RealType getBondPotential();
166 >    RealType getBendPotential();
167 >    RealType getTorsionPotential();
168 >    RealType getInversionPotential();
169  
170 <    void setTime(RealType time) {
136 <      currentTime_ =time;
137 <      //time at statData is redundant
138 <      statData[Stats::TIME] = currentTime_;
139 <    }
170 >    RealType getShortRangePotential();
171  
172 <    RealType getChi() {
173 <      return chi_;
174 <    }
172 >    void     setLongRangePotential(const potVec lrPot);
173 >    RealType getLongRangePotential();
174 >    potVec   getLongRangePotentials();
175  
176 <    void setChi(RealType chi) {
177 <      chi_ = chi;
178 <    }
176 >    void     setExcludedPotentials(const potVec exPot);
177 >    potVec   getExcludedPotentials();
178 >  
179 >    void     setRestraintPotential(const RealType rp);
180 >    RealType getRestraintPotential();
181  
182 <    RealType getIntegralOfChiDt() {
183 <      return integralOfChiDt_;
151 <    }
182 >    void     setRawPotential(const RealType rp);
183 >    RealType getRawPotential();
184  
185 <    void setIntegralOfChiDt(RealType integralOfChiDt) {
186 <      integralOfChiDt_ = integralOfChiDt;
187 <    }
188 <            
189 <    Mat3x3d getEta() {
190 <      return eta_;
191 <    }
185 >    RealType getPotentialEnergy();
186 >    RealType getKineticEnergy();
187 >    RealType getTranslationalKineticEnergy();
188 >    RealType getRotationalKineticEnergy();
189 >    void     setKineticEnergy(const RealType ke);
190 >    void     setTranslationalKineticEnergy(const RealType tke);
191 >    void     setRotationalKineticEnergy(const RealType rke);
192 >    RealType getTotalEnergy();
193 >    void     setTotalEnergy(const RealType te);
194 >    RealType getConservedQuantity();
195 >    void     setConservedQuantity(const RealType cq);
196 >    RealType getTemperature();
197 >    void     setTemperature(const RealType temp);
198 >    RealType getElectronicTemperature();
199 >    void     setElectronicTemperature(const RealType eTemp);
200 >    RealType getPressure();
201 >    void     setPressure(const RealType pressure);
202  
203 <    void setEta(const Mat3x3d& eta) {
204 <      eta_ = eta;
205 <    }
203 >    Mat3x3d  getPressureTensor();
204 >    void     setPressureTensor(const Mat3x3d& pressureTensor);
205 >
206 >    Mat3x3d  getStressTensor();
207 >    void     setStressTensor(const Mat3x3d& stressTensor);
208 >
209 >    Vector3d getConductiveHeatFlux();
210 >    void     setConductiveHeatFlux(const Vector3d& chf);
211 >
212 >    Vector3d getConvectiveHeatFlux();
213 >    void     setConvectiveHeatFlux(const Vector3d& chf);
214 >
215 >    Vector3d getHeatFlux();
216 >    
217 >    Vector3d getSystemDipole();
218 >    void     setSystemDipole(const Vector3d& bd);
219 >
220 >    pair<RealType, RealType> getThermostat();
221 >    void setThermostat(const pair<RealType, RealType>& thermostat);
222 >
223 >    pair<RealType, RealType> getElectronicThermostat();
224 >    void setElectronicThermostat(const pair<RealType, RealType>& eThermostat);
225              
226 +    Mat3x3d  getBarostat();
227 +    void     setBarostat(const Mat3x3d& barostat);
228 +
229 +    Mat3x3d  getInertiaTensor();
230 +    void     setInertiaTensor(const Mat3x3d& inertiaTensor);
231 +
232 +    RealType getGyrationalVolume();
233 +    void     setGyrationalVolume(const RealType gv);
234 +
235 +    RealType getHullVolume();
236 +    void     setHullVolume(const RealType hv);
237 +    
238 +    void     setOrthoTolerance(RealType orthoTolerance);
239 +
240      DataStorage atomData;
241      DataStorage rigidbodyData;
242 <    Stats statData;
243 <            
169 <  private:
170 <    RealType currentTime_;
242 >    DataStorage cgData;
243 >    FrameData   frameData;
244  
245 <    Mat3x3d hmat_;
246 <    Mat3x3d invHmat_;
247 <    int orthoRhombic_;
245 >    bool hasTotalEnergy;        
246 >    bool hasTranslationalKineticEnergy;    
247 >    bool hasRotationalKineticEnergy;    
248 >    bool hasKineticEnergy;    
249 >    bool hasShortRangePotential;
250 >    bool hasLongRangePotential;
251 >    bool hasPotentialEnergy;    
252 >    bool hasXYarea;
253 >    bool hasVolume;        
254 >    bool hasPressure;      
255 >    bool hasTemperature;    
256 >    bool hasElectronicTemperature;
257 >    bool hasCOM;            
258 >    bool hasCOMvel;
259 >    bool hasCOMw;
260 >    bool hasPressureTensor;    
261 >    bool hasSystemDipole;    
262 >    bool hasConvectiveHeatFlux;
263 >    bool hasInertiaTensor;
264 >    bool hasGyrationalVolume;
265 >    bool hasHullVolume;
266 >    bool hasConservedQuantity;
267  
268 <    RealType chi_;
269 <    RealType integralOfChiDt_;
270 <    Mat3x3d eta_;
179 <            
180 <    int id_; /**< identification number of the snapshot */
268 >  private:
269 >    RealType orthoTolerance_;
270 >    
271    };
272  
273    typedef DataStorage (Snapshot::*DataStoragePointer);

Comparing:
trunk/src/brains/Snapshot.hpp (property svn:keywords), Revision 963 by tim, Wed May 17 21:51:42 2006 UTC vs.
branches/development/src/brains/Snapshot.hpp (property svn:keywords), Revision 1774 by gezelter, Wed Aug 8 16:03:02 2012 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines