SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
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, unsigned int wid)
 Constructor. More...
 
size_t GetIteration () const
 Get the current iteration. More...
 
size_t 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 (size_t 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< 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...
 
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, bool mass_weight=true) const
 Compute center of mass of a group of atoms based on index. More...
 
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::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. More...
 
std::vector< double > SerializeVelocities ()
 Return the serialized velocities across all local cores. More...
 
std::vector< int > SerializeIDs ()
 Return the serialized IDs across all local cores. More...
 
 ~Snapshot ()
 Destructor.
 

Private Attributes

mxx::comm comm_
 Local snapshot (walker) communicator.
 
unsigned int wid_
 Walker ID.
 
unsigned int 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_
 Virial tensor.
 
Vector3 origin_
 Box origin.
 
Bool3 isperiodic_
 Periodicity of box.
 
std::vector< Vector3positions_
 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.
 
size_t iteration_
 Iteration of Simulation.
 
size_t 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 47 of file Snapshot.h.

Constructor & Destructor Documentation

◆ Snapshot()

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

Constructor.

Parameters
commMPI communicator
widWalker ID

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

Definition at line 94 of file Snapshot.h.

94  :
95  comm_(comm), wid_(wid), H_(), Hinv_(), virial_(Matrix3::Zero()),
96  origin_({0,0,0}), isperiodic_({true, true, true}), positions_(0),
97  velocities_(0), forces_(0), masses_(0), atomids_(0), types_(0),
98  iteration_(0), temperature_(0), energy_(0), kb_(0)
99  {}
double energy_
Average per-particle energy.
Definition: Snapshot.h:78
std::vector< double > masses_
Masses.
Definition: Snapshot.h:70
Vector3 origin_
Box origin.
Definition: Snapshot.h:63
Label atomids_
List of Atom IDs.
Definition: Snapshot.h:72
double temperature_
Current temperature.
Definition: Snapshot.h:77
Bool3 isperiodic_
Periodicity of box.
Definition: Snapshot.h:65
mxx::comm comm_
Local snapshot (walker) communicator.
Definition: Snapshot.h:51
std::vector< Vector3 > positions_
Positions.
Definition: Snapshot.h:67
std::vector< Vector3 > velocities_
Velocities.
Definition: Snapshot.h:68
Label types_
List of Atom types.
Definition: Snapshot.h:73
Matrix3 H_
Parrinello-Rahman box H-matrix.
Definition: Snapshot.h:58
Matrix3 virial_
Virial tensor.
Definition: Snapshot.h:61
unsigned int wid_
Walker ID.
Definition: Snapshot.h:53
Matrix3 Hinv_
Parinello-Rahman box inverse.
Definition: Snapshot.h:59
double kb_
Kb from the MD driver.
Definition: Snapshot.h:79
std::vector< Vector3 > forces_
Forces.
Definition: Snapshot.h:69
size_t iteration_
Iteration of Simulation.
Definition: Snapshot.h:75

Member Function Documentation

◆ ApplyMinimumImage() [1/2]

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 405 of file Snapshot.h.

406  {
407  Vector3 scaled = Hinv_*v;
408 
409  for(int i = 0; i < 3; ++i)
410  scaled[i] -= isperiodic_[i]*roundf(scaled[i]);
411 
412  return H_*scaled;
413  }
double roundf(double x)
Quick helper function to round a double.
Definition: Snapshot.h:35
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33

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

Here is the call graph for this function:

◆ ApplyMinimumImage() [2/2]

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

Apply minimum image to a vector.

Parameters
vVector of interest

Definition at line 390 of file Snapshot.h.

391  {
392  Vector3 scaled = Hinv_*(*v);
393 
394  for(int i = 0; i < 3; ++i)
395  scaled[i] -= isperiodic_[i]*roundf(scaled[i]);
396 
397  *v = H_*scaled;
398  }

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

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CenterOfMass()

Vector3 SSAGES::Snapshot::CenterOfMass ( const Label indices,
bool  mass_weight = true 
) const
inline

Compute center of mass of a group of atoms based on index.

Parameters
indicesIDs of particles of interest.
mass_weightMass-weighting for COM calculation.
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.
Default behavior is to use a mass-weighted average.

