SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
Public Member Functions | Private Attributes | List of all members
SSAGES::Snapshot Class Reference

Class containing a snapshot of the current simulation in time. More...

#include <Snapshot.h>

Public Member Functions

 Snapshot (const MPI_Comm &comm, uint wid)
 Constructor. More...
 
int GetIteration () const
 Get the current iteration. More...
 
int GetTargetIterations () const
 Get target iterations. More...
 
double GetTemperature () const
 Get current temperature. More...
 
double GetEnergy () const
 Get per-particle energy. More...
 
const Matrix3GetHMatrix () const
 Get system H-matrix. More...
 
const Matrix3GetVirial () const
 Get box virial. More...
 
Matrix3GetVirial ()
 Get box virial. More...
 
const Vector3GetOrigin () const
 Get origin of the system. More...
 
const Bool3IsPeriodic () const
 Get periodicity of three dimensions. More...
 
double GetVolume () const
 Get system volume. More...
 
double GetKb () const
 Get system Kb. More...
 
double GetDielectric () const
 Get dielectric constant. More...
 
double Getqqrd2e () const
 Get qqrd2e value. More...
 
const mxx::comm & GetCommunicator () const
 Get communicator for walker. More...
 
unsigned GetWalkerID () const
 Get walker ID. More...
 
unsigned GetNumAtoms () const
 Get number of atoms in this snapshot. More...
 
void SetIteration (int iteration)
 Set the iteration. More...
 
void SetTargetIterations (int target)
 Set target iterations. More...
 
void SetTemperature (double temperature)
 Change the temperature. More...
 
void SetEnergy (double energy)
 Change the energy. More...
 
void SetHMatrix (const Matrix3 &hmat)
 Change the Box H-matrix. More...
 
void SetVirial (const Matrix3 &virial)
 Change the box virial. More...
 
void SetOrigin (const Vector3 &origin)
 Change the box origin. More...
 
void SetPeriodicity (const Bool3 &isperiodic)
 Change the periodicity of the system. More...
 
void SetKb (double kb)
 Change the kb. More...
 
void SetDielectric (double dielectric)
 Set the dielectric constant. More...
 
void Setqqrd2e (double qqrd2e)
 Set the value for qqrd2e. More...
 
void SetNumAtoms (unsigned int natoms)
 Set number of atoms in this snapshot. More...
 
const std::vector< Vector3 > & GetPositions () const
 Access the particle positions. More...
 
std::vector< Vector3 > & GetPositions ()
 Access the particle positions. More...
 
const std::vector< Integer3 > & GetImageFlags () const
 Access the particles image flags. More...
 
std::vector< Integer3 > & GetImageFlags ()
 Access the particles image flags. More...
 
const std::vector< Vector3 > & GetVelocities () const
 Access the particle velocities. More...
 
std::vector< Vector3 > & GetVelocities ()
 Access the particle velocities. More...
 
const std::vector< Vector3 > & GetForces () const
 Access the per-particle forces. More...
 
std::vector< Vector3 > & GetForces ()
 Access the per-particle forces. More...
 
const std::vector< double > & GetMasses () const
 Const access to the particle masses. More...
 
std::vector< double > & GetMasses ()
 Const access to the particle masses. More...
 
Vector3 ScaleVector (const Vector3 &v) const
 Scale a vector into fractional coordinates. More...
 
Vector3 UnwrapVector (const Vector3 &v, const Integer3 &image) const
 Unwrap a vector's real coordinates according to its image replica count. More...
 
void ApplyMinimumImage (Vector3 *v) const
 Apply minimum image to a vector. More...
 
Vector3 ApplyMinimumImage (const Vector3 &v) const
 Apply minimum image to a vector. More...
 
double TotalMass (const Label &indices) const
 Compute the total mass of a group of particles based on index. More...
 
Vector3 CenterOfMass (const Label &indices) const
 
Vector3 CenterOfMass (const Label &indices, double mtot) const
 
const LabelGetAtomIDs () const
 Access the atom IDs. More...
 
LabelGetAtomIDs ()
 Access the atom IDs. More...
 
int GetLocalIndex (int id) const
 Gets the local atom index corresponding to an atom ID. More...
 
void GetLocalIndices (const Label &ids, Label *indices) const
 
const std::vector< double > & GetCharges () const
 Access the atom charges. More...
 
std::vector< double > & GetCharges ()
 Access the atom charges. More...
 
const LabelGetAtomTypes () const
 Access the atom types. More...
 
LabelGetAtomTypes ()
 Access the atom types. More...
 
const std::vector< std::vector
< double > > & 
GetSigmas () const
 Access the atom sigmas. More...
 
std::vector< std::vector
< double > > & 
GetSigmas ()
 Access the atom sigmas. More...
 
const std::string & GetSnapshotID () const
 Access the snapshot ID. More...
 
std::string & GetSnapshotID ()
 Access the snapshot ID. More...
 
bool HasChanged () const
 Query if Snapshot was modified. More...
 
void Changed (bool state)
 Set the "changed" flag of the Snapshot. More...
 
std::vector< double > SerializePositions ()
 Return the serialized positions across all local cores.
 
