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

A collection of one or more AtomGraph objects, with similar components arranged in contiguous arrays (often padded by the GPU warp size to prevent one system from flowing into another). Work units covering all systems are laid out in additional, contiguous arrays. While individual topologies (AtomGraphs) have Hybrid POINTER-kind objects for details such as charges and atom type indices and all of the POINTER-kind objects targeted a single ARRAY-kind object of the correct memory type, the synthesis, which may have many topologies and a great deal more overall information, stores most of its data in a series of ARRAY-kind objects, one for each member variable. More...

#include <atomgraph_synthesis.h>

Public Member Functions

int getUniqueTopologyCount () const
 Get the number of unique topologies described by the synthesis.
 
const std::vector< AtomGraph * > & getUniqueTopologies () const
 Get a const reference to the list of all pointers for unique topologies making up this synthesis.
 
std::vector< int > getTopologyIndices () const
 Get a const reference to the list of all topology indices, indicating the topologies describing each systems contained within this synthesis.
 
RestraintApparatusgetSystemRestraintPointer (int system_index) const
 Get a restraint apparatus pointer for a sepcific system contained within the synthesis.
 
int getSystemCount () const
 Get the number of systems that this synthesis describes.
 
int getPaddedAtomCount () const
 Get the entire padded number of atoms covering all systems.
 
int getAtomOffset (int system_index) const
 Get the starting point for atoms of a specific system in the lineup of all topologies.
 
int getVirtualSiteCount () const
 Get the total number of virtual sites across all systems, including replicas.
 
int getBondTermCount () const
 Get the total number of bond terms across all systems, including replicas.
 
int getAngleTermCount () const
 Get the total number of bond angle terms across all systems, including replicas.
 
int getDihedralTermCount () const
 Get the total number of cosine-based dihedral terms across all systems and replicas.
 
int getUreyBradleyTermCount () const
 Get the total number of Urey-Bradley terms across all systems and replicas.
 
int getCharmmImproperTermCount () const
 Get the total number of CHARMM improper terms across all systems and replicas.
 
int getCmapTermCount () const
 Get the total number of CMAP terms across all systems and replicas.
 
int getLJTypeCount () const
 Get the number of unique atom types (a parameter, not an extensive quantity dependent on the number of systems).
 
int getChargeTypeCount () const
 Get the number of unique charge parameters.
 
int getBondParameterCount () const
 Get the number of unique harmonic bond parameter sets.
 
int getAngleParameterCount () const
 Get the number of unique harmonic bond angle parameter sets.
 
int getDihedralParameterCount () const
 Get the number of unique cosine-based dihedral parameter sets.
 
int getUreyBradleyParameterCount () const
 Get the number of unique Urey-Bradley harmonic angle parameter sets.
 
int getCharmmImproperParameterCount () const
 Get the number of unique CHARMM improper parameter sets.
 
int getCmapSurfaceCount () const
 Get the number of unique CMAP surfaces.
 
const Hybrid< int > & getSystemAtomCounts () const
 Get the sizes of all individual systems as a const reference to the Hybrid array member variable.
 
const Hybrid< int > & getSystemAtomOffsets () const
 Get the starting locations of all individual systems as a const reference to the Hybrid array member variable.
 
const Hybrid< int > & getDegreesOfFreedom () const
 Get the numbers of unconstrained degrees of freedom for all individual systems in the.
 
const Hybrid< int > & getConstrainedDegreesOfFreedom () const
 Get the numbers of constrained degrees of freedom for all individual systems in the.
 
UnitCellType getUnitCellType () const
 Get the unit cell type that will be taken for all systems (TRICLINIC subsumes ORTHORHOMBIC in a sort of type promotion).
 
ImplicitSolventModel getImplicitSolventModel () const
 Get the implicit solvent model in use across all systems.
 
double getDielectricConstant () const
 Get the dielectric constant (supporting the implicit solvent model) for all systems.
 
double getSaltConcentration () const
 Get the salt concentration (supporting the implicit solvent model) for all systems.
 
double getCoulombConstant () const
 Get the fundamental Coulomb constant defining the electrostatics of all systems.
 
int getValenceWorkUnitCount () const
 Get the overall number of valence work units needed to account for interactions in all systems.
 