Definition at line 442 of file Snapshot.h.

443  {
444  // Store coordinates and masses in vectors to gather.
445  std::vector<double> pos, mass, gpos, gmass;
446  std::vector<int> pcounts(comm_.size(), 0), mcounts(comm_.size(), 0);
447  std::vector<int> pdispls(comm_.size()+1, 0), mdispls(comm_.size()+1, 0);
448 
449  // Calculate total mass
450  double mtot = mass_weight ? TotalMass(indices) : indices.size();
451 
452  pcounts[comm_.rank()] = 3*indices.size();
453  mcounts[comm_.rank()] = indices.size();
454 
455  // Reduce counts.
456  MPI_Allreduce(MPI_IN_PLACE, pcounts.data(), pcounts.size(), MPI_INT, MPI_SUM, comm_);
457  MPI_Allreduce(MPI_IN_PLACE, mcounts.data(), mcounts.size(), MPI_INT, MPI_SUM, comm_);
458 
459  // Compute displacements.
460  std::partial_sum(pcounts.begin(), pcounts.end(), pdispls.begin() + 1);
461  std::partial_sum(mcounts.begin(), mcounts.end(), mdispls.begin() + 1);
462 
463  // Fill up mass and position vectors.
464  for(auto& idx : indices)
465  {
466  auto& p = positions_[idx];
467  pos.push_back(p[0]);
468  pos.push_back(p[1]);
469  pos.push_back(p[2]);
470  mass.push_back(masses_[idx]);
471  }
472 
473  // Re-size receiving vectors.
474  gpos.resize(pdispls.back(), 0);
475  gmass.resize(mdispls.back(), 1.0); // fictitious mass of 1
476 
477  // All-gather data.
478  MPI_Allgatherv(pos.data(), pos.size(), MPI_DOUBLE, gpos.data(), pcounts.data(), pdispls.data(), MPI_DOUBLE, comm_);
479  if(mass_weight) { // overwrite gmass with true mass values
480  MPI_Allgatherv(mass.data(), mass.size(), MPI_DOUBLE, gmass.data(), mcounts.data(), mdispls.data(), MPI_DOUBLE, comm_);
481  }
482 
483  // Loop through atoms and compute mass weighted sum.
484  // We march linearly through list and find nearest image
485  // to each successive particle to properly unwrap object.
486  Vector3 ppos = {gpos[0], gpos[1], gpos[2]}; // Previous unwrapped position.
487  Vector3 cpos = ppos;
488  Vector3 xcm = gmass[0]*cpos;
489 
490  for(size_t i = 1, j = 3; i < gmass.size(); ++i, j += 3)
491  {
492  cpos = {gpos[j], gpos[j+1], gpos[j+2]};
493  cpos = ApplyMinimumImage(cpos - ppos) + ppos;
494  xcm += gmass[i]*cpos;
495  ppos = cpos;
496  }
497 
498  return xcm/mtot;
499  }
double TotalMass(const Label &indices) const
Compute the total mass of a group of particles based on index.
Definition: Snapshot.h:423
void ApplyMinimumImage(Vector3 *v) const
Apply minimum image to a vector.
Definition: Snapshot.h:390

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

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Changed()

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 596 of file Snapshot.h.

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

References changed_.

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

Here is the caller graph for this function:

◆ GetAtomIDs() [1/2]

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

Access the atom IDs.

Returns
List of atom IDs

Definition at line 508 of file Snapshot.h.

509  {
510  changed_ = true;
511  return atomids_;
512  }

References atomids_, and changed_.

◆ GetAtomIDs() [2/2]

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

Access the atom IDs.

Returns
List of atom IDs

Definition at line 505 of file Snapshot.h.

505 { return atomids_; }

References atomids_.

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

Here is the caller graph for this function:

◆ GetAtomTypes() [1/2]

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

Access the atom types.

Returns
List of atom types

Definition at line 567 of file Snapshot.h.

568  {
569  changed_ = true;
570  return types_;
571  }

References changed_, and types_.

◆ GetAtomTypes() [2/2]

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

Access the atom types.

Returns
List of atom types

Definition at line 564 of file Snapshot.h.

564 { return types_; }

References types_.