std::vector< double > SerializeVelocities ()
 Return the serialized velocities across all local cores.
 
std::vector< int > SerializeIDs ()
 Return the serialized positions across all local cores.
 
 ~Snapshot ()
 Destructor.
 

Private Attributes

mxx::comm comm_
 Local snapshot (walker) communicator.
 
uint wid_
 Walker ID.
 
uint nlocal_
 Number of atoms located on this snapshot.
 
std::string ID_
 ID string.
 
Matrix3 H_
 Parrinello-Rahman box H-matrix.
 
Matrix3 Hinv_
 Parinello-Rahman box inverse.
 
Matrix3 virial_
 
Vector3 origin_
 < Virial tensor. More...
 
Bool3 isperiodic_
 Periodicity of box.
 
std::vector< Vector3positions_
 Positions.
 
std::vector< Integer3images_
 Unwrapped positions.
 
std::vector< Vector3velocities_
 Velocities.
 
std::vector< Vector3forces_
 Forces.
 
std::vector< double > masses_
 Masses.
 
std::vector< double > charges_
 Charges.
 
Label atomids_
 List of Atom IDs.
 
Label types_
 List of Atom types.
 
std::vector< std::vector
< double > > 
sigma_
 Sigma.
 
int iteration_
 Iteration of Simulation.
 
int targetiter_
 Iteration target of simulation.
 
double temperature_
 Current temperature.
 
double energy_
 Average per-particle energy.
 
double kb_
 Kb from the MD driver.
 
double dielectric_
 Dielectric.
 
double qqrd2e_
 qqrd2e
 
bool changed_
 TRUE is Simulation state changed
 

Detailed Description

Class containing a snapshot of the current simulation in time.

This contains information on the particle positions, velocities, etc... and additional information on the state of the system.

Definition at line 43 of file Snapshot.h.

Constructor & Destructor Documentation

SSAGES::Snapshot::Snapshot ( const MPI_Comm &  comm,
uint  wid 
)
inline

Constructor.

Parameters
commMPI communicator
widWalker ID

Initialize a snapshot with MPI communicator and a correpsonding walker ID.

Definition at line 92 of file Snapshot.h.

92  :
93  comm_(comm), wid_(wid), H_(), Hinv_(), virial_(Matrix3::Zero()),
94  origin_({0,0,0}), isperiodic_({true, true, true}), positions_(0),
95  images_(0), velocities_(0), forces_(0), masses_(0), atomids_(0),
96  types_(0), iteration_(0), temperature_(0), energy_(0), kb_(0)
97  {}
std::vector< Integer3 > images_
Unwrapped positions.
Definition: Snapshot.h:64
int iteration_
Iteration of Simulation.
Definition: Snapshot.h:73
std::vector< Vector3 > forces_
Forces.
Definition: Snapshot.h:66
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54
double kb_
Kb from the MD driver.
Definition: Snapshot.h:77
std::vector< Vector3 > positions_
Positions.
Definition: Snapshot.h:63
Matrix3 Hinv_
Parinello-Rahman box inverse.
Definition: Snapshot.h:55
double energy_
Average per-particle energy.
Definition: Snapshot.h:76
double temperature_
Current temperature.
Definition: Snapshot.h:75
Bool3 isperiodic_
Periodicity of box.
Definition: Snapshot.h:61
mxx::comm comm_
Local snapshot (walker) communicator.
Definition: Snapshot.h:47
std::vector< double > masses_
Masses.
Definition: Snapshot.h:67
Label types_
List of Atom types.
Definition: Snapshot.h:70
std::vector< Vector3 > velocities_
Velocities.
Definition: Snapshot.h:65
uint wid_
Walker ID.
Definition: Snapshot.h:49
Vector3 origin_
< Virial tensor.
Definition: Snapshot.h:59
Label atomids_
List of Atom IDs.
Definition: Snapshot.h:69

Member Function Documentation

void SSAGES::Snapshot::ApplyMinimumImage ( Vector3 v) const
inline

Apply minimum image to a vector.

Parameters
vVector of interest

Definition at line 418 of file Snapshot.h.

References H_, Hinv_, isperiodic_, and SSAGES::roundf().

Referenced by CenterOfMass(), SSAGES::AngleCV::Evaluate(), SSAGES::TorsionalCV::Evaluate(), SSAGES::GyrationTensorCV::Evaluate(), SSAGES::ParticlePositionCV::Evaluate(), SSAGES::PairwiseCV::Evaluate(), SSAGES::RouseModeCV::Evaluate(), and SSAGES::ParticleSeparationCV::Evaluate().

419  {
420  Vector3 scaled = Hinv_*(*v);
421 
422  for(int i = 0; i < 3; ++i)
423  scaled[i] -= isperiodic_[i]*roundf(scaled[i]);
424 
425  *v = H_*scaled;
426  }
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54
Matrix3 Hinv_
Parinello-Rahman box inverse.
Definition: Snapshot.h:55
Bool3 isperiodic_
Periodicity of box.
Definition: Snapshot.h:61
double roundf(double x)
Quick helper function to round a double.
Definition: Snapshot.h:31
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33

Here is the call graph for this function:

