STORMM Source Documentation
Loading...
Searching...
No Matches
stormm::synthesis::PhaseSpaceSynthesis Class Reference

A fixed-precision representation of coordinates, velocities, and forces to manage a set of simulations. The time steps are stored in units of femtoseconds. Coordinates are stored as long long integers to utilize the 32-bit int pipeline for difference computations, bypassing fp64 computations wherever possible to get the best throughput on many visualization, AI-oriented, and gaming cards. More...

#include <phasespace_synthesis.h>

Public Member Functions

HybridFormat getFormat () const
 Get the memory format of the object.
 
int getSystemCount () const
 Get the number of systems in the object.
 
int getUniqueTopologyCount () const
 Get the number of unique topologies in the object.
 
UnitCellType getUnitCellType () const
 Get the unit cell type.
 
CoordinateCycle getCyclePosition () const
 Get the current position in the coordinate cycle.
 
int getGlobalPositionBits () const
 Get the global position scaling bit count.
 
int getLocalPositionBits () const
 Get the local position scaling bit count.
 
int getVelocityBits () const
 Get the velocity scaling bit count.
 
int getForceAccumulationBits () const
 Get the force accumulation bit count.
 
float getForceScalingFactor () const
 Get the force scaling factor used in fixed-precision accumulation.
 
float getInverseForceScalingFactor () const
 Get the inverse of the force scaling factor used in fixed-precision accumulation.
 
const std::vector< AtomGraph * > & getUniqueTopologies () const
 Get a list of unique topology pointers from this coordinate synthesis, const-casted for accessibility to other functions.
 
int getAtomOffset (int system_index) const
 Get the starting index of atoms for one of the systems, using its index.
 
int getAtomCount (int system_index) const
 Get the number of atoms in one of the systems, using its index.
 
int getPaddedAtomCount () const
 Get the total number of atoms across all systems in the synthesis, including padding between systems and after the end of the final system.
 
int getUniqueTopologyIndex (int system_index) const
 Get the index of the unique topology used by a particular system in the synthesis.
 
std::vector< int > getUniqueTopologyIndices () const
 Get the unique topology index of all systems in the synthesis.
 
std::vector< int > getUniqueTopologyExampleIndices () const
 Get a list of system indices from within this coordinate synthesis, providing examples of each unique topology as presented in the order the topology pointers appear in the output of getUniqueTopologies() above.
 
int getTopologyInstanceCount (int topology_index) const
 Get the number of systems that make use of a particular topology.
 
std::vector< int > getSystemIndicesByTopology (int topology_index) const
 Get a list of system indices which all reference the same unique topology (the unique topology index is determined by the order in which each unique topology appears in the overall list of systems in the object).
 
const PhaseSpaceSynthesisgetSelfPointer () const
 Get a const pointer to the object itself.
 
