| 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. |
| 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 |
> |
namespace OpenMD{ |
| 53 |
|
|
| 54 |
|
/** |
| 55 |
< |
* @class Snapshot Snapshot.hpp "brains/Snapshot.hpp" |
| 56 |
< |
* @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. |
| 55 |
> |
* FrameData is a structure for holding system-wide dynamic data |
| 56 |
> |
* about the simulation. |
| 57 |
|
*/ |
| 58 |
+ |
|
| 59 |
+ |
struct FrameData { |
| 60 |
+ |
int id; /**< identification number of the snapshot */ |
| 61 |
+ |
RealType currentTime; /**< current time */ |
| 62 |
+ |
Mat3x3d hmat; /**< axes of the periodic box in matrix form */ |
| 63 |
+ |
Mat3x3d invHmat; /**< the inverse of the Hmat matrix */ |
| 64 |
+ |
bool orthoRhombic; /**< is this an orthorhombic periodic box? */ |
| 65 |
+ |
RealType volume; /**< total volume of this frame */ |
| 66 |
+ |
RealType pressure; /**< pressure of this frame */ |
| 67 |
+ |
RealType totalEnergy; /**< total energy of this frame */ |
| 68 |
+ |
RealType kineticEnergy; /**< kinetic energy of this frame */ |
| 69 |
+ |
RealType potentialEnergy; /**< potential energy of this frame */ |
| 70 |
+ |
RealType shortRangePotential; /**< short-range contributions to the potential*/ |
| 71 |
+ |
RealType longRangePotential; /**< long-range contributions to the potential */ |
| 72 |
+ |
RealType bondPotential; /**< bonded contribution to the potential */ |
| 73 |
+ |
RealType bendPotential; /**< angle-bending contribution to the potential */ |
| 74 |
+ |
RealType torsionPotential; /**< dihedral (torsion angle) contribution to the potential */ |
| 75 |
+ |
RealType inversionPotential; /**< inversion (planarity) contribution to the potential */ |
| 76 |
+ |
potVec lrPotentials; /**< breakdown of long-range potentials by family */ |
| 77 |
+ |
potVec excludedPotentials; /**< breakdown of excluded potentials by family */ |
| 78 |
+ |
RealType restraintPotential; /**< potential energy of restraints */ |
| 79 |
+ |
RealType rawPotential; /**< unrestrained potential energy (when restraints are applied) */ |
| 80 |
+ |
RealType temperature; /**< temperature of this frame */ |
| 81 |
+ |
RealType chi; /**< thermostat velocity */ |
| 82 |
+ |
RealType integralOfChiDt; /**< thermostat position */ |
| 83 |
+ |
RealType electronicTemperature; /**< temperature of the electronic degrees of freedom */ |
| 84 |
+ |
RealType chiQ; /**< fluctuating charge thermostat velocity */ |
| 85 |
+ |
RealType integralOfChiQDt; /**< fluctuating charge thermostat position */ |
| 86 |
+ |
Mat3x3d eta; /**< barostat matrix */ |
| 87 |
+ |
Vector3d COM; /**< location of system center of mass */ |
| 88 |
+ |
Vector3d COMvel; /**< system center of mass velocity */ |
| 89 |
+ |
Vector3d COMw; /**< system center of mass angular velocity */ |
| 90 |
+ |
Mat3x3d stressTensor; /**< stress tensor */ |
| 91 |
+ |
Mat3x3d pressureTensor; /**< pressure tensor */ |
| 92 |
+ |
Vector3d systemDipole; /**< total system dipole moment */ |
| 93 |
+ |
Vector3d conductiveHeatFlux; /**< heat flux vector (conductive only) */ |
| 94 |
+ |
}; |
| 95 |
+ |
|
| 96 |
+ |
|
| 97 |
+ |
/** |
| 98 |
+ |
* @class Snapshot |
| 99 |
+ |
* @brief The Snapshot class is a repository storing dynamic data during a |
| 100 |
+ |
* Simulation. Every Snapshot contains FrameData (for global information) |
| 101 |
+ |
* as well as DataStorage (one for Atoms, one for RigidBodies, and one for |
| 102 |
+ |
* CutoffGroups). |
| 103 |
+ |
*/ |
| 104 |
|
class Snapshot { |
| 69 |
– |
public: |
| 70 |
– |
|
| 71 |
– |
Snapshot(int nAtoms, int nRigidbodies) : atomData(nAtoms), rigidbodyData(nRigidbodies), |
| 72 |
– |
currentTime_(0), orthoRhombic_(0), chi_(0.0), integralOfChiDt_(0.0), eta_(0.0), id_(-1) { |
| 105 |
|
|
| 106 |
+ |
public: |
| 107 |
+ |
Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups) : |
| 108 |
+ |
atomData(nAtoms), rigidbodyData(nRigidbodies), |
| 109 |
+ |
cgData(nCutoffGroups, DataStorage::dslPosition), |
| 110 |
+ |
orthoTolerance_(1e-6), hasCOM_(false), hasVolume_(false), |
| 111 |
+ |
hasShortRangePotential_(false), |
| 112 |
+ |
hasBondPotential_(false), hasBendPotential_(false), |
| 113 |
+ |
hasTorsionPotential_(false), hasInversionPotential_(false), |
| 114 |
+ |
hasLongRangePotential_(false), hasLongRangePotentialFamilies_(false), |
| 115 |
+ |
hasRestraintPotential_(false), hasRawPotential_(false), |
| 116 |
+ |
hasExcludedPotentials_(false) |
| 117 |
+ |
{ |
| 118 |
+ |
|
| 119 |
+ |
frameData.id = -1; |
| 120 |
+ |
frameData.currentTime = 0; |
| 121 |
+ |
frameData.hmat = Mat3x3d(0.0); |
| 122 |
+ |
frameData.invHmat = Mat3x3d(0.0); |
| 123 |
+ |
frameData.orthoRhombic = false; |
| 124 |
+ |
frameData.volume = 0.0; |
| 125 |
+ |
frameData.pressure = 0.0; |
| 126 |
+ |
frameData.totalEnergy = 0.0; |
| 127 |
+ |
frameData.kineticEnergy = 0.0; |
| 128 |
+ |
frameData.potentialEnergy = 0.0; |
| 129 |
+ |
frameData.temperature = 0.0; |
| 130 |
+ |
frameData.chi = 0.0; |
| 131 |
+ |
frameData.integralOfChiDt = 0.0; |
| 132 |
+ |
frameData.electronicTemperature = 0.0; |
| 133 |
+ |
frameData.chiQ = 0.0; |
| 134 |
+ |
frameData.integralOfChiQDt = 0.0; |
| 135 |
+ |
frameData.eta = Mat3x3d(0.0); |
| 136 |
+ |
frameData.COM = V3Zero; |
| 137 |
+ |
frameData.COMvel = V3Zero; |
| 138 |
+ |
frameData.COMw = V3Zero; |
| 139 |
+ |
frameData.stressTensor = Mat3x3d(0.0); |
| 140 |
+ |
frameData.pressureTensor = Mat3x3d(0.0); |
| 141 |
+ |
frameData.systemDipole = Vector3d(0.0); |
| 142 |
+ |
frameData.conductiveHeatFlux = Vector3d(0.0, 0.0, 0.0); |
| 143 |
|
} |
| 144 |
|
|
| 145 |
< |
Snapshot(int nAtoms, int nRigidbodies, int storageLayout) |
| 146 |
< |
: atomData(nAtoms, storageLayout), rigidbodyData(nRigidbodies, storageLayout), |
| 147 |
< |
currentTime_(0), orthoRhombic_(0), chi_(0.0), integralOfChiDt_(0.0), eta_(0.0), id_(-1) { |
| 145 |
> |
Snapshot(int nAtoms, int nRigidbodies, int nCutoffGroups, int storageLayout) : |
| 146 |
> |
atomData(nAtoms, storageLayout), |
| 147 |
> |
rigidbodyData(nRigidbodies, storageLayout), |
| 148 |
> |
cgData(nCutoffGroups, DataStorage::dslPosition), |
| 149 |
> |
orthoTolerance_(1e-6), hasCOM_(false), hasVolume_(false), |
| 150 |
> |
hasShortRangePotential_(false), |
| 151 |
> |
hasBondPotential_(false), hasBendPotential_(false), |
| 152 |
> |
hasTorsionPotential_(false), hasInversionPotential_(false), |
| 153 |
> |
hasLongRangePotential_(false), hasLongRangePotentialFamilies_(false), |
| 154 |
> |
hasRestraintPotential_(false), hasRawPotential_(false), |
| 155 |
> |
hasExcludedPotentials_(false) |
| 156 |
> |
{ |
| 157 |
|
|
| 158 |
< |
} |
| 159 |
< |
|
| 158 |
> |
frameData.id = -1; |
| 159 |
> |
frameData.currentTime = 0; |
| 160 |
> |
frameData.hmat = Mat3x3d(0.0); |
| 161 |
> |
frameData.invHmat = Mat3x3d(0.0); |
| 162 |
> |
frameData.orthoRhombic = false; |
| 163 |
> |
frameData.volume = 0.0; |
| 164 |
> |
frameData.pressure = 0.0; |
| 165 |
> |
frameData.totalEnergy = 0.0; |
| 166 |
> |
frameData.kineticEnergy = 0.0; |
| 167 |
> |
frameData.potentialEnergy = 0.0; |
| 168 |
> |
frameData.temperature = 0.0; |
| 169 |
> |
frameData.chi = 0.0; |
| 170 |
> |
frameData.integralOfChiDt = 0.0; |
| 171 |
> |
frameData.electronicTemperature = 0.0; |
| 172 |
> |
frameData.chiQ = 0.0; |
| 173 |
> |
frameData.integralOfChiQDt = 0.0; |
| 174 |
> |
frameData.eta = Mat3x3d(0.0); |
| 175 |
> |
frameData.COM = V3Zero; |
| 176 |
> |
frameData.COMvel = V3Zero; |
| 177 |
> |
frameData.COMw = V3Zero; |
| 178 |
> |
frameData.stressTensor = Mat3x3d(0.0); |
| 179 |
> |
frameData.pressureTensor = Mat3x3d(0.0); |
| 180 |
> |
frameData.systemDipole = V3Zero; |
| 181 |
> |
frameData.conductiveHeatFlux = Vector3d(0.0, 0.0, 0.0); |
| 182 |
> |
} |
| 183 |
> |
|
| 184 |
|
/** Returns the id of this Snapshot */ |
| 185 |
|
int getID() { |
| 186 |
< |
return id_; |
| 186 |
> |
return frameData.id; |
| 187 |
|
} |
| 188 |
|
|
| 189 |
|
/** Sets the id of this Snapshot */ |
| 190 |
|
void setID(int id) { |
| 191 |
< |
id_ = id; |
| 191 |
> |
frameData.id = id; |
| 192 |
|
} |
| 193 |
|
|
| 194 |
|
int getSize() { |
| 205 |
|
return rigidbodyData.getSize(); |
| 206 |
|
} |
| 207 |
|
|
| 208 |
+ |
/** Returns the number of rigid bodies */ |
| 209 |
+ |
int getNumberOfCutoffGroups() { |
| 210 |
+ |
return cgData.getSize(); |
| 211 |
+ |
} |
| 212 |
+ |
|
| 213 |
|
/** Returns the H-Matrix */ |
| 214 |
|
Mat3x3d getHmat() { |
| 215 |
< |
return hmat_; |
| 215 |
> |
return frameData.hmat; |
| 216 |
|
} |
| 217 |
|
|
| 218 |
|
/** Sets the H-Matrix */ |
| 219 |
|
void setHmat(const Mat3x3d& m); |
| 220 |
|
|
| 221 |
|
RealType getVolume() { |
| 222 |
< |
return hmat_.determinant(); |
| 222 |
> |
if (!hasVolume_) { |
| 223 |
> |
frameData.volume = frameData.hmat.determinant(); |
| 224 |
> |
hasVolume_ = true; |
| 225 |
> |
} |
| 226 |
> |
return frameData.volume; |
| 227 |
|
} |
| 228 |
|
|
| 229 |
+ |
void setVolume(RealType volume){ |
| 230 |
+ |
hasVolume_=true; |
| 231 |
+ |
frameData.volume = volume; |
| 232 |
+ |
} |
| 233 |
+ |
|
| 234 |
|
/** Returns the inverse H-Matrix */ |
| 235 |
|
Mat3x3d getInvHmat() { |
| 236 |
< |
return invHmat_; |
| 236 |
> |
return frameData.invHmat; |
| 237 |
|
} |
| 238 |
|
|
| 239 |
|
/** Wrapping the vector according to periodic boundary condition*/ |
| 240 |
|
void wrapVector(Vector3d& v); |
| 241 |
+ |
/** Scaling a vector to multiples of the periodic box */ |
| 242 |
+ |
Vector3d scaleVector(Vector3d &v); |
| 243 |
|
|
| 244 |
+ |
|
| 245 |
+ |
Vector3d getCOM(); |
| 246 |
+ |
Vector3d getCOMvel(); |
| 247 |
+ |
Vector3d getCOMw(); |
| 248 |
|
|
| 249 |
|
RealType getTime() { |
| 250 |
< |
return currentTime_; |
| 250 |
> |
return frameData.currentTime; |
| 251 |
|
} |
| 252 |
|
|
| 253 |
|
void increaseTime(RealType dt) { |
| 255 |
|
} |
| 256 |
|
|
| 257 |
|
void setTime(RealType time) { |
| 258 |
< |
currentTime_ =time; |
| 258 |
> |
frameData.currentTime =time; |
| 259 |
|
//time at statData is redundant |
| 260 |
< |
statData[Stats::TIME] = currentTime_; |
| 260 |
> |
statData[Stats::TIME] = frameData.currentTime; |
| 261 |
|
} |
| 262 |
|
|
| 263 |
+ |
void setShortRangePotential(RealType srp) { |
| 264 |
+ |
frameData.shortRangePotential = srp; |
| 265 |
+ |
hasShortRangePotential_ = true; |
| 266 |
+ |
statData[Stats::SHORT_RANGE_POTENTIAL] = frameData.shortRangePotential; |
| 267 |
+ |
} |
| 268 |
+ |
|
| 269 |
+ |
RealType getShortRangePotential() { |
| 270 |
+ |
return frameData.shortRangePotential; |
| 271 |
+ |
} |
| 272 |
+ |
|
| 273 |
+ |
void setBondPotential(RealType bp) { |
| 274 |
+ |
frameData.bondPotential = bp; |
| 275 |
+ |
hasBondPotential_ = true; |
| 276 |
+ |
statData[Stats::BOND_POTENTIAL] = frameData.bondPotential; |
| 277 |
+ |
} |
| 278 |
+ |
|
| 279 |
+ |
void setBendPotential(RealType bp) { |
| 280 |
+ |
frameData.bendPotential = bp; |
| 281 |
+ |
hasBendPotential_ = true; |
| 282 |
+ |
statData[Stats::BEND_POTENTIAL] = frameData.bendPotential; |
| 283 |
+ |
} |
| 284 |
+ |
|
| 285 |
+ |
void setTorsionPotential(RealType tp) { |
| 286 |
+ |
frameData.torsionPotential = tp; |
| 287 |
+ |
hasTorsionPotential_ = true; |
| 288 |
+ |
statData[Stats::DIHEDRAL_POTENTIAL] = frameData.torsionPotential; |
| 289 |
+ |
} |
| 290 |
+ |
|
| 291 |
+ |
void setInversionPotential(RealType ip) { |
| 292 |
+ |
frameData.inversionPotential = ip; |
| 293 |
+ |
hasInversionPotential_ = true; |
| 294 |
+ |
statData[Stats::INVERSION_POTENTIAL] = frameData.inversionPotential; |
| 295 |
+ |
} |
| 296 |
+ |
|
| 297 |
+ |
void setLongRangePotential(RealType lrp) { |
| 298 |
+ |
frameData.longRangePotential = lrp; |
| 299 |
+ |
hasLongRangePotential_ = true; |
| 300 |
+ |
statData[Stats::LONG_RANGE_POTENTIAL] = frameData.longRangePotential; |
| 301 |
+ |
} |
| 302 |
+ |
|
| 303 |
+ |
RealType getLongRangePotential() { |
| 304 |
+ |
return frameData.longRangePotential; |
| 305 |
+ |
} |
| 306 |
+ |
|
| 307 |
+ |
void setLongRangePotentialFamilies(potVec lrPot) { |
| 308 |
+ |
frameData.lrPotentials = lrPot; |
| 309 |
+ |
hasLongRangePotentialFamilies_ = true; |
| 310 |
+ |
statData[Stats::VANDERWAALS_POTENTIAL] = frameData.lrPotentials[VANDERWAALS_FAMILY]; |
| 311 |
+ |
statData[Stats::ELECTROSTATIC_POTENTIAL] = frameData.lrPotentials[ELECTROSTATIC_FAMILY]; |
| 312 |
+ |
statData[Stats::METALLIC_POTENTIAL] = frameData.lrPotentials[METALLIC_FAMILY]; |
| 313 |
+ |
statData[Stats::HYDROGENBONDING_POTENTIAL] = frameData.lrPotentials[HYDROGENBONDING_FAMILY]; |
| 314 |
+ |
} |
| 315 |
+ |
|
| 316 |
+ |
potVec getLongRangePotentials() { |
| 317 |
+ |
return frameData.lrPotentials; |
| 318 |
+ |
} |
| 319 |
+ |
|
| 320 |
+ |
void setExcludedPotentials(potVec exPot) { |
| 321 |
+ |
frameData.excludedPotentials = exPot; |
| 322 |
+ |
hasExcludedPotentials_ = true; |
| 323 |
+ |
} |
| 324 |
+ |
|
| 325 |
+ |
potVec getExcludedPotentials() { |
| 326 |
+ |
return frameData.excludedPotentials; |
| 327 |
+ |
} |
| 328 |
+ |
|
| 329 |
+ |
|
| 330 |
+ |
void setRestraintPotential(RealType rp) { |
| 331 |
+ |
frameData.restraintPotential = rp; |
| 332 |
+ |
hasRestraintPotential_ = true; |
| 333 |
+ |
statData[Stats::RESTRAINT_POTENTIAL] = frameData.restraintPotential; |
| 334 |
+ |
} |
| 335 |
+ |
|
| 336 |
+ |
RealType getRestraintPotential() { |
| 337 |
+ |
return frameData.restraintPotential; |
| 338 |
+ |
} |
| 339 |
+ |
|
| 340 |
+ |
void setRawPotential(RealType rp) { |
| 341 |
+ |
frameData.rawPotential = rp; |
| 342 |
+ |
hasRawPotential_ = true; |
| 343 |
+ |
statData[Stats::RAW_POTENTIAL] = frameData.rawPotential; |
| 344 |
+ |
} |
| 345 |
+ |
|
| 346 |
+ |
RealType getRawPotential() { |
| 347 |
+ |
return frameData.rawPotential; |
| 348 |
+ |
} |
| 349 |
+ |
|
| 350 |
|
RealType getChi() { |
| 351 |
< |
return chi_; |
| 351 |
> |
return frameData.chi; |
| 352 |
|
} |
| 353 |
|
|
| 354 |
|
void setChi(RealType chi) { |
| 355 |
< |
chi_ = chi; |
| 355 |
> |
frameData.chi = chi; |
| 356 |
|
} |
| 357 |
|
|
| 358 |
|
RealType getIntegralOfChiDt() { |
| 359 |
< |
return integralOfChiDt_; |
| 359 |
> |
return frameData.integralOfChiDt; |
| 360 |
|
} |
| 361 |
|
|
| 362 |
|
void setIntegralOfChiDt(RealType integralOfChiDt) { |
| 363 |
< |
integralOfChiDt_ = integralOfChiDt; |
| 363 |
> |
frameData.integralOfChiDt = integralOfChiDt; |
| 364 |
|
} |
| 365 |
|
|
| 366 |
+ |
RealType getChiElectronic() { |
| 367 |
+ |
return frameData.chiQ; |
| 368 |
+ |
} |
| 369 |
+ |
|
| 370 |
+ |
void setChiElectronic(RealType chiQ) { |
| 371 |
+ |
frameData.chiQ = chiQ; |
| 372 |
+ |
} |
| 373 |
+ |
|
| 374 |
+ |
RealType getIntegralOfChiElectronicDt() { |
| 375 |
+ |
return frameData.integralOfChiQDt; |
| 376 |
+ |
} |
| 377 |
+ |
|
| 378 |
+ |
void setIntegralOfChiElectronicDt(RealType integralOfChiQDt) { |
| 379 |
+ |
frameData.integralOfChiQDt = integralOfChiQDt; |
| 380 |
+ |
} |
| 381 |
+ |
|
| 382 |
+ |
void setOrthoTolerance(RealType orthoTolerance) { |
| 383 |
+ |
orthoTolerance_ = orthoTolerance; |
| 384 |
+ |
} |
| 385 |
+ |
|
| 386 |
|
Mat3x3d getEta() { |
| 387 |
< |
return eta_; |
| 387 |
> |
return frameData.eta; |
| 388 |
|
} |
| 389 |
|
|
| 390 |
|
void setEta(const Mat3x3d& eta) { |
| 391 |
< |
eta_ = eta; |
| 391 |
> |
frameData.eta = eta; |
| 392 |
|
} |
| 393 |
< |
|
| 393 |
> |
|
| 394 |
> |
Mat3x3d getStressTensor() { |
| 395 |
> |
return frameData.stressTensor; |
| 396 |
> |
} |
| 397 |
> |
|
| 398 |
> |
void setStressTensor(const Mat3x3d& stressTensor) { |
| 399 |
> |
frameData.stressTensor = stressTensor; |
| 400 |
> |
} |
| 401 |
> |
|
| 402 |
> |
Vector3d getConductiveHeatFlux() { |
| 403 |
> |
return frameData.conductiveHeatFlux; |
| 404 |
> |
} |
| 405 |
> |
|
| 406 |
> |
void setConductiveHeatFlux(const Vector3d& heatFlux) { |
| 407 |
> |
frameData.conductiveHeatFlux = heatFlux; |
| 408 |
> |
} |
| 409 |
> |
|
| 410 |
> |
bool hasCOM() { |
| 411 |
> |
return hasCOM_; |
| 412 |
> |
} |
| 413 |
> |
|
| 414 |
> |
void setCOMprops(const Vector3d& COM, const Vector3d& COMvel, const Vector3d& COMw) { |
| 415 |
> |
frameData.COM = COM; |
| 416 |
> |
frameData.COMvel = COMvel; |
| 417 |
> |
frameData.COMw = COMw; |
| 418 |
> |
hasCOM_ = true; |
| 419 |
> |
} |
| 420 |
> |
|
| 421 |
|
DataStorage atomData; |
| 422 |
|
DataStorage rigidbodyData; |
| 423 |
+ |
DataStorage cgData; |
| 424 |
+ |
FrameData frameData; |
| 425 |
|
Stats statData; |
| 168 |
– |
|
| 169 |
– |
private: |
| 170 |
– |
RealType currentTime_; |
| 426 |
|
|
| 427 |
< |
Mat3x3d hmat_; |
| 428 |
< |
Mat3x3d invHmat_; |
| 429 |
< |
int orthoRhombic_; |
| 430 |
< |
|
| 431 |
< |
RealType chi_; |
| 432 |
< |
RealType integralOfChiDt_; |
| 433 |
< |
Mat3x3d eta_; |
| 434 |
< |
|
| 435 |
< |
int id_; /**< identification number of the snapshot */ |
| 427 |
> |
private: |
| 428 |
> |
RealType orthoTolerance_; |
| 429 |
> |
bool hasCOM_; |
| 430 |
> |
bool hasVolume_; |
| 431 |
> |
bool hasShortRangePotential_; |
| 432 |
> |
bool hasBondPotential_; |
| 433 |
> |
bool hasBendPotential_; |
| 434 |
> |
bool hasTorsionPotential_; |
| 435 |
> |
bool hasInversionPotential_; |
| 436 |
> |
bool hasLongRangePotential_; |
| 437 |
> |
bool hasLongRangePotentialFamilies_; |
| 438 |
> |
bool hasRestraintPotential_; |
| 439 |
> |
bool hasRawPotential_; |
| 440 |
> |
bool hasExcludedPotentials_; |
| 441 |
|
}; |
| 442 |
|
|
| 443 |
|
typedef DataStorage (Snapshot::*DataStoragePointer); |