Here is the caller graph for this function:

Vector3 SSAGES::Snapshot::ApplyMinimumImage ( const Vector3 v) const
inline

Apply minimum image to a vector.

Parameters
vVector of interest
Returns
Return new vector.

Definition at line 433 of file Snapshot.h.

References H_, Hinv_, isperiodic_, and SSAGES::roundf().

434  {
435  Vector3 scaled = Hinv_*v;
436 
437  for(int i = 0; i < 3; ++i)
438  scaled[i] -= isperiodic_[i]*roundf(scaled[i]);
439 
440  return H_*scaled;
441  }
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54
Matrix3 Hinv_
Parinello-Rahman box inverse.
Definition: Snapshot.h:55
Bool3 isperiodic_
Periodicity of box.
Definition: Snapshot.h:61
double roundf(double x)
Quick helper function to round a double.
Definition: Snapshot.h:31
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33

Here is the call graph for this function:

Vector3 SSAGES::Snapshot::CenterOfMass ( const Label indices) const
inline

Compute center of mass of a group of atoms based on idex with provided Total mass.

Parameters
indicesIDs of particles of interest.
Returns
Vector3 Center of mass of particles.

Definition at line 466 of file Snapshot.h.

References TotalMass().

Referenced by SSAGES::ParticleCoordinateCV::Evaluate(), SSAGES::GyrationTensorCV::Evaluate(), SSAGES::ParticlePositionCV::Evaluate(), SSAGES::RouseModeCV::Evaluate(), and SSAGES::ParticleSeparationCV::Evaluate().

467  {
468  // Get total mass.
469  auto mtot = TotalMass(indices);
470 
471  return CenterOfMass(indices, mtot);
472  }
double TotalMass(const Label &indices) const
Compute the total mass of a group of particles based on index.
Definition: Snapshot.h:451
Vector3 CenterOfMass(const Label &indices) const
Definition: Snapshot.h:466

Here is the call graph for this function:

Here is the caller graph for this function:

Vector3 SSAGES::Snapshot::CenterOfMass ( const Label indices,
double  mtot 
) const
inline

Compute center of mass of a group of atoms based on idex with provided Total mass.

Parameters
indicesIDs of particles of interest.
mtotTotal mass of particle group.
Returns
Vector3 Center of mass of particles.
Note
Each processor passes in the local indices of the atoms of interest and this function will collect the data and compute the center of mass.

Definition at line 483 of file Snapshot.h.

References ApplyMinimumImage(), comm_, masses_, and positions_.

484  {
485  // Store coorinates and masses in vectors to gather.
486  std::vector<double> pos, mass, gpos, gmass;
487  std::vector<int> pcounts(comm_.size(), 0), mcounts(comm_.size(), 0);
488  std::vector<int> pdispls(comm_.size()+1, 0), mdispls(comm_.size()+1, 0);
489 
490  pcounts[comm_.rank()] = 3*indices.size();
491  mcounts[comm_.rank()] = indices.size();
492 
493  // Reduce counts.
494  MPI_Allreduce(MPI_IN_PLACE, pcounts.data(), pcounts.size(), MPI_INT, MPI_SUM, comm_);
495  MPI_Allreduce(MPI_IN_PLACE, mcounts.data(), mcounts.size(), MPI_INT, MPI_SUM, comm_);
496 
497  // Compute displacements.
498  std::partial_sum(pcounts.begin(), pcounts.end(), pdispls.begin() + 1);
499  std::partial_sum(mcounts.begin(), mcounts.end(), mdispls.begin() + 1);
500 
501  // Fill up mass and position vectors.
502  for(auto& idx : indices)
503  {
504  auto& p = positions_[idx];
505  pos.push_back(p[0]);
506  pos.push_back(p[1]);
507  pos.push_back(p[2]);
508  mass.push_back(masses_[idx]);
509  }
510 
511  // Re-size receiving vectors.
512  gpos.resize(pdispls.back(), 0);
513  gmass.resize(mdispls.back(), 0);
514 
515  // All-gather data.
516  MPI_Allgatherv(pos.data(), pos.size(), MPI_DOUBLE, gpos.data(), pcounts.data(), pdispls.data(), MPI_DOUBLE, comm_);
517  MPI_Allgatherv(mass.data(), mass.size(), MPI_DOUBLE, gmass.data(), mcounts.data(), mdispls.data(), MPI_DOUBLE, comm_);
518 
519  // Loop through atoms and compute mass weighted sum.
520  // We march linearly through list and find nearest image
521  // to each successive particle to properly unwrap object.
522  Vector3 ppos = {gpos[0], gpos[1], gpos[2]}; // Previous unwrapped position.
523  Vector3 cpos = ppos;
524  Vector3 xcm = gmass[0]*cpos;
525 
526  for(size_t i = 1, j = 3; i < gmass.size(); ++i, j += 3)
527  {
528  cpos = {gpos[j], gpos[j+1], gpos[j+2]};
529  cpos = ApplyMinimumImage(cpos - ppos) + ppos;
530  xcm += gmass[i]*cpos;
531  ppos = cpos;
532  }
533 
534  return xcm/mtot;
535  }
std::vector< Vector3 > positions_
Positions.
Definition: Snapshot.h:63
mxx::comm comm_
Local snapshot (walker) communicator.
Definition: Snapshot.h:47
std::vector< double > masses_
Masses.
Definition: Snapshot.h:67
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33
void ApplyMinimumImage(Vector3 *v) const
Apply minimum image to a vector.
Definition: Snapshot.h:418