int2 getValenceWorkUnitSize () const
 Get the maximum size of valence work units, whether for organic molecules or solvent (rigid water) residues. This will have been either set by the user or tailored by the automated heuristics to produce the best saturation.
 
ValenceKernelSize getValenceThreadBlockSize () const
 Get the necessary thread block size for evaluating the valence work units.
 
const Hybrid< int2 > & getValenceWorkUnitAbstracts () const
 Get the abstracts (condensed lists of import and instruction set limits) for the valence work units spanning all systems in this synthesis.
 
NbwuKind getNonbondedWorkType () const
 Get the type of non-bonded work required by systems in this synthesis.
 
int getNonbondedWorkUnitCount () const
 Get the number of non-bonded work units serving systems in this synthesis.
 
int getRandomCacheDepth () const
 Get the depth of the random numbers cache that the non-bonded work units of this topology synthesis will try to fill. This will raise a warning and return zero if no cache depth has yet been set.
 
int getReductionWorkUnitCount () const
 Get the number of reduction work units spanning all systems.
 
RdwuPerSystem getRdwuPerSystem () const
 Get a qualitative assessment of the number of reduction work units assigned to any one system.
 
const Hybrid< int > & getReductionWorkUnitAbstracts () const
 Get the reduction work unit abstracts (this is not required for valence or non-bonded work units as the abstracts needed for implementing these work units are produced by member functions of this class (see below). The reduction work, however, is managed by an abstract splicing together elements of this class, PhaseSpaceSynthesis, and a small class allocating temporary storage space.
 
SyValenceKit< double > getDoublePrecisionValenceKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with double-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.
 
SyValenceKit< float > getSinglePrecisionValenceKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with single-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.
 
SyRestraintKit< double, double2, double4getDoublePrecisionRestraintKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with double-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.
 
SyRestraintKit< float, float2, float4getSinglePrecisionRestraintKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with single-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.
 
SyNonbondedKit< double, double2getDoublePrecisionNonbondedKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with double-precision real numbers for computing non-bonded interactions for all systems based on work units stored in this object.
 
SyNonbondedKit< float, float2getSinglePrecisionNonbondedKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with single-precision real numbers for computing non-bonded interactions for all systems based on work units stored in this object.
 
SyAtomUpdateKit< double, double2, double4getDoublePrecisionAtomUpdateKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with double-precision real numbers for updating atom and virtual site positions based on work units stored in this object.
 
SyAtomUpdateKit< float, float2, float4getSinglePrecisionAtomUpdateKit (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get a minimal kit with single-precision real numbers for updating atom and virtual site positions based on work units stored in this object.
 
const AtomGraphSynthesisgetSelfPointer () const
 Get a const pointer to the object itself, useful for retrieving a valid pointer when the object is available as a const reference.
 
void loadNonbondedWorkUnits (const StaticExclusionMaskSynthesis &poly_se, InitializationTask init_request=InitializationTask::NONE, int random_cache_depth=0, const GpuDetails &gpu=null_gpu)
 Construct non-bonded work units for all unique topologies (there are no restraints for non-bonded interactions that might distinguish systems with the same topology, as was a consideration when developing the valence work units). Load the instructions into the topology synthesis for availability on the GPU.
 
 AtomGraphSynthesis (const std::vector< AtomGraph * > &topologies_in, const std::vector< RestraintApparatus * > &restraints_in, const std::vector< int > &topology_indices_in, const std::vector< int > &restraint_indices_in, ExceptionResponse policy_in=ExceptionResponse::WARN, const GpuDetails &gpu=null_gpu, StopWatch *timer_in=nullptr)
 The constructor takes a series of topologies and NMR restraints. The NMR restraints point to specific topologies and thereby apply to any coordinate sets that also point to those topologies.
 
 AtomGraphSynthesis (const std::vector< AtomGraph * > &topologies_in, const std::vector< RestraintApparatus * > &restraints_in, ExceptionResponse policy_in=ExceptionResponse::WARN, const GpuDetails &gpu=null_gpu, StopWatch *timer_in=nullptr)
 
 AtomGraphSynthesis (const std::vector< AtomGraph * > &topologies_in, const std::vector< int > &topology_indices_in, ExceptionResponse policy_in=ExceptionResponse::WARN, const GpuDetails &gpu=null_gpu, StopWatch *timer_in=nullptr)
 
 AtomGraphSynthesis (const std::vector< AtomGraph * > &topologies_in, ExceptionResponse policy_in=ExceptionResponse::WARN, const GpuDetails &gpu=null_gpu, StopWatch *timer_in=nullptr)
 
 AtomGraphSynthesis (const AtomGraphSynthesis &original)
 The presence of POINTER-kind Hybrid objects necessitates manually written copy and move constructors as well as copy and move assignment operators, but all are valid.
 
 AtomGraphSynthesis (AtomGraphSynthesis &&original)
 
AtomGraphSynthesisoperator= (const AtomGraphSynthesis &original)
 
AtomGraphSynthesisoperator= (AtomGraphSynthesis &&original)
 
const AtomGraphgetSystemTopologyPointer (int system_index) const
 Get a topology pointer for a specific system contained within the synthesis.
 
const std::vector< AtomGraph * > getSystemTopologyPointer () const
 
int getAtomCount () const
 Get the number of atoms in one or more systems.
 
int getAtomCount (int system_index) const
 
std::vector< AtomicRadiusSet > getPBRadiiSet () const
 Get the name of the PB radii set for one or more systems.
 
std::vector< AtomicRadiusSet > getPBRadiiSet (int low_limit, int high_limit) const
 
AtomicRadiusSet getPBRadiiSet (int index) const
 
template<typename T>
std::vector< T > getPartialCharges (HybridTargetLevel tier=HybridTargetLevel::HOST) const
 Get partial charges stored within the synthesis.
 
template<typename T>
std::vector< T > getPartialCharges (int system_index, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
template<typename T>
std::vector< T > getPartialCharges (int system_index, int low_index, int high_index, HybridTargetLevel tier=HybridTargetLevel::HOST) const
 
template<typename T>
std::vector< T > getAtomicMasses (HybridTargetLevel tier) const
 Get the masses (or inverse masses) of atoms in the synthesis.
 
template<typename T>
std::vector< T > getAtomicMasses (HybridTargetLevel tier, int system_index) const
 
template<typename T>
std::vector< T > getAtomicMasses (HybridTargetLevel tier, int system_index, int low_index, int high_index) const
 
void setImplicitSolventModel ()
 Apply an implicit solvent model to the synthesis. Any mode of operation other than taking the original topologies' parameters as is (calling this member function with no arguments) will cause changes ot the underlying topologies.
 
void setImplicitSolventModel (ImplicitSolventModel igb_in, const NeckGeneralizedBornKit< double > &ngbk, const std::vector< AtomicRadiusSet > &radii_sets_in, double dielectric_in=80.0, double saltcon_in=0.0, ExceptionResponse policy=ExceptionResponse::WARN)
 
void setImplicitSolventModel (ImplicitSolventModel igb_in, const NeckGeneralizedBornKit< double > &ngbk, AtomicRadiusSet radii_set=AtomicRadiusSet::NONE, double dielectric_in=80.0, double saltcon_in=0.0, ExceptionResponse policy=ExceptionResponse::WARN)
 
void setImplicitSolventModel (ImplicitSolventModel igb_in, const NeckGeneralizedBornTable &ngb_tab, const std::vector< AtomicRadiusSet > &radii_sets_in, double dielectric_in=80.0, double saltcon_in=0.0, ExceptionResponse policy=ExceptionResponse::WARN)
 
void setImplicitSolventModel (ImplicitSolventModel igb_in, const NeckGeneralizedBornTable &ngb_tab, AtomicRadiusSet radii_set=AtomicRadiusSet::NONE, double dielectric_in=80.0, double saltcon_in=0.0, ExceptionResponse policy=ExceptionResponse::WARN)
 

Detailed Description

A collection of one or more AtomGraph objects, with similar components arranged in contiguous arrays (often padded by the GPU warp size to prevent one system from flowing into another). Work units covering all systems are laid out in additional, contiguous arrays. While individual topologies (AtomGraphs) have Hybrid POINTER-kind objects for details such as charges and atom type indices and all of the POINTER-kind objects targeted a single ARRAY-kind object of the correct memory type, the synthesis, which may have many topologies and a great deal more overall information, stores most of its data in a series of ARRAY-kind objects, one for each member variable.

Constructor & Destructor Documentation

◆ AtomGraphSynthesis() [1/2]

stormm::synthesis::AtomGraphSynthesis::AtomGraphSynthesis ( const std::vector< AtomGraph * > & topologies_in,
const std::vector< RestraintApparatus * > & restraints_in,
const std::vector< int > & topology_indices_in,
const std::vector< int > & restraint_indices_in,
ExceptionResponse policy_in = ExceptionResponse::WARN,
const GpuDetails & gpu = null_gpu,
StopWatch * timer_in = nullptr )

The constructor takes a series of topologies and NMR restraints. The NMR restraints point to specific topologies and thereby apply to any coordinate sets that also point to those topologies.

Overloaded:

  • Prepare the energy surface with or without restraints
  • Accept lists of unique topologies and restraints with auxiliary lists of indices for replicating them, or accept explicit lists of topologies and restraints with one index for each system.
Parameters
topologies_inList of input topology pointers
restraints_inList of input restraint collections
topology_indices_inIndicators of which topologies describe each system in this synthesis. This allows a synthesis to describe many copies of each system in its work units while storing relatively small amounts of data.
restraint_indices_inIndicators of which collections of restraints supplement the energy surfaces of each system within this synthesis. Different sets of restraints allow many systems referencing the same topology to evolve on different energy surfaces.
policyInstruction on what to do if questionable input is encountered
vwu_atom_limitMaximum number of atoms to include in each valence work unit, taken from user input or a pre-computed value based on known launch bounds of the kernel that will evaluate the work units
timer_inTimer to track the wall time of this setup

◆ AtomGraphSynthesis() [2/2]

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

The presence of POINTER-kind Hybrid objects necessitates manually written copy and move constructors as well as copy and move assignment operators, but all are valid.

Parameters
originalThe original object to copy or move
otherAnother object placed on the right hand side of the assignment statement

Member Function Documentation

◆ getAtomCount()

int stormm::synthesis::AtomGraphSynthesis::getAtomCount ( ) const

Get the number of atoms in one or more systems.

Overloaded:

  • Get the total number of atoms, summed over all systems, including replicas. This does not include padding.
  • Get the total number of atoms in a specific sytem, without padding.
Parameters
system_indexindex of the system of interest

◆ getAtomicMasses()

template<typename T>
std::vector< T > stormm::synthesis::AtomGraphSynthesis::getAtomicMasses ( HybridTargetLevel tier) const

Get the masses (or inverse masses) of atoms in the synthesis.

Overloaded:

  • Get masses for all systems, padded by the warp size between systems
  • Get masses for a single system
  • Get masses for a specific range of atoms in a single system
Parameters
tierObtain the data from arrays on the host or the device
system_indexIndex of the specific system to query
low_indexLower limit of atoms in the system to query (will be validated)
high_indexUpper limit of atoms in the system to query (will be validated)

◆ getAtomOffset()

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

Get the starting point for atoms of a specific system in the lineup of all topologies.

Parameters
system_indexindex of the system of interest

◆ getDoublePrecisionAtomUpdateKit()

SyAtomUpdateKit< double, double2, double4 > stormm::synthesis::AtomGraphSynthesis::getDoublePrecisionAtomUpdateKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with double-precision real numbers for updating atom and virtual site positions based on work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getDoublePrecisionNonbondedKit()

SyNonbondedKit< double, double2 > stormm::synthesis::AtomGraphSynthesis::getDoublePrecisionNonbondedKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with double-precision real numbers for computing non-bonded interactions for all systems based on work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getDoublePrecisionRestraintKit()

SyRestraintKit< double, double2, double4 > stormm::synthesis::AtomGraphSynthesis::getDoublePrecisionRestraintKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with double-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getDoublePrecisionValenceKit()

SyValenceKit< double > stormm::synthesis::AtomGraphSynthesis::getDoublePrecisionValenceKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with double-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getPartialCharges()

template<typename T>
std::vector< T > stormm::synthesis::AtomGraphSynthesis::getPartialCharges ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get partial charges stored within the synthesis.

Overloaded:

  • Get partial charges for all systems, padded by the warp size between systems
  • Get partial charges for a single system
  • Get partial charges for a specific range of atoms in a single system
Parameters
tierObtain the data from arrays on the host or the device
system_indexIndex of the specific system to query
low_indexLower limit of atoms in the system to query (will be validated)
high_indexUpper limit of atoms in the system to query (will be validated)

◆ getPBRadiiSet()

std::vector< AtomicRadiusSet > stormm::synthesis::AtomGraphSynthesis::getPBRadiiSet ( ) const

Get the name of the PB radii set for one or more systems.

Overloaded:

  • Get the PB radii set for all systems
  • Get the PB radii set for a series of systems between low and high limits
  • Get the PB radii set for a specific system

◆ getSinglePrecisionAtomUpdateKit()

SyAtomUpdateKit< float, float2, float4 > stormm::synthesis::AtomGraphSynthesis::getSinglePrecisionAtomUpdateKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with single-precision real numbers for updating atom and virtual site positions based on work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getSinglePrecisionNonbondedKit()

SyNonbondedKit< float, float2 > stormm::synthesis::AtomGraphSynthesis::getSinglePrecisionNonbondedKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with single-precision real numbers for computing non-bonded interactions for all systems based on work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getSinglePrecisionRestraintKit()

SyRestraintKit< float, float2, float4 > stormm::synthesis::AtomGraphSynthesis::getSinglePrecisionRestraintKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with single-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getSinglePrecisionValenceKit()

SyValenceKit< float > stormm::synthesis::AtomGraphSynthesis::getSinglePrecisionValenceKit ( HybridTargetLevel tier = HybridTargetLevel::HOST) const

Get a minimal kit with single-precision parameter detail for computing valence interactions for all systems based on the work units stored in this object.

Parameters
tierLevel at which to obtain pointers for the abstract

◆ getSystemRestraintPointer()

RestraintApparatus * stormm::synthesis::AtomGraphSynthesis::getSystemRestraintPointer ( int system_index) const

Get a restraint apparatus pointer for a sepcific system contained within the synthesis.

Parameters
system_indexIndex of the system of interest

◆ getSystemTopologyPointer()

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

Get a topology pointer for a specific system contained 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

◆ loadNonbondedWorkUnits()

void stormm::synthesis::AtomGraphSynthesis::loadNonbondedWorkUnits ( const StaticExclusionMaskSynthesis & poly_se,
InitializationTask init_request = InitializationTask::NONE,
int random_cache_depth = 0,
const GpuDetails & gpu = null_gpu )

Construct non-bonded work units for all unique topologies (there are no restraints for non-bonded interactions that might distinguish systems with the same topology, as was a consideration when developing the valence work units). Load the instructions into the topology synthesis for availability on the GPU.

Parameters
poly_seSynthesis of static exclusion masks for a compilation of systems in isolated boundary conditions.
init_requestIndicate a pattern for the non-bonded work units to initialize accumulators for subsequent force and energy calculations.
random_cache_depthNumber of random values to store for each atom in the synthesis (the actual number of values stored is multiplied by three, for Cartesian X, Y, and Z contributions)
gpuDetails of the GPU in use (this may change the profile of the workload)

◆ setImplicitSolventModel()

void stormm::synthesis::AtomGraphSynthesis::setImplicitSolventModel ( )

Apply an implicit solvent model to the synthesis. Any mode of operation other than taking the original topologies' parameters as is (calling this member function with no arguments) will cause changes ot the underlying topologies.

Overloaded:

  • Take the implicit solvent models as described in the original topologies
  • Apply specific atomic radius sets to each system
  • Apply a common atomic radius set ot all systems
  • Accept a GB neck model, its double-precision abstract, or no such model
Parameters
igb_inA specific, established implicit solvent model (from literature and hard-coded into STORMM) to apply
dielectric_inThe desired dielectric constant
saltcon_inThe intended salt concentration (affects the GB decay parameter Kappa)
radii_sets_inRadii to impart to the topology (this is often coupled to the choice of implicit solvent model, but for purposes of experimentation or new model development might be flexible)
policyIndicator of what to do if the topology's PB radii to not meet the implicit solvent model requirements, or there is some other problem

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