◆ GetCharges() [1/2]

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

Access the atom charges.

Returns
List of atom charges

Definition at line 554 of file Snapshot.h.

555  {
556  changed_ = true;
557  return charges_;
558  }
std::vector< double > charges_
Charges.
Definition: Snapshot.h:71

References changed_, and charges_.

◆ GetCharges() [2/2]

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

Access the atom charges.

Returns
List of atom charges

Definition at line 551 of file Snapshot.h.

551 { return charges_; }

References charges_.

◆ GetCommunicator()

const mxx::comm& SSAGES::Snapshot::GetCommunicator ( ) const
inline

◆ GetDielectric()

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

Get dielectric constant.

Returns
dielectric constant.

Definition at line 171 of file Snapshot.h.

171 { return dielectric_; }
double dielectric_
Dielectric.
Definition: Snapshot.h:81

References dielectric_.

◆ GetEnergy()

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

Get per-particle energy.

Returns
Current average per-particle energy

Definition at line 123 of file Snapshot.h.

123 { return energy_; }

References energy_.

◆ GetForces() [1/2]

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

Access the per-particle forces.

Returns
List of per-particle forces

Definition at line 354 of file Snapshot.h.

355  {
356  changed_ = true;
357  return forces_;
358  }

References changed_, and forces_.

◆ GetForces() [2/2]

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

◆ GetHMatrix()

const Matrix3& SSAGES::Snapshot::GetHMatrix ( ) const
inline

Get system H-matrix.

Returns
Parrinello-Rahman H-matrix of simulation box

Definition at line 129 of file Snapshot.h.

129 { return H_; }

References H_.

◆ GetIteration()

size_t SSAGES::Snapshot::GetIteration ( ) const
inline

Get the current iteration.

Returns
Current Iteration

Definition at line 105 of file Snapshot.h.

105 {return iteration_; }

References iteration_.

Referenced by SSAGES::Logger::PostIntegration(), SSAGES::ANN::PostIntegration(), SSAGES::BFS::PostIntegration(), SSAGES::CFF::PostIntegration(), SSAGES::Meta::PostIntegration(), SSAGES::Umbrella::PostIntegration(), and SSAGES::Hook::PostIntegrationHook().

Here is the caller graph for this function:

◆ GetKb()

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

Get system Kb.

Returns
Kb of the current simulation box

Definition at line 165 of file Snapshot.h.

165 { return kb_; }

References kb_.

Referenced by SSAGES::BFS::PostIntegration(), SSAGES::ANN::PreSimulation(), and SSAGES::CFF::PreSimulation().

Here is the caller graph for this function:

◆ GetLocalIndex()

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 519 of file Snapshot.h.

520  {
521 
522  auto s = std::find(atomids_.begin(), atomids_.end(), id);
523  if(s == atomids_.end())
524  return -1;
525  else
526  return s - atomids_.begin();
527  }

References atomids_.

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

Here is the caller graph for this function:

◆ GetLocalIndices()

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 537 of file Snapshot.h.

538  {
539  for(auto& id : ids)
540  {
541  auto idx = GetLocalIndex(id);
542  if(idx != -1)
543  indices->push_back(idx);
544  }
545  }
int GetLocalIndex(int id) const
Gets the local atom index corresponding to an atom ID.
Definition: Snapshot.h:519

References GetLocalIndex().

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

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMasses() [1/2]

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 370 of file Snapshot.h.

371  {
372  changed_ = true;
373  return masses_;
374  }

References changed_, and masses_.

◆ GetMasses() [2/2]

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 367 of file Snapshot.h.

367 { return masses_; }

References masses_.

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

Here is the caller graph for this function:

◆ GetNumAtoms()

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

◆ GetOrigin()

const Vector3& SSAGES::Snapshot::GetOrigin ( ) const
inline

Get origin of the system.

Returns
Vector containing coordinates of box origin.

Definition at line 147 of file Snapshot.h.

147 { return origin_; }

References origin_.

◆ GetPositions() [1/2]

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

Access the particle positions.

Returns
List of particle positions

Definition at line 328 of file Snapshot.h.

329  {
330  changed_ = true;
331  return positions_;
332  }

References changed_, and positions_.

◆ GetPositions() [2/2]

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

◆ Getqqrd2e()

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

Get qqrd2e value.

Returns
Value of qqrd2e.

Definition at line 177 of file Snapshot.h.

177 { return qqrd2e_; }
double qqrd2e_
qqrd2e
Definition: Snapshot.h:82

References qqrd2e_.

◆ GetSnapshotID() [1/2]

std::string& SSAGES::Snapshot::GetSnapshotID ( )
inline

Access the snapshot ID.

Returns
Snapshot ID

Definition at line 580 of file Snapshot.h.

581  {
582  changed_ = true;
583  return ID_;
584  }
std::string ID_
ID string.
Definition: Snapshot.h:56

References changed_, and ID_.

◆ GetSnapshotID() [2/2]

const std::string& SSAGES::Snapshot::GetSnapshotID ( ) const
inline

Access the snapshot ID.

Returns
Snapshot ID

Definition at line 577 of file Snapshot.h.

577 { return ID_; }

References ID_.

◆ GetTargetIterations()

size_t SSAGES::Snapshot::GetTargetIterations ( ) const
inline

Get target iterations.

Returns
Target Iterations

Definition at line 111 of file Snapshot.h.

111 {return targetiter_; }
size_t targetiter_
Iteration target of simulation.
Definition: Snapshot.h:76

References targetiter_.

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

Here is the caller graph for this function:

◆ GetTemperature()

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

Get current temperature.

Returns
System temperature

Definition at line 117 of file Snapshot.h.

117 {return temperature_; }

References temperature_.

◆ GetVelocities() [1/2]

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

Access the particle velocities.

Returns
List of particle velocities

Definition at line 341 of file Snapshot.h.

342  {
343  changed_ = true;
344  return velocities_;
345  }

References changed_, and velocities_.

◆ GetVelocities() [2/2]

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

Access the particle velocities.

Returns
List of particle velocities

Definition at line 338 of file Snapshot.h.

338 { return velocities_; }

References velocities_.

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

Here is the caller graph for this function:

◆ GetVirial() [1/2]

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

Get box virial.

Returns
Virial tensor of the simulation box.

Definition at line 141 of file Snapshot.h.

141 { return virial_; }

References virial_.

◆ GetVirial() [2/2]

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

Get box virial.

Returns
Virial tensor of the simulation box.

Definition at line 135 of file Snapshot.h.

135 { return virial_; }

References virial_.

Referenced by SSAGES::ANN::PostIntegration(), SSAGES::BFS::PostIntegration(), SSAGES::CFF::PostIntegration(), SSAGES::Meta::PostIntegration(), and SSAGES::Umbrella::PostIntegration().

Here is the caller graph for this function:

◆ GetVolume()

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

Get system volume.

Returns
Volume of the current simulation box

Definition at line 159 of file Snapshot.h.

159 { return H_.determinant(); }

References H_.

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

Here is the caller graph for this function:

◆ GetWalkerID()

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

Get walker ID.

Returns
ID of the Walker

Definition at line 195 of file Snapshot.h.

195 { return wid_; }

References wid_.

Referenced by SSAGES::BFS::PreSimulation(), and SSAGES::StringMethod::PreSimulation().

Here is the caller graph for this function:

◆ HasChanged()

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

Query if Snapshot was modified.

Returns
True if Snapshot was modified, else return False

Definition at line 590 of file Snapshot.h.

590 { return changed_; }

References changed_.

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

Here is the caller graph for this function:

◆ IsPeriodic()

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

Get periodicity of three dimensions.

Returns
Three dimensional boolean containing periodicity of each dimension.

Definition at line 153 of file Snapshot.h.

153 {return isperiodic_; }

References isperiodic_.

◆ ScaleVector()

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 381 of file Snapshot.h.

382  {
383  return Hinv_*(v-origin_);
384  }

References Hinv_, and origin_.

◆ SerializeIDs()

std::vector<int> SSAGES::Snapshot::SerializeIDs ( )
inline

Return the serialized IDs across all local cores.

Returns
IDs across all local cores.

Definition at line 673 of file Snapshot.h.

674  {
675  std::vector<int> mcounts(comm_.size(), 0);
676  std::vector<int> mdispls(comm_.size()+1, 0);
677 
678  mcounts[comm_.rank()] = nlocal_;
679 
680  // Reduce counts.
681  MPI_Allreduce(MPI_IN_PLACE, mcounts.data(), mcounts.size(), MPI_INT, MPI_SUM, comm_);
682 
683  // Compute displacements.
684  std::partial_sum(mcounts.begin(), mcounts.end(), mdispls.begin() + 1);
685 
686  // Re-size receiving vectors.
687  std::vector<int> IDs;
688  IDs.resize(mdispls.back(), 0);
689 
690  // All-gather data.
691  MPI_Allgatherv(atomids_.data(), atomids_.size(), MPI_INT, IDs.data(), mcounts.data(), mdispls.data(), MPI_INT, comm_);
692  return IDs;
693  }

References atomids_, comm_, and nlocal_.

Referenced by SSAGES::FiniteTempString::PostIntegration(), and SSAGES::Swarm::PostIntegration().

Here is the caller graph for this function:

◆ SerializePositions()

std::vector<double> SSAGES::Snapshot::SerializePositions ( )
inline

Return the serialized positions across all local cores.

Returns
Positions across all local cores.

Definition at line 602 of file Snapshot.h.

603  {
604 
605  std::vector<int> pcounts(comm_.size(), 0);
606  std::vector<int> pdispls(comm_.size()+1, 0);
607 
608  pcounts[comm_.rank()] = 3*nlocal_;
609 
610  // Reduce counts.
611  MPI_Allreduce(MPI_IN_PLACE, pcounts.data(), pcounts.size(), MPI_INT, MPI_SUM, comm_);
612 
613  // Compute displacements.
614  std::partial_sum(pcounts.begin(), pcounts.end(), pdispls.begin() + 1);
615 
616  // Re-size receiving vectors.
617  std::vector<double> positions;
618  positions.resize(pdispls.back(), 0);
619 
620  std::vector<double> ptemp;
621 
622  for(auto& p : positions_)
623  {
624  ptemp.push_back(p[0]);
625  ptemp.push_back(p[1]);
626  ptemp.push_back(p[2]);
627  }
628 
629  // All-gather data.
630  MPI_Allgatherv(ptemp.data(), ptemp.size(), MPI_DOUBLE, positions.data(), pcounts.data(), pdispls.data(), MPI_DOUBLE, comm_);
631  return positions;
632  }

References comm_, nlocal_, and positions_.

Referenced by SSAGES::FiniteTempString::PostIntegration(), and SSAGES::Swarm::PostIntegration().

Here is the caller graph for this function:

◆ SerializeVelocities()

std::vector<double> SSAGES::Snapshot::SerializeVelocities ( )
inline

Return the serialized velocities across all local cores.

Returns
Velocities across all local cores.

Definition at line 638 of file Snapshot.h.

639  {
640  std::vector<int> vcounts(comm_.size(), 0);
641  std::vector<int> vdispls(comm_.size()+1, 0);
642 
643  vcounts[comm_.rank()] = 3*nlocal_;
644 
645  // Reduce counts.
646  MPI_Allreduce(MPI_IN_PLACE, vcounts.data(), vcounts.size(), MPI_INT, MPI_SUM, comm_);
647 
648  // Compute displacements.
649  std::partial_sum(vcounts.begin(), vcounts.end(), vdispls.begin() + 1);
650 
651  // Re-size receiving vectors.
652  std::vector<double> velocities;
653  velocities.resize(vdispls.back(), 0);
654 
655  std::vector<double> vtemp;
656 
657  for(auto& v : velocities_)
658  {
659  vtemp.push_back(v[0]);
660  vtemp.push_back(v[1]);
661  vtemp.push_back(v[2]);
662  }
663 
664  // All-gather data.
665  MPI_Allgatherv(vtemp.data(), vtemp.size(), MPI_DOUBLE, velocities.data(), vcounts.data(), vdispls.data(), MPI_DOUBLE, comm_);
666  return velocities;
667  }

References comm_, nlocal_, and velocities_.

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

Here is the caller graph for this function:

◆ SetDielectric()

void SSAGES::Snapshot::SetDielectric ( double  dielectric)
inline

Set the dielectric constant.

Parameters
dielectricValue for the dielectric constant.

Definition at line 299 of file Snapshot.h.

300  {
301  dielectric_ = dielectric;
302  changed_ = true;
303  }

References changed_, and dielectric_.

◆ SetEnergy()

void SSAGES::Snapshot::SetEnergy ( double  energy)
inline

Change the energy.

Parameters
energyNew value for the energy

Definition at line 238 of file Snapshot.h.

239  {
240  energy_ = energy;
241  changed_ = true;
242  }

References changed_, and energy_.

◆ SetHMatrix()

void SSAGES::Snapshot::SetHMatrix ( const Matrix3 hmat)
inline

Change the Box H-matrix.

Parameters
hmatNew H-matrix for the system

Definition at line 248 of file Snapshot.h.

249  {
250  H_ = hmat;
251  Hinv_ = hmat.inverse();
252  changed_ = true;
253  }

References changed_, H_, and Hinv_.

◆ SetIteration()

void SSAGES::Snapshot::SetIteration ( size_t  iteration)
inline

Set the iteration.

Parameters
iterationNew value for the iteration

Definition at line 208 of file Snapshot.h.

209  {
210  iteration_ = iteration;
211  changed_ = true;
212  }

References changed_, and iteration_.

◆ SetKb()

void SSAGES::Snapshot::SetKb ( double  kb)
inline

Change the kb.

Parameters
kbNew value for the kb

Definition at line 289 of file Snapshot.h.

290  {
291  kb_ = kb;
292  changed_ = true;
293  }

References changed_, and kb_.

◆ SetNumAtoms()

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 319 of file Snapshot.h.

319 { nlocal_ = natoms; }

References nlocal_.

◆ SetOrigin()

void SSAGES::Snapshot::SetOrigin ( const Vector3 origin)
inline

Change the box origin.

Parameters
originNew origin for the system

Definition at line 269 of file Snapshot.h.

270  {
271  origin_ = origin;
272  changed_ = true;
273  }

References changed_, and origin_.

◆ SetPeriodicity()

void SSAGES::Snapshot::SetPeriodicity ( const Bool3 isperiodic)
inline

Change the periodicity of the system.

Parameters
isperiodicPeriodicity of three dimensions

Definition at line 279 of file Snapshot.h.

280  {
281  isperiodic_ = isperiodic;
282  changed_ = true;
283  }

References changed_, and isperiodic_.

◆ Setqqrd2e()

void SSAGES::Snapshot::Setqqrd2e ( double  qqrd2e)
inline

Set the value for qqrd2e.

Parameters
qqrd2eValue for qqrd2e.

Definition at line 309 of file Snapshot.h.

310  {
311  qqrd2e_ = qqrd2e;
312  changed_ = true;
313  }

References changed_, and qqrd2e_.

◆ SetTargetIterations()

void SSAGES::Snapshot::SetTargetIterations ( int  target)
inline

Set target iterations.

Parameters
targetNew value for target iterations

Definition at line 218 of file Snapshot.h.

219  {
220  targetiter_ = target;
221  changed_ = true;
222  }

References changed_, and targetiter_.

◆ SetTemperature()

void SSAGES::Snapshot::SetTemperature ( double  temperature)
inline

Change the temperature.

Parameters
temperatureNew value for the temperature

Definition at line 228 of file Snapshot.h.

229  {
230  temperature_ = temperature;
231  changed_ = true;
232  }

References changed_, and temperature_.

◆ SetVirial()

void SSAGES::Snapshot::SetVirial ( const Matrix3 virial)
inline

Change the box virial.

Parameters
virialNew virial tensor for the system.

Definition at line 259 of file Snapshot.h.

260  {
261  virial_ = virial;
262  changed_ = true;
263  }

References changed_, and virial_.

◆ TotalMass()

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 423 of file Snapshot.h.

424  {
425  auto mtot = 0.;
426  for(auto& i : indices)
427  mtot += masses_[i];
428  MPI_Allreduce(MPI_IN_PLACE, &mtot, 1, MPI_DOUBLE, MPI_SUM, comm_);
429  return mtot;
430  }

References comm_, and masses_.

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

Here is the caller graph for this function:

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