Here is the call graph for this function:

void SSAGES::Snapshot::Changed ( bool  state)
inline

Set the "changed" flag of the Snapshot.

Parameters
stateState to which the "changed" flag is set

Definition at line 645 of file Snapshot.h.

References changed_.

Referenced by SSAGES::Hook::PostIntegrationHook(), SSAGES::Hook::PostSimulationHook(), and SSAGES::Hook::PreSimulationHook().

645 { changed_ = state; }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81

Here is the caller graph for this function:

const Label& SSAGES::Snapshot::GetAtomIDs ( ) const
inline

Access the atom IDs.

Returns
List of atom IDs

Definition at line 541 of file Snapshot.h.

References atomids_.

Referenced by SSAGES::ForwardFlux::AppendTrajectoryFile(), SSAGES::PairwiseCV::Evaluate(), SSAGES::FiniteTempString::PostIntegration(), SSAGES::Swarm::PostIntegration(), SSAGES::ForwardFlux::ReadFFSConfiguration(), and SSAGES::ForwardFlux::WriteFFSConfiguration().

541 { return atomids_; }
Label atomids_
List of Atom IDs.
Definition: Snapshot.h:69

Here is the caller graph for this function:

Label& SSAGES::Snapshot::GetAtomIDs ( )
inline

Access the atom IDs.

Returns
List of atom IDs

Definition at line 544 of file Snapshot.h.

References atomids_, and changed_.