PhaseSpace exportSystem (int index, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Export a system's coordinates, velocities, and forces to a PhaseSpace object.
 
void primeConjugateGradientCalculation ()
 Set the alternate positions to the current forces, and initialize velocities to zero, as part of the first step of conjugate gradient optimization. This primes the system so that the alternate coordinates and velocities arrays can hold the prior forces and temporary conjugate gradient memory.
 
void printTrajectory (const std::vector< int > &system_indices, const std::string &file_name, double current_time, CoordinateFileKind output_kind, PrintSituation expectation) const
 Print a list of structures to a trajectory file. Download will be performed when calling this function, over the subset of relevant frames and data.
 
 PhaseSpaceSynthesis (const std::vector< PhaseSpace > &ps_list, const std::vector< AtomGraph * > &ag_list, int globalpos_scale_bits_in=default_globalpos_scale_bits, int localpos_scale_bits_in=default_localpos_scale_bits, int velocity_scale_bits_in=default_velocity_scale_bits, int force_scale_bits_in=default_force_scale_bits, HybridFormat format_in=default_hpc_format, const GpuDetails &gpu=null_gpu)
 The constructor works from a series PhaseSpace object, importing The LabFrame is not a POINTER-kind object of any sort.
 
 PhaseSpaceSynthesis (const std::vector< PhaseSpace > &ps_list, const std::vector< AtomGraph * > &ag_list, const std::vector< int > &index_key, int globalpos_scale_bits_in=default_globalpos_scale_bits, int localpos_scale_bits_in=default_localpos_scale_bits, int velocity_scale_bits_in=default_velocity_scale_bits, int force_scale_bits_in=default_force_scale_bits, HybridFormat format_in=default_hpc_format, const GpuDetails &gpu=null_gpu)
 
 PhaseSpaceSynthesis (const std::vector< PhaseSpace > &ps_list, const std::vector< int > &ps_index_key, const std::vector< AtomGraph * > &ag_list, const std::vector< int > &ag_index_key, int globalpos_scale_bits_in=default_globalpos_scale_bits, int localpos_scale_bits_in=default_localpos_scale_bits, int velocity_scale_bits_in=default_velocity_scale_bits, int force_scale_bits_in=default_force_scale_bits, HybridFormat format_in=default_hpc_format, const GpuDetails &gpu=null_gpu)
 
 PhaseSpaceSynthesis (const PhaseSpaceSynthesis &original)
 Copy and move constructors work much like their counterparts in the smaller PhaseSpace object.
 
 PhaseSpaceSynthesis (const PhaseSpaceSynthesis &original, HybridFormat format_in)
 
 PhaseSpaceSynthesis (PhaseSpaceSynthesis &&original)
 
PhaseSpaceSynthesisoperator= (const PhaseSpaceSynthesis &other)
 Copy and move assignment operators must also be explicitly written out.
 
PhaseSpaceSynthesisoperator= (PhaseSpaceSynthesis &&other)
 
const AtomGraphgetSystemTopologyPointer (int system_index) const
 Get the topology pointer for a particular system within the synthesis.
 
const std::vector< AtomGraph * > & getSystemTopologyPointer () const
 
const Hybrid< llint > * getCoordinateHandle (CartesianDimension dim, TrajectoryKind kind, CoordinateCycle orientation) const
 Get a pointer to one of the Hybrid objects holding system coordinates. If the synthesis uses extended precision data, the result of this function must be supplemented with the equivalent call to getCoordinateOverflowHandle().
 
const Hybrid< llint > * getCoordinateHandle (CartesianDimension dim, TrajectoryKind kind=TrajectoryKind::POSITIONS) const
 
Hybrid< llint > * getCoordinateHandle (CartesianDimension dim, TrajectoryKind kind, CoordinateCycle orientation)
 
Hybrid< llint > * getCoordinateHandle (CartesianDimension dim, TrajectoryKind kind=TrajectoryKind::POSITIONS)
 
const Hybrid< int > * getCoordinateOverflowHandle (CartesianDimension dim, TrajectoryKind kind, CoordinateCycle orientation) const
 Get a pointer to one of the Hybrid objects holding system coordinates. Overloading and descriptions of input parameters follow from getCoordinateHandle(), above.
 
const Hybrid< int > * getCoordinateOverflowHandle (CartesianDimension dim, TrajectoryKind kind=TrajectoryKind::POSITIONS) const
 
Hybrid< int > * getCoordinateOverflowHandle (CartesianDimension dim, TrajectoryKind kind, CoordinateCycle orientation)
 
Hybrid< int > * getCoordinateOverflowHandle (CartesianDimension dim, TrajectoryKind kind=TrajectoryKind::POSITIONS)
 
const Hybrid< double > * getBoxTransformsHandle (CoordinateCycle orientation) const
 Get a pointer to the box space transforms. Overloading and description of the optional input parameter follows from getCoordinateHandle(), above.
 
const Hybrid< double > * getBoxTransformsHandle () const
 
Hybrid< double > * getBoxTransformsHandle (CoordinateCycle orientation)
 
Hybrid< double > * getBoxTransformsHandle ()
 
const Hybrid< double > * getInverseTransformsHandle (CoordinateCycle orientation) const
 Get a pointer to the inverse transform that takes coordinates back into real space. Overloading and description of the optional input parameter follows from getCoordinateHandle(), above.
 
const Hybrid< double > * getInverseTransformsHandle () const
 
Hybrid< double > * getInverseTransformsHandle (CoordinateCycle orientation)
 
Hybrid< double > * getInverseTransformsHandle ()
 
const Hybrid< double > * getBoxDimensionsHandle (CoordinateCycle orientation) const
 Get a pointer to the transform that takes coordinates back into real space. Overloading and description of the optional input parameter follows from getCoordinateHandle(), above.
 
const Hybrid< double > * getBoxDimensionsHandle () const
 
Hybrid< double > * getBoxDimensionsHandle (CoordinateCycle orientation)
 
Hybrid< double > * getBoxDimensionsHandle ()
 
const Hybrid< llint > * getBoxVectorsHandle (CoordinateCycle orientation) const
 Get a pointer to the box space transforms. Overloading and description of the optional input parameter follows from getCoordinateHandle(), above.
 
const Hybrid< llint > * getBoxVectorsHandle () const
 
Hybrid< llint > * getBoxVectorsHandle (CoordinateCycle orientation)
 
Hybrid< llint > * getBoxVectorsHandle ()
 
const Hybrid< int > * getBoxVectorsOverflowHandle (CoordinateCycle orientation) const
 Get a pointer to the box space transforms. Overloading and description of the optional input parameter follows from getCoordinateHandle(), above.
 
const Hybrid< int > * getBoxVectorsOverflowHandle () const
 
Hybrid< int > * getBoxVectorsOverflowHandle (CoordinateCycle orientation)
 
Hybrid< int > * getBoxVectorsOverflowHandle ()
 
const PsSynthesisReader data (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get the reader or writer, as appropriate based on the const-ness of this object.
 
const PsSynthesisReader data (CoordinateCycle orientation, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
PsSynthesisWriter data (HybridTargetLevel tier=HybridTargetLevel::HOST)
 
PsSynthesisWriter data (CoordinateCycle orientation, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
const PsSynthesisBorders borders (CoordinateCycle orientation, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get the read-only summary of the system sizes. The data is never needed in as a collection of device viewable pointers to host-side data, as this information is constant as of the creation of the object and therefore consistent on both the CPU host and GPU device.
 
const PsSynthesisBorders borders (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
void extractSystem (PhaseSpace *ps, int index, HybridTargetLevel origin=HybridTargetLevel::HOST, HybridTargetLevel destination=HybridTargetLevel::HOST, const GpuDetails &gpu=null_gpu) const
 Extract the phase space (plus forces) of a specific system within the synthesis.
 
CoordinateFrame exportCoordinates (int index, CoordinateCycle orientation, HybridFormat format_out, TrajectoryKind trajkind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Export a system's coordinates, velocities, or forces to a compact CoordinateFrame object.
 
CoordinateFrame exportCoordinates (int index, CoordinateCycle orientation, TrajectoryKind trajkind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
CoordinateFrame exportCoordinates (int index, HybridFormat format_out, TrajectoryKind trajkind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
CoordinateFrame exportCoordinates (int index, TrajectoryKind trajkind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
std::vector< double > getInterlacedCoordinates (int index, CoordinateCycle orientation, TrajectoryKind trajkind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get the interlaced X, Y, and Z coordinates of one system in particular. The result is comparable to the getInterlacedCoordinates() functions held by CoordinateFrame and PhaseSpace objects. The result will be provided in internal units of Angstroms (for positions), Angstroms per femtosecond (for velocities), or kcal/mol-A (for forces).
 
std::vector< double > getInterlacedCoordinates (int index, TrajectoryKind trajkind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
void updateCyclePosition ()
 Modify the object's cycle position, as in the course of propagating dynamics.
 
void updateCyclePosition (CoordinateCycle time_point)
 
void initializeForces (CoordinateCycle orientation, int index=-1)
 Initialize the forces within a PhaseSpaceSynthesis object. This is the analog of the eponymous function in the PhaseSpace object.
 
void initializeForces (int index=-1)
 
void importSystem (const PhaseSpaceReader &psr, int system_index, CoordinateCycle orientation, HybridTargetLevel tier=HybridTargetLevel::HOST)
 Import a system from one of the other coordinate objects, or from a series of C-style arrays with trusted lengths. The imported system's size must correspond to that expected by the atom count of the system it will replace.
 
void importSystem (const PhaseSpaceReader &psr, int system_index, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const PhaseSpaceWriter &psw, int system_index, CoordinateCycle orientation, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const PhaseSpaceWriter &psw, int system_index, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const PhaseSpace &ps, int system_index, CoordinateCycle orientation, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const PhaseSpace &ps, int system_index, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const CoordinateFrameReader &cfr, int system_index, CoordinateCycle orientation, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const CoordinateFrameReader &cfr, int system_index, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const CoordinateFrameWriter &cfw, int system_index, CoordinateCycle orientation, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const CoordinateFrameWriter &cfw, int system_index, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const CoordinateFrame &cf, int system_index, CoordinateCycle orientation, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
void importSystem (const CoordinateFrame &cf, int system_index, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const T *x_import, const T *y_import, const T *z_import, const double *box_xform_in, const double *inverse_xform_in, const double *box_dimensions_in, int system_index, CoordinateCycle orientation, double inverse_scaling_factor=1.0, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const CoordinateSeriesReader< T > &csr, int frame_index, int system_index, CoordinateCycle orientation, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const CoordinateSeriesReader< T > &csr, int frame_index, int system_index, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const CoordinateSeriesWriter< T > &csw, int frame_index, int system_index, CoordinateCycle orientation, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const CoordinateSeriesWriter< T > &csw, int frame_index, int system_index, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const CoordinateSeries< T > &cs, int frame_index, int system_index, CoordinateCycle orientation, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 
template<typename T>
void importSystem (const CoordinateSeries< T > &cs, int frame_index, int system_index, TrajectoryKind kind=TrajectoryKind::POSITIONS, HybridTargetLevel tier=HybridTargetLevel::HOST)
 

Detailed Description

A fixed-precision representation of coordinates, velocities, and forces to manage a set of simulations. The time steps are stored in units of femtoseconds. Coordinates are stored as long long integers to utilize the 32-bit int pipeline for difference computations, bypassing fp64 computations wherever possible to get the best throughput on many visualization, AI-oriented, and gaming cards.

Constructor & Destructor Documentation

◆ PhaseSpaceSynthesis() [1/2]

stormm::synthesis::PhaseSpaceSynthesis::PhaseSpaceSynthesis ( const std::vector< PhaseSpace > & ps_list,
const std::vector< AtomGraph * > & ag_list,
int globalpos_scale_bits_in = default_globalpos_scale_bits,
int localpos_scale_bits_in = default_localpos_scale_bits,
int velocity_scale_bits_in = default_velocity_scale_bits,
int force_scale_bits_in = default_force_scale_bits,
HybridFormat format_in = default_hpc_format,
const GpuDetails & gpu = null_gpu )

The constructor works from a series PhaseSpace object, importing The LabFrame is not a POINTER-kind object of any sort.

Overloaded:

  • Take arrays of PhaseSpace objects and AtomGraph pointers
  • Take a SystemCache object and unpack it
  • Skip the specification of thermostats, barostats, and the time step (just use defaults), but explicitly specify at least the global position scaling bit count and possibly other bit counts
Parameters
ps_listArray of input coordinates, velocities, and forces objects
ag_listArray of pointers to input topologies
index_keyIndices of the given topology and coordinate objects to assemble into a larger list of systems to be held within the resulting PhaseSpaceSynthesis object.
ps_index_keyIndices of the given coordinate objects to assemble into a larger list of systems to be held within the resulting PhaseSpaceSynthesis object.
ag_index_keyIndices of the given topology objects to assemble into a larger list of systems to be held within the resulting PhaseSpaceSynthesis object.

◆ PhaseSpaceSynthesis() [2/2]

stormm::synthesis::PhaseSpaceSynthesis::PhaseSpaceSynthesis ( const PhaseSpaceSynthesis & original)

Copy and move constructors work much like their counterparts in the smaller PhaseSpace object.

Parameters
originalThe original PhaseSpaceSynthesis object
format_inMake the copy with an arbitrary memory layout, copying data based on the priority rules set forth in the deepCopy() function for Hybrid objects.

Member Function Documentation

◆ borders()

const PsSynthesisBorders stormm::synthesis::PhaseSpaceSynthesis::borders ( CoordinateCycle orientation,
HybridTargetLevel tier = HybridTargetLevel::HOST ) const

Get the read-only summary of the system sizes. The data is never needed in as a collection of device viewable pointers to host-side data, as this information is constant as of the creation of the object and therefore consistent on both the CPU host and GPU device.

Parameters
orientationIndicate whether to retrieve box information at the WHITE or BLACK stage of the coordinate cycle. If unspecified, the object's internal coordinate cycle stage will guide the selection.
tierThe level (host or device) at which to get the set of pointers

◆ data()

const PsSynthesisReader stormm::synthesis::PhaseSpaceSynthesis::data ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get the reader or writer, as appropriate based on the const-ness of this object.

Overloaded:

  • Get a read-only abstract for a const object
  • Get a writeable abstract for a mutable object
  • Choose the stage in the time cycle that shall be deemed the relevant coordinates
Parameters
tierThe level (host or device) at which to get the set of pointers
orientationStage in the object's time cycle to take as the current coordinates

◆ exportCoordinates()

CoordinateFrame stormm::synthesis::PhaseSpaceSynthesis::exportCoordinates ( int index,
CoordinateCycle orientation,
HybridFormat format_out,
TrajectoryKind trajkind = TrajectoryKind::POSITIONS,
HybridTargetLevel tier = HybridTargetLevel::HOST ) const

Export a system's coordinates, velocities, or forces to a compact CoordinateFrame object.

Overloaded:

  • Provide the stage of the coordinate cycle from which to obtain coordinates, or make the current point in the time cycle the basis of the new frame, or specify the point in the time cycle to draw from
  • Build the new frame based on the current format of the PhaseSpaceSynthesis, or specify a format
Parameters
indexIndex of the system of interest within the synthesis
trajkindType of coordinates to copy
orientationStage of the coordinate cycle from which to obtain coordinates
format_outThe format of the new frame
tierThe level (host or device) at which to get the data

◆ exportSystem()

PhaseSpace stormm::synthesis::PhaseSpaceSynthesis::exportSystem ( int index,
HybridTargetLevel tier = HybridTargetLevel::HOST ) const

Export a system's coordinates, velocities, and forces to a PhaseSpace object.

Parameters
indexIndex of the system of interest within the synthesis
tierThe level (host or device) at which to get the data

◆ extractSystem()

void stormm::synthesis::PhaseSpaceSynthesis::extractSystem ( PhaseSpace * ps,
int index,
HybridTargetLevel origin = HybridTargetLevel::HOST,
HybridTargetLevel destination = HybridTargetLevel::HOST,
const GpuDetails & gpu = null_gpu ) const

Extract the phase space (plus forces) of a specific system within the synthesis.

Parameters
psPointer to an allocated PhaseSpace object (i.e. the original) ready to accept data from the synthesis (which may have evolved since it was first loaded)
trajkindType of coordinates to copy
indexIndex of the system of interest within the synthesis
originThe level (host or device) at which to get the data
destinationThe level (host or device) at which to get the data
gpuDetails of the GPU in use

◆ getAtomCount()

int stormm::synthesis::PhaseSpaceSynthesis::getAtomCount ( int system_index) const

Get the number of atoms in one of the systems, using its index.

Parameters
system_indexIndex of the system of interest

◆ getAtomOffset()

int stormm::synthesis::PhaseSpaceSynthesis::getAtomOffset ( int system_index) const

Get the starting index of atoms for one of the systems, using its index.

Parameters
system_indexIndex of the system of interest

◆ getCoordinateHandle()

const Hybrid< llint > * stormm::synthesis::PhaseSpaceSynthesis::getCoordinateHandle ( CartesianDimension dim,
TrajectoryKind kind,
CoordinateCycle orientation ) const

Get a pointer to one of the Hybrid objects holding system coordinates. If the synthesis uses extended precision data, the result of this function must be supplemented with the equivalent call to getCoordinateOverflowHandle().

Overloaded:

  • Get a const pointer to a const coordinate synthesis
  • Get a non-const pointer to a mutable coordinate synthesis
  • Indicate the specific stage in the time cycle, or accept the object's current stage
Parameters
dimIndicate Cartesian X, Y, or Z values
kindIndicate positions, velocities, or forces
orientationThe point in the time cycle (WHITE or BLACK) at which to set the pointer

◆ getInterlacedCoordinates()

std::vector< double > stormm::synthesis::PhaseSpaceSynthesis::getInterlacedCoordinates ( int index,
CoordinateCycle orientation,
TrajectoryKind trajkind = TrajectoryKind::POSITIONS,
HybridTargetLevel tier = HybridTargetLevel::HOST ) const

Get the interlaced X, Y, and Z coordinates of one system in particular. The result is comparable to the getInterlacedCoordinates() functions held by CoordinateFrame and PhaseSpace objects. The result will be provided in internal units of Angstroms (for positions), Angstroms per femtosecond (for velocities), or kcal/mol-A (for forces).

Overloaded:

  • Provide the stage of the coordinate cycle from which to obtain coordinates
  • Assume that the WHITE stage of the coordinate cycle is desired
Parameters
indexIndex of the system of interest within the synthesis
trajkindType of coordinates to copy
orientationStage of the coordinate cycle from which to obtain coordinates
tierThe level (host or device) at which to get the data

◆ getSystemIndicesByTopology()

std::vector< int > stormm::synthesis::PhaseSpaceSynthesis::getSystemIndicesByTopology ( int topology_index) const

Get a list of system indices which all reference the same unique topology (the unique topology index is determined by the order in which each unique topology appears in the overall list of systems in the object).

Parameters
topology_indexIndex of the unique topology of interest

◆ getSystemTopologyPointer()

const AtomGraph * stormm::synthesis::PhaseSpaceSynthesis::getSystemTopologyPointer ( int system_index) const

Get the topology pointer for a particular system within the synthesis.

Overloaded:

  • Get the topology pointer for a specific system
  • Get topology pointers for all systems
Parameters
system_indexIndex of the system of interest within the synthesis

◆ getTopologyInstanceCount()

int stormm::synthesis::PhaseSpaceSynthesis::getTopologyInstanceCount ( int topology_index) const

Get the number of systems that make use of a particular topology.

Parameters
topology_indexIndex of the unique topology of interest

◆ getUniqueTopologyIndex()

int stormm::synthesis::PhaseSpaceSynthesis::getUniqueTopologyIndex ( int system_index) const

Get the index of the unique topology used by a particular system in the synthesis.

Parameters
system_indexIndex of the system of interest

◆ importSystem()

void stormm::synthesis::PhaseSpaceSynthesis::importSystem ( const PhaseSpaceReader & psr,
int system_index,
CoordinateCycle orientation,
HybridTargetLevel tier = HybridTargetLevel::HOST )

Import a system from one of the other coordinate objects, or from a series of C-style arrays with trusted lengths. The imported system's size must correspond to that expected by the atom count of the system it will replace.

Overloaded:

  • Provide a PhaseSpace object (all coordinates, velocities, and forces of the input object will be transferred from the WHITE stage of the time cycle in the PhaseSpace object, into the specified stage of the time cycle in this synthesis, and other stages of the time cycle will be transferred accordingly).
  • Provide three arrays of Cartesian X, Y, and Z coordinates, a scaling factor if the data type is fixed-precision integral, plus indications of whether the data is for positions, velocities, or forces, and at what stage of the time cycle the data is to enter the PhaseSpace object.
  • Provide a CoordinateFrame or CoordinateSeries object with a frame number, plus indications of whether the object truly contains positions, velocities, or forces, and what stage of the time cycle the data is to enter the PhaseSpaceSynthesis.
  • Provide an abstract of any of the major coordinate objects, plus the other information to target the import to the correct system, HPC tier, and place in the time cycle.
Parameters
psComplete phase space and time cycle data intended to replace one of the systems in the synthesis
system_indexIndex of the system within this synthesis that the imported coordinates shall replace
orientationStage of the time cycle at which the WHITE stage of an input PhaseSpace object is to enter the synthesis, or at which the data in raw arrays, a CoordinateFrame, or a CoordinateSeries object is to enter the synthesis
tierThe level (host or device) at which to perform the transfer
x_importInput Cartesian X coordinates (these could be positions, velocities, or forces)
y_importInput Cartesian Y coordinates
z_importInput Cartesian Z coordinates
box_xform_inTransformation matrix to take coordinates into fractional space (for positions only–provide nullptr for velocities or forces)
inverse_xform_inTransformation matrix to take coordinates back to real space. The units of elements in this matrix are Angstroms.
box_dimensions_inDimensions of the box (redundant with the information stored in either of the transformation matrices, but convenient and perhaps best able to preserve bitwise information to pass it directly). The units of this array are Angstroms.
kindSpecifies whether the Cartesian X, Y, and Z data are positions, velocities, or forces
scaling_factorScaling factor to take the input X, Y, and Z data into internal units of Angstroms, Angstroms per femtosecond, or kcal/mol-A^2
cfInput coordinate frame object containing X, Y, and Z data as double-precision objects
csInput coordinate series object with X, Y, and Z data for many frames, one of which will be copied over. This object contains its own scaling factor.
frame_indexIndex of a CoordinateSeries object to be transferred

◆ initializeForces()

void stormm::synthesis::PhaseSpaceSynthesis::initializeForces ( CoordinateCycle orientation,
int index = -1 )

Initialize the forces within a PhaseSpaceSynthesis object. This is the analog of the eponymous function in the PhaseSpace object.

Overloaded:

  • Initialize one or more systems' forces on the host if STORMM is compiled for CPU only
  • Initialize one or more systems' forces on either the host or the device if STORMM is compiled for HPC
  • Initialize forces for a specific point in the time cycle, or the current point
Parameters
orientationPoint in the time cycle at which to clear / initialize the forces
gpuDetails of the GPU in use
tierThe level (host or device) at which to initialize forces
indexIndex of the system of interest within the synthesis–if negative, all systems will have their forces initialized

◆ operator=()

PhaseSpaceSynthesis & stormm::synthesis::PhaseSpaceSynthesis::operator= ( const PhaseSpaceSynthesis & other)

Copy and move assignment operators must also be explicitly written out.

Parameters
otherThe other PhaseSpaceSynthesis object

◆ primeConjugateGradientCalculation()

void stormm::synthesis::PhaseSpaceSynthesis::primeConjugateGradientCalculation ( )

Set the alternate positions to the current forces, and initialize velocities to zero, as part of the first step of conjugate gradient optimization. This primes the system so that the alternate coordinates and velocities arrays can hold the prior forces and temporary conjugate gradient memory.

Parameters
gpuDetails of the GPU in use
tierThe level (host or device) at which to initialize vectors

◆ printTrajectory()

void stormm::synthesis::PhaseSpaceSynthesis::printTrajectory ( const std::vector< int > & system_indices,
const std::string & file_name,
double current_time,
CoordinateFileKind output_kind,
PrintSituation expectation ) const

Print a list of structures to a trajectory file. Download will be performed when calling this function, over the subset of relevant frames and data.

Parameters
system_indicesList of system coordinates / velocities / forces to print
file_nameName of the file to write, or base name of a set of files to write
current_timeCurrent time progress of the group of simulations
output_kindThe type of trajectory file to write
expectationThe state that the output trajectory file is expected to be found in

◆ updateCyclePosition()

void stormm::synthesis::PhaseSpaceSynthesis::updateCyclePosition ( )

Modify the object's cycle position, as in the course of propagating dynamics.

Overloaded:

  • Increment the cycle position one step forward (this toggles between the WHITE and BLACK states, and is the same as decrementing the cycle position unless the time cycles takes on a third possible stage)
  • Set the time cycle to a specific point
Parameters
time_pointThe point in the time cycle which the object is to take as "current"

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