545  {
546  changed_ = true;
547  return atomids_;
548  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
Label atomids_
List of Atom IDs.
Definition: Snapshot.h:69
const Label& SSAGES::Snapshot::GetAtomTypes ( ) const
inline

Access the atom types.

Returns
List of atom types

Definition at line 600 of file Snapshot.h.

References types_.

Referenced by SSAGES::COPSSImage::force_pol(), and SSAGES::COPSSImage::PostIntegration().

600 { return types_; }
Label types_
List of Atom types.
Definition: Snapshot.h:70

Here is the caller graph for this function:

Label& SSAGES::Snapshot::GetAtomTypes ( )
inline

Access the atom types.

Returns
List of atom types

Definition at line 603 of file Snapshot.h.

References changed_, and types_.

604  {
605  changed_ = true;
606  return types_;
607  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
Label types_
List of Atom types.
Definition: Snapshot.h:70
const std::vector<double>& SSAGES::Snapshot::GetCharges ( ) const
inline

Access the atom charges.

Returns
List of atom charges

Definition at line 587 of file Snapshot.h.

References charges_.

Referenced by SSAGES::COPSSImage::force_pol(), and SSAGES::COPSSImage::PostIntegration().

587 { return charges_; }
std::vector< double > charges_
Charges.
Definition: Snapshot.h:68

Here is the caller graph for this function:

std::vector<double>& SSAGES::Snapshot::GetCharges ( )
inline

Access the atom charges.

Returns
List of atom charges

Definition at line 590 of file Snapshot.h.

References changed_, and charges_.

591  {
592  changed_ = true;
593  return charges_;
594  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< double > charges_
Charges.
Definition: Snapshot.h:68
const mxx::comm& SSAGES::Snapshot::GetCommunicator ( ) const
inline

Get communicator for walker.

Returns
Communicator containing processors for simulation instance (walker).

Access the communicator containing the set of processors that belong to a single simulation box.

Definition at line 184 of file Snapshot.h.

References comm_.

Referenced by SSAGES::BoxVolumeCV::Evaluate(), SSAGES::AngleCV::Evaluate(), SSAGES::TorsionalCV::Evaluate(), SSAGES::GyrationTensorCV::Evaluate(), SSAGES::PairwiseCV::Evaluate(), SSAGES::AngleCV::Initialize(), SSAGES::PairwiseCV::Initialize(), SSAGES::ParticleCoordinateCV::Initialize(), SSAGES::TorsionalCV::Initialize(), SSAGES::GyrationTensorCV::Initialize(), SSAGES::ParticlePositionCV::Initialize(), SSAGES::RouseModeCV::Initialize(), SSAGES::ParticleSeparationCV::Initialize(), and SSAGES::COPSS::PostIntegration().

185  {
186  return comm_;
187  }
mxx::comm comm_
Local snapshot (walker) communicator.
Definition: Snapshot.h:47

Here is the caller graph for this function:

double SSAGES::Snapshot::GetDielectric ( ) const
inline

Get dielectric constant.

Returns
dielectric constant.

Definition at line 169 of file Snapshot.h.

References dielectric_.

Referenced by SSAGES::COPSSImage::force_pol().

169 { return dielectric_; }
double dielectric_
Dielectric.
Definition: Snapshot.h:79

Here is the caller graph for this function:

double SSAGES::Snapshot::GetEnergy ( ) const
inline

Get per-particle energy.

Returns
Current average per-particle energy

Definition at line 121 of file Snapshot.h.

References energy_.

121 { return energy_; }
double energy_
Average per-particle energy.
Definition: Snapshot.h:76
const std::vector<Vector3>& SSAGES::Snapshot::GetForces ( ) const
inline
std::vector<Vector3>& SSAGES::Snapshot::GetForces ( )
inline

Access the per-particle forces.

Returns
List of per-particle forces

Definition at line 365 of file Snapshot.h.

References changed_, and forces_.

366  {
367  changed_ = true;
368  return forces_;
369  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< Vector3 > forces_
Forces.
Definition: Snapshot.h:66
const Matrix3& SSAGES::Snapshot::GetHMatrix ( ) const
inline

Get system H-matrix.

Returns
Parrinello-Rahman H-matrix of simulation box

Definition at line 127 of file Snapshot.h.

References H_.

Referenced by SSAGES::Umbrella::PostIntegration().

127 { return H_; }
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54

Here is the caller graph for this function:

const std::vector<Integer3>& SSAGES::Snapshot::GetImageFlags ( ) const
inline

Access the particles image flags.

Returns
List of particle image flags

Definition at line 336 of file Snapshot.h.

References images_.

336 { return images_; }
std::vector< Integer3 > images_
Unwrapped positions.
Definition: Snapshot.h:64
std::vector<Integer3>& SSAGES::Snapshot::GetImageFlags ( )
inline

Access the particles image flags.

Returns
List of particle image flags

Definition at line 339 of file Snapshot.h.

References changed_, and images_.

340  {
341  changed_ = true;
342  return images_;
343  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< Integer3 > images_
Unwrapped positions.
Definition: Snapshot.h:64
int SSAGES::Snapshot::GetIteration ( ) const
inline

Get the current iteration.

Returns
Current Iteration

Definition at line 103 of file Snapshot.h.

References iteration_.

Referenced by SSAGES::Logger::PostIntegration(), SSAGES::Umbrella::PostIntegration(), SSAGES::Meta::PostIntegration(), SSAGES::Basis::PostIntegration(), and SSAGES::Hook::PostIntegrationHook().

103 {return iteration_; }
int iteration_
Iteration of Simulation.
Definition: Snapshot.h:73

Here is the caller graph for this function:

double SSAGES::Snapshot::GetKb ( ) const
inline

Get system Kb.

Returns
Kb of the current simulation box

Definition at line 163 of file Snapshot.h.

References kb_.

Referenced by SSAGES::Basis::PostIntegration().

163 { return kb_; }
double kb_
Kb from the MD driver.
Definition: Snapshot.h:77

Here is the caller graph for this function:

int SSAGES::Snapshot::GetLocalIndex ( int  id) const
inline

Gets the local atom index corresponding to an atom ID.

Parameters
idAtom ID.
Returns
Local atom index or -1 if not found.

Definition at line 555 of file Snapshot.h.

References atomids_.

Referenced by SSAGES::AngleCV::Evaluate(), SSAGES::TorsionalCV::Evaluate(), SSAGES::PairwiseCV::Evaluate(), GetLocalIndices(), SSAGES::PairwiseCV::Initialize(), SSAGES::ParticleCoordinateCV::Initialize(), SSAGES::GyrationTensorCV::Initialize(), SSAGES::ParticlePositionCV::Initialize(), SSAGES::RouseModeCV::Initialize(), SSAGES::ParticleSeparationCV::Initialize(), SSAGES::FiniteTempString::PostIntegration(), and SSAGES::Swarm::PostIntegration().

556  {
557 
558  auto s = std::find(atomids_.begin(), atomids_.end(), id);
559  if(s == atomids_.end())
560  return -1;
561  else
562  return s - atomids_.begin();
563  }
Label atomids_
List of Atom IDs.
Definition: Snapshot.h:69

Here is the caller graph for this function:

void SSAGES::Snapshot::GetLocalIndices ( const Label ids,
Label indices 
) const
inline

Gets the local atom indices corresponding to atom IDs in the vector.

Parameters
idsVector of atom ID's.
indicesPointer to container for local atom indices.
Note
If atom does not exist on processor, it will be ignored.

Definition at line 573 of file Snapshot.h.

References GetLocalIndex().

Referenced by SSAGES::ParticleCoordinateCV::Evaluate(), SSAGES::GyrationTensorCV::Evaluate(), SSAGES::ParticlePositionCV::Evaluate(), SSAGES::PairwiseCV::Evaluate(), SSAGES::RouseModeCV::Evaluate(), SSAGES::ParticleSeparationCV::Evaluate(), SSAGES::AlphaRMSDCV::Evaluate(), SSAGES::AntiBetaRMSDCV::Evaluate(), SSAGES::ParallelBetaRMSDCV::Evaluate(), SSAGES::AngleCV::Initialize(), SSAGES::TorsionalCV::Initialize(), and SSAGES::RouseModeCV::setMasses().

574  {
575  for(auto& id : ids)
576  {
577  auto idx = GetLocalIndex(id);
578  if(idx != -1)
579  indices->push_back(idx);
580  }
581  }
int GetLocalIndex(int id) const
Gets the local atom index corresponding to an atom ID.
Definition: Snapshot.h:555

Here is the call graph for this function:

Here is the caller graph for this function:

const std::vector<double>& SSAGES::Snapshot::GetMasses ( ) const
inline

Const access to the particle masses.

Returns
List of Masses

Note that the Masses can be either stored as per-atom or per-type depending on the Lammps Atom type used.

Definition at line 378 of file Snapshot.h.

References masses_.

Referenced by SSAGES::ParticleCoordinateCV::Evaluate(), SSAGES::GyrationTensorCV::Evaluate(), SSAGES::ParticlePositionCV::Evaluate(), SSAGES::RouseModeCV::Evaluate(), SSAGES::ParticleSeparationCV::Evaluate(), and SSAGES::ABF::PostIntegration().

378 { return masses_; }
std::vector< double > masses_
Masses.
Definition: Snapshot.h:67

Here is the caller graph for this function:

std::vector<double>& SSAGES::Snapshot::GetMasses ( )
inline

Const access to the particle masses.

Returns
List of Masses

Note that the Masses can be either stored as per-atom or per-type depending on the Lammps Atom type used.

Definition at line 381 of file Snapshot.h.

References changed_, and masses_.

382  {
383  changed_ = true;
384  return masses_;
385  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< double > masses_
Masses.
Definition: Snapshot.h:67
unsigned SSAGES::Snapshot::GetNumAtoms ( ) const
inline
const Vector3& SSAGES::Snapshot::GetOrigin ( ) const
inline

Get origin of the system.

Returns
Vector containing coordinates of box origin.

Definition at line 145 of file Snapshot.h.

References origin_.

145 { return origin_; }
Vector3 origin_
< Virial tensor.
Definition: Snapshot.h:59
const std::vector<Vector3>& SSAGES::Snapshot::GetPositions ( ) const
inline
std::vector<Vector3>& SSAGES::Snapshot::GetPositions ( )
inline

Access the particle positions.

Returns
List of particle positions

Definition at line 326 of file Snapshot.h.

References changed_, and positions_.

327  {
328  changed_ = true;
329  return positions_;
330  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< Vector3 > positions_
Positions.
Definition: Snapshot.h:63
double SSAGES::Snapshot::Getqqrd2e ( ) const
inline

Get qqrd2e value.

Returns
Value of qqrd2e.

Definition at line 175 of file Snapshot.h.

References qqrd2e_.

Referenced by SSAGES::COPSSImage::force_pol().

175 { return qqrd2e_; }
double qqrd2e_
qqrd2e
Definition: Snapshot.h:80

Here is the caller graph for this function:

const std::vector<std::vector<double> >& SSAGES::Snapshot::GetSigmas ( ) const
inline

Access the atom sigmas.

Returns
List of atom sigmas

Definition at line 613 of file Snapshot.h.

References sigma_.

613 { return sigma_; }
std::vector< std::vector< double > > sigma_
Sigma.
Definition: Snapshot.h:71
std::vector<std::vector<double> >& SSAGES::Snapshot::GetSigmas ( )
inline

Access the atom sigmas.

Returns
List of atom sigmas

Definition at line 616 of file Snapshot.h.

References changed_, and sigma_.

617  {
618  changed_ = true;
619  return sigma_;
620  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< std::vector< double > > sigma_
Sigma.
Definition: Snapshot.h:71
const std::string& SSAGES::Snapshot::GetSnapshotID ( ) const
inline

Access the snapshot ID.

Returns
Snapshot ID

Definition at line 626 of file Snapshot.h.

References ID_.

626 { return ID_; }
std::string ID_
ID string.
Definition: Snapshot.h:52
std::string& SSAGES::Snapshot::GetSnapshotID ( )
inline

Access the snapshot ID.

Returns
Snapshot ID

Definition at line 629 of file Snapshot.h.

References changed_, and ID_.

630  {
631  changed_ = true;
632  return ID_;
633  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::string ID_
ID string.
Definition: Snapshot.h:52
int SSAGES::Snapshot::GetTargetIterations ( ) const
inline

Get target iterations.

Returns
Target Iterations

Definition at line 109 of file Snapshot.h.

References targetiter_.

Referenced by SSAGES::Meta::PreSimulation().

109 {return targetiter_; }
int targetiter_
Iteration target of simulation.
Definition: Snapshot.h:74

Here is the caller graph for this function:

double SSAGES::Snapshot::GetTemperature ( ) const
inline

Get current temperature.

Returns
System temperature

Definition at line 115 of file Snapshot.h.

References temperature_.

Referenced by SSAGES::Basis::PostIntegration().

115 {return temperature_; }
double temperature_
Current temperature.
Definition: Snapshot.h:75

Here is the caller graph for this function:

const std::vector<Vector3>& SSAGES::Snapshot::GetVelocities ( ) const
inline

Access the particle velocities.

Returns
List of particle velocities

Definition at line 349 of file Snapshot.h.

References velocities_.

Referenced by SSAGES::ForwardFlux::AppendTrajectoryFile(), SSAGES::Swarm::PostIntegration(), SSAGES::ABF::PostIntegration(), SSAGES::ForwardFlux::ReadFFSConfiguration(), and SSAGES::ForwardFlux::WriteFFSConfiguration().

349 { return velocities_; }
std::vector< Vector3 > velocities_
Velocities.
Definition: Snapshot.h:65

Here is the caller graph for this function:

std::vector<Vector3>& SSAGES::Snapshot::GetVelocities ( )
inline

Access the particle velocities.

Returns
List of particle velocities

Definition at line 352 of file Snapshot.h.

References changed_, and velocities_.

353  {
354  changed_ = true;
355  return velocities_;
356  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
std::vector< Vector3 > velocities_
Velocities.
Definition: Snapshot.h:65
const Matrix3& SSAGES::Snapshot::GetVirial ( ) const
inline

Get box virial.

Virial tensor of the simulation box.

Definition at line 133 of file Snapshot.h.

Referenced by SSAGES::Umbrella::PostIntegration(), SSAGES::Meta::PostIntegration(), and SSAGES::Basis::PostIntegration().

133 { return virial_; }

Here is the caller graph for this function:

Matrix3& SSAGES::Snapshot::GetVirial ( )
inline

Get box virial.

Virial tensor of the simulation box.

Definition at line 139 of file Snapshot.h.

139 { return virial_; }
double SSAGES::Snapshot::GetVolume ( ) const
inline

Get system volume.

Returns
Volume of the current simulation box

Definition at line 157 of file Snapshot.h.

References H_.

Referenced by SSAGES::BoxVolumeCV::Evaluate().

157 { return H_.determinant(); }
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54

Here is the caller graph for this function:

unsigned SSAGES::Snapshot::GetWalkerID ( ) const
inline

Get walker ID.

Returns
ID of the Walker

Definition at line 193 of file Snapshot.h.

References wid_.

Referenced by SSAGES::COPSS::PostIntegration(), SSAGES::StringMethod::PreSimulation(), and SSAGES::Basis::PreSimulation().

193 { return wid_; }
uint wid_
Walker ID.
Definition: Snapshot.h:49

Here is the caller graph for this function:

bool SSAGES::Snapshot::HasChanged ( ) const
inline

Query if Snapshot was modified.

Returns
True if Snapshot was modified, else return False

Definition at line 639 of file Snapshot.h.

References changed_.

Referenced by SSAGES::Hook::PostIntegrationHook(), SSAGES::Hook::PostSimulationHook(), and SSAGES::Hook::PreSimulationHook().

639 { return changed_; }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81

Here is the caller graph for this function:

const Bool3& SSAGES::Snapshot::IsPeriodic ( ) const
inline

Get periodicity of three dimensions.

Returns
Three dimensional boolean containing periodicity of each dimension.

Definition at line 151 of file Snapshot.h.

References isperiodic_.

151 {return isperiodic_; }
Bool3 isperiodic_
Periodicity of box.
Definition: Snapshot.h:61
Vector3 SSAGES::Snapshot::ScaleVector ( const Vector3 v) const
inline

Scale a vector into fractional coordinates.

Parameters
vVector of interest
Returns
Scaled vector in fractional coordinates

Definition at line 392 of file Snapshot.h.

References Hinv_, and origin_.

393  {
394  return Hinv_*(v-origin_);
395  }
Matrix3 Hinv_
Parinello-Rahman box inverse.
Definition: Snapshot.h:55
Vector3 origin_
< Virial tensor.
Definition: Snapshot.h:59
void SSAGES::Snapshot::SetDielectric ( double  dielectric)
inline

Set the dielectric constant.

Parameters
dielectricValue for the dielectric constant.

Definition at line 297 of file Snapshot.h.

References changed_, and dielectric_.

298  {
299  dielectric_ = dielectric;
300  changed_ = true;
301  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
double dielectric_
Dielectric.
Definition: Snapshot.h:79
void SSAGES::Snapshot::SetEnergy ( double  energy)
inline

Change the energy.

Parameters
energyNew value for the energy

Definition at line 236 of file Snapshot.h.

References changed_, and energy_.

237  {
238  energy_ = energy;
239  changed_ = true;
240  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
double energy_
Average per-particle energy.
Definition: Snapshot.h:76
void SSAGES::Snapshot::SetHMatrix ( const Matrix3 hmat)
inline

Change the Box H-matrix.

Parameters
hmatNew H-matrix for the system

Definition at line 246 of file Snapshot.h.

References changed_, H_, and Hinv_.

247  {
248  H_ = hmat;
249  Hinv_ = hmat.inverse();
250  changed_ = true;
251  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54
Matrix3 Hinv_
Parinello-Rahman box inverse.
Definition: Snapshot.h:55
void SSAGES::Snapshot::SetIteration ( int  iteration)
inline

Set the iteration.

Parameters
iterationNew value for the iteration

Definition at line 206 of file Snapshot.h.

References changed_, and iteration_.

207  {
208  iteration_ = iteration;
209  changed_ = true;
210  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
int iteration_
Iteration of Simulation.
Definition: Snapshot.h:73
void SSAGES::Snapshot::SetKb ( double  kb)
inline

Change the kb.

Parameters
kbNew value for the kb

Definition at line 287 of file Snapshot.h.

References changed_, and kb_.

288  {
289  kb_ = kb;
290  changed_ = true;
291  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
double kb_
Kb from the MD driver.
Definition: Snapshot.h:77
void SSAGES::Snapshot::SetNumAtoms ( unsigned int  natoms)
inline

Set number of atoms in this snapshot.

Parameters
natomsNumber of atoms in this snapshot

Definition at line 317 of file Snapshot.h.

References nlocal_.

317 { nlocal_ = natoms; }
uint nlocal_
Number of atoms located on this snapshot.
Definition: Snapshot.h:50
void SSAGES::Snapshot::SetOrigin ( const Vector3 origin)
inline

Change the box origin.

Parameters
originNew origin for the system

Definition at line 267 of file Snapshot.h.

References changed_, and origin_.

268  {
269  origin_ = origin;
270  changed_ = true;
271  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
Vector3 origin_
< Virial tensor.
Definition: Snapshot.h:59
void SSAGES::Snapshot::SetPeriodicity ( const Bool3 isperiodic)
inline

Change the periodicity of the system.

Parameters
isperiodicPeriodicity of three dimensions

Definition at line 277 of file Snapshot.h.

References changed_, and isperiodic_.

278  {
279  isperiodic_ = isperiodic;
280  changed_ = true;
281  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
Bool3 isperiodic_
Periodicity of box.
Definition: Snapshot.h:61
void SSAGES::Snapshot::Setqqrd2e ( double  qqrd2e)
inline

Set the value for qqrd2e.

Parameters
qqrd2eValue for qqrd2e.

Definition at line 307 of file Snapshot.h.

References changed_, and qqrd2e_.

308  {
309  qqrd2e_ = qqrd2e;
310  changed_ = true;
311  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
double qqrd2e_
qqrd2e
Definition: Snapshot.h:80
void SSAGES::Snapshot::SetTargetIterations ( int  target)
inline

Set target iterations.

Parameters
targetNew value for target iterations

Definition at line 216 of file Snapshot.h.

References changed_, and targetiter_.

217  {
218  targetiter_ = target;
219  changed_ = true;
220  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
int targetiter_
Iteration target of simulation.
Definition: Snapshot.h:74
void SSAGES::Snapshot::SetTemperature ( double  temperature)
inline

Change the temperature.

Parameters
temperatureNew value for the temperature

Definition at line 226 of file Snapshot.h.

References changed_, and temperature_.

227  {
228  temperature_ = temperature;
229  changed_ = true;
230  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
double temperature_
Current temperature.
Definition: Snapshot.h:75
void SSAGES::Snapshot::SetVirial ( const Matrix3 virial)
inline

Change the box virial.

Parameters
virialNew virial tensor for the system.

Definition at line 257 of file Snapshot.h.

References changed_.

258  {
259  virial_ = virial;
260  changed_ = true;
261  }
bool changed_
TRUE is Simulation state changed
Definition: Snapshot.h:81
double SSAGES::Snapshot::TotalMass ( const Label indices) const
inline

Compute the total mass of a group of particles based on index.

Parameters
indicesIDs of particles of interest.
Returns
double Total mass of particles.
Note
This function reduces the mass across the communicator associated with the snapshot.

Definition at line 451 of file Snapshot.h.

References comm_, and masses_.

Referenced by CenterOfMass(), SSAGES::ParticleCoordinateCV::Evaluate(), SSAGES::GyrationTensorCV::Evaluate(), SSAGES::ParticlePositionCV::Evaluate(), SSAGES::ParticleSeparationCV::Evaluate(), and SSAGES::RouseModeCV::setMasses().

452  {
453  auto mtot = 0.;
454  for(auto& i : indices)
455  mtot += masses_[i];
456  MPI_Allreduce(MPI_IN_PLACE, &mtot, 1, MPI_DOUBLE, MPI_SUM, comm_);
457  return mtot;
458  }
mxx::comm comm_
Local snapshot (walker) communicator.
Definition: Snapshot.h:47
std::vector< double > masses_
Masses.
Definition: Snapshot.h:67

Here is the caller graph for this function:

Vector3 SSAGES::Snapshot::UnwrapVector ( const Vector3 v,
const Integer3 image 
) const
inline

Unwrap a vector's real coordinates according to its image replica count.

Parameters
vVector of interest
imageInteger vector representing mirror images in the three dimensions.
Returns
Vector3 Unwrapped vector in real coordinates. This function takes a set of (wrapped) coordinates and returns the unwrapped coordinates.
Note
This function does not require the initial coordinates to be within the simulation box.

Definition at line 409 of file Snapshot.h.

References H_.

410  {
411  return H_*image.cast<double>()+v;
412  }
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:54

Member Data Documentation

Vector3 SSAGES::Snapshot::origin_
private

< Virial tensor.

Box origin.

Definition at line 59 of file Snapshot.h.

Referenced by GetOrigin(), ScaleVector(), and SetOrigin().


The documentation for this class was generated from the following file: