STORMM Source Documentation
Loading...
Searching...
No Matches
stormm::topology::AtomGraphStage Class Reference

The topology stage is a CPU–bound object composed of Standard Template Library vectors for its data arrays, with each member corresponding to a (probably eponymous) member of the AtomGraph class. This object allows a topology to be built atom by atom, deconstructed, rebuilt, and even exported to a (much less mutable) AtomGraph. More...

#include <atomgraph_stage.h>

Public Member Functions

const std::string & getTitle () const
 Get the title of the topology as it is staged.
 
const std::string & getOutputFile () const
 Get the name of the output file to be referenced by the resulting topology. Printed results will go to this file.
 
int getAtomCount () const
 Get the total number of atoms.
 
int getResidueCount () const
 Get the total number of residues.
 
int getMoleculeCount () const
 Get the total number of molecules.
 
double getElectrostatic14Screening () const
 Get the general 1:4 screening factor on electrostatic interactions.
 
double getVanDerWaals14Screening () const
 Get the general 1:4 screening factor on van-der Waals interactions.
 
AtomGraph exportTopology (int minimum_solute_size=default_minimum_solute_atoms, const std::vector< char4 > &solute_included_residue_names={}, double hmass_repartition_factor=1.0) const
 Export the topology based on the internal data arrays. A final list of checks and fixes will be applied to ensure consistent parameterization.
 
void setTitle (const std::string &title_in)
 Set the title of the topology.
 
void setOutputFile (const std::string &output_file_in)
 Set the output file.
 
void addUbrdParameters (double k_eq_in, double l_eq_in, const char4 atom_i_type=unknown_atom_type, const char4 atom_k_type=unknown_atom_type)
 Add new Urey-Bradley stretching parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.
 
void addCImpParameters (double k_eq_in, double phi_in, const char4 atom_i_type=unknown_atom_type, const char4 atom_j_type=unknown_atom_type, const char4 atom_k_type=unknown_atom_type, const char4 atom_l_type=unknown_atom_type)
 Add new parameters for a CHARMM improper dihedral, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.
 
void addCmapParameters (const int dimension_in, const std::vector< double > &surf_in, const char4 atom_i_type=unknown_atom_type, const char4 atom_j_type=unknown_atom_type, const char4 atom_k_type=unknown_atom_type, const char4 atom_l_type=unknown_atom_type, const char4 atom_m_type=unknown_atom_type)
 Add new parameters for a CHARMM correction map two-dimensional surface term, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.
 
void addBondParameters (double k_eq_in, double l_eq_in, const char4 atom_i_type=unknown_atom_type, const char4 atom_j_type=unknown_atom_type)
 Add new bond parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.
 
void addAnglParameters (double k_eq_in, double th_eq_in, const char4 atom_i_type=unknown_atom_type, const char4 atom_j_type=unknown_atom_type, const char4 atom_k_type=unknown_atom_type)
 Add new angle parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.
 
void addDiheParameters (double ampl_in, double period_in, double phi_in, const char4 atom_i_type=unknown_atom_type, const char4 atom_j_type=unknown_atom_type, const char4 atom_k_type=unknown_atom_type, const char4 atom_l_type=unknown_atom_type)
 Add new dihedral parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.
 
std::vector< int > addAtoms (const std::vector< int > &z_numbers, int placement=-1, const std::vector< int > &added_residue_limits={}, const std::vector< char4 > &added_atom_names={}, const std::vector< char4 > &added_atom_types={}, const std::vector< char4 > &added_residue_names={}, const std::vector< int > &added_charge_indices={}, const std::vector< int > &added_lennard_jones_indices={}, const std::vector< double > &added_atomic_charges={})
 Add atoms to the system. This will create space for new atoms at the end of each array or at the position specified, moving all existing atoms back and promoting their topological indices as needed. Residue and molecule limits will be adjust such that the added atoms become part of the residue containing the placement point (including placement points at the tail of the residue). If the residue limits are [0, 5), [5, 14), and [14, 18) and six atoms are added at placement 5, the first residue gains the atoms and the rsidue limits become [0, 11), [11, 20), and [20, 24). An array of indices indicating the positions at which the atoms were added will be returned.
 
int addAtom (int z_number, int placement=-1, const char4 atom_name={ 'U', 'N', 'K', ' ' }, const char4 atom_type={ 'U', 'N', 'K', ' ' }, const char4 res_name={ 'U', 'N', 'K', ' ' }, int charge_index=-1, int lennard_jones_index=-1, double atomic_charge=0.0)
 Add an atom to the topology. This is a single-atom variant of the addAtoms() function above, which it calls internally, and returns a single index where the atom was placed.
 
std::vector< int > addVirtualSites (const std::vector< int > &parameter_indices, const std::vector< int > &frame1_atoms, const std::vector< int > &frame2_atoms, const std::vector< int > &frame3_atoms, const std::vector< int > &frame4_atoms)
 Add a series of virtual sites to a topology. A frame parameter set of each specified index must be added in advance, and any frame atoms referenced by the sites must also exist in advance. If other atoms or other virtual sites are added later, the atom indexing of existing virtual sites will update along with the growing topology. The indices of the added virtual sites in the topology are returned so that future additions, whether of virtual sites or other atoms, can be adjusted accordingly.
 
int addVirtualSite (int parameter_index, int frame1_atom, int frame2_atom, int frame3_atom=-1, int frame4_atom=-1)
 Add a virtual site to the topology. This is a single-particle variant of the addVirtualSites() function above, which it calls internally, and returns a single index where the virtual site was placed.
 
void setBonds (const std::vector< int > &atom_i, const std::vector< int > &atom_j, const std::vector< int > &parameter_index)
 Set multiple bonds throughout the structure. The bond term lists, bond parameter lists, and number of bonds will be updated as appropriate.
 
void setBond (int atom_i, int atom_j, int parameter_index)
 Set a bond between two atoms. The bond term lists, bond parameter lists, and number of bonds will be updated as appropriate.
 
void setElectrostatic14Screening (double screening_factor_in)
 Set the general screening factor on 1:4 electrostatic interactions.
 
void setVanDerWaals14Screening (double screening_factor_in)
 Set the general screening factor on 1:4 van-der Waals interactions.
 
 AtomGraphStage (int atom_count_in=0, const std::vector< int > &residue_limits_in={}, ExceptionResponse policy_in=ExceptionResponse::DIE)
 The constructor can accept a number of atoms, residues, and molecules, an existing AtomGraph, or a ChemicalFeatures object.
 
 AtomGraphStage (const AtomGraph *ag_in, ExceptionResponse policy_in=ExceptionResponse::DIE)
 
 AtomGraphStage (const AtomGraph &ag_in, ExceptionResponse policy_in=ExceptionResponse::DIE)
 
 AtomGraphStage (const MdlMol *cmpd_in, ExceptionResponse policy_in=ExceptionResponse::DIE)
 
 AtomGraphStage (const MdlMol &cmpd_in, ExceptionResponse policy_in=ExceptionResponse::DIE)
 
 AtomGraphStage (const AtomGraphStage &original)=default
 The utility of this object lies in its simple construction and independent arrays, which can be resized as needed. This implies no pointers to repair and no POINTER-kind Hybrid objects which could require special copy or move constructors and assignment operators. All of the defaults apply.
 
 AtomGraphStage (AtomGraphStage &&original)=default
 
AtomGraphStageoperator= (const AtomGraphStage &other)=default
 
AtomGraphStageoperator= (AtomGraphStage &&other)=default
 

Detailed Description

The topology stage is a CPU–bound object composed of Standard Template Library vectors for its data arrays, with each member corresponding to a (probably eponymous) member of the AtomGraph class. This object allows a topology to be built atom by atom, deconstructed, rebuilt, and even exported to a (much less mutable) AtomGraph.

Constructor & Destructor Documentation

◆ AtomGraphStage()

stormm::topology::AtomGraphStage::AtomGraphStage ( const AtomGraphStage & original)
default

The utility of this object lies in its simple construction and independent arrays, which can be resized as needed. This implies no pointers to repair and no POINTER-kind Hybrid objects which could require special copy or move constructors and assignment operators. All of the defaults apply.

Parameters
originalThe object to copy or move
otherAn existing object placed on the left hand side of the assignment statement

Member Function Documentation

◆ addAnglParameters()

void stormm::topology::AtomGraphStage::addAnglParameters ( double k_eq_in,
double th_eq_in,
const char4 atom_i_type = unknown_atom_type,
const char4 atom_j_type = unknown_atom_type,
const char4 atom_k_type = unknown_atom_type )

Add new angle parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.

Parameters
k_eq_inForce constant (stiffness) for the harmonic bond angle interaction, in units of kcal/mol-radian^2
th_eq_inEquilibrium angle, in radians
atom_i_typeNon-bonded atom type of the first atom in the angle. The atom types are reflexive, an angle spanning types tI, tJ, and tK being equivalent to an angle between types tK, tJ, and tI.
atom_j_typeNon-bonded atom type of the second atom in the angle
atom_k_typeNon-bonded atom type of the third atom in the angle

◆ addAtom()

int stormm::topology::AtomGraphStage::addAtom ( int z_number,
int placement = -1,
const char4 atom_name = { 'U', 'N', 'K', ' ' },
const char4 atom_type = { 'U', 'N', 'K', ' ' },
const char4 res_name = { 'U', 'N', 'K', ' ' },
int charge_index = -1,
int lennard_jones_index = -1,
double atomic_charge = 0.0 )

Add an atom to the topology. This is a single-atom variant of the addAtoms() function above, which it calls internally, and returns a single index where the atom was placed.

Parameters
z_numberAtomic number of the atom to add, required to specify the element
placementTopological index at which to insert the atom into the topology as it exists when the function is called. If negative or beyond the end of the list, the atom will be palced at the end of the list.
atom_nameName of the atom
atom_typeName of the non-bonded atom type
res_nameName of the residue that the atom is to be part of, if the atom is seeding a new residue
charge_indexPartial charge parameter index of the atom
lennard_jones_indexLennard-Jones parameter index of the atom
atomic_chargePartial charge of the atom

◆ addAtoms()

std::vector< int > stormm::topology::AtomGraphStage::addAtoms ( const std::vector< int > & z_numbers,
int placement = -1,
const std::vector< int > & added_residue_limits = {},
const std::vector< char4 > & added_atom_names = {},
const std::vector< char4 > & added_atom_types = {},
const std::vector< char4 > & added_residue_names = {},
const std::vector< int > & added_charge_indices = {},
const std::vector< int > & added_lennard_jones_indices = {},
const std::vector< double > & added_atomic_charges = {} )

Add atoms to the system. This will create space for new atoms at the end of each array or at the position specified, moving all existing atoms back and promoting their topological indices as needed. Residue and molecule limits will be adjust such that the added atoms become part of the residue containing the placement point (including placement points at the tail of the residue). If the residue limits are [0, 5), [5, 14), and [14, 18) and six atoms are added at placement 5, the first residue gains the atoms and the rsidue limits become [0, 11), [11, 20), and [20, 24). An array of indices indicating the positions at which the atoms were added will be returned.

Parameters
z_numbersAtomic numbers of the atoms to add
placementTopological index to give the first added atom. Negative values will add atoms to the back of the list. Other added atoms will be inserted in sequence.
added_residue_limitsOptional definition of the boundaries of residues within the added atom set. Indices apply to the list of added atoms, not the topology itself, although they will affect how the residue bounds of the topology develop after the additions.
added_atom_namesOptional names of the added atoms
added_atom_typesOptional charater codes for added atom non-bonded types
added_residue_namesOptional names of the residues composed by the added atoms. Must meet the dimension of the added residue limits array.
added_charge_indicesOptional charge parameter indices for added atoms
added_atomic_chargesOptional partial charges for added atoms
added_lennard_jones_indicesOptional Lennard-Jones parameter indices to impart to the added atoms

◆ addBondParameters()

void stormm::topology::AtomGraphStage::addBondParameters ( double k_eq_in,
double l_eq_in,
const char4 atom_i_type = unknown_atom_type,
const char4 atom_j_type = unknown_atom_type )

Add new bond parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.

Parameters
k_eq_inForce constant (stiffness) for the harmonic bond interaction, in units of kcal/mol-A^2
l_eq_inEquilibrium length of the bond, in Angstroms
atom_i_typeNon-bonded atom type of the first atom in the bond. The atom types are reflexive, a bond between types tI and tJ being equivalent to a bond between types tJ and tI.
atom_j_typeNon-bonded atom type of the second atom in the bond

◆ addCImpParameters()

void stormm::topology::AtomGraphStage::addCImpParameters ( double k_eq_in,
double phi_in,
const char4 atom_i_type = unknown_atom_type,
const char4 atom_j_type = unknown_atom_type,
const char4 atom_k_type = unknown_atom_type,
const char4 atom_l_type = unknown_atom_type )

Add new parameters for a CHARMM improper dihedral, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.

Parameters
ampl_inForce constant (amplitude) for the harmonic improper torsion term, in units of kcal/mol-rad^2
phi_inPhase angle for the improper term, in unis of radians
atom_i_typeNon-bonded atom type of the first atom in the improper torsion. The atom types are again reflexive, a harmonic improper spanning types tI, tJ, tK, and tL being equivalent to an angle between types tL, tK, tJ, and tI.
atom_j_typeNon-bonded atom type of the second atom in the improper torsion
atom_k_typeNon-bonded atom type of the third atom in the improper torsion
atom_l_typeNon-bonded atom type of the fourth atom in the improper torsion

◆ addCmapParameters()

void stormm::topology::AtomGraphStage::addCmapParameters ( const int dimension_in,
const std::vector< double > & surf_in,
const char4 atom_i_type = unknown_atom_type,
const char4 atom_j_type = unknown_atom_type,
const char4 atom_k_type = unknown_atom_type,
const char4 atom_l_type = unknown_atom_type,
const char4 atom_m_type = unknown_atom_type )

Add new parameters for a CHARMM correction map two-dimensional surface term, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.

Parameters
dimension_inDimension of the two-dimensional surface, assumed to be square and periodic on the range [-pi, pi)
surf_inDetails of the two-dimensional energy surface, which will be solved by Catmull-Rom splines wherein the first and second cross-derivatives at each grid point are set to be equal between segements which exactly interpolate the potential values as provided. Units of kcal/mol.
atom_i_typeNon-bonded atom type of the first atom in the CMAP. The atom types are again reflexive, a harmonic improper spanning types tI, tJ, tK, tL, and tM being equivalent to an angle between types tM, tL, tK, tJ, and tI.
atom_j_typeNon-bonded atom type of the second atom in the CMAP
atom_k_typeNon-bonded atom type of the third atom in the CMAP
atom_l_typeNon-bonded atom type of the fourth atom in the CMAP
atom_m_typeNon-bonded atom type of the fifth atom in the CMAP

◆ addDiheParameters()

void stormm::topology::AtomGraphStage::addDiheParameters ( double ampl_in,
double period_in,
double phi_in,
const char4 atom_i_type = unknown_atom_type,
const char4 atom_j_type = unknown_atom_type,
const char4 atom_k_type = unknown_atom_type,
const char4 atom_l_type = unknown_atom_type )

Add new dihedral parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.

Parameters
ampl_inForce constant (amplitude) for the cosine-based torsion term, in units of kcal/mol
period_inPeriodicity of the cosine argument (unitless)
phi_inPhase angle for the cosine term, in unis of radians
atom_i_typeNon-bonded atom type of the first atom in the torsion. The atom types are again reflexive, a torsion spanning types tI, tJ, tK, and tL being equivalent to an angle between types tL, tK, tJ, and tI.
atom_j_typeNon-bonded atom type of the second atom in the torsion
atom_k_typeNon-bonded atom type of the third atom in the torsion
atom_l_typeNon-bonded atom type of the fourth atom in the torsion

◆ addUbrdParameters()

void stormm::topology::AtomGraphStage::addUbrdParameters ( double k_eq_in,
double l_eq_in,
const char4 atom_i_type = unknown_atom_type,
const char4 atom_k_type = unknown_atom_type )

Add new Urey-Bradley stretching parameters, accessible either by an index determined by the order of addition or non-bonded atom types of the participating atoms.

Parameters
k_eq_inForce constant (stiffness) for the harmonic bond interaction, in units of kcal/mol-A^2
l_eq_inEquilibrium length of the bond, in Angstroms
atom_i_typeNon-bonded atom type of the first atom in the term. The atom types are reflexive, a Urey-Bradley stretching term between types tI and tK being equivalent to a term between types tK and tI.
atom_k_typeNon-bonded atom type of the second atom in the bond

◆ addVirtualSite()

int stormm::topology::AtomGraphStage::addVirtualSite ( int parameter_index,
int frame1_atom,
int frame2_atom,
int frame3_atom = -1,
int frame4_atom = -1 )

Add a virtual site to the topology. This is a single-particle variant of the addVirtualSites() function above, which it calls internally, and returns a single index where the virtual site was placed.

Parameters
parameter_indexIndex of the virtual site parameter set to apply (this must match the number of frame atoms provided)
frame1_atomThe parent atom of the virtual site
frame2_atomSecond atom in the frame (all site types have at least two atoms)
frame3_atomOptional third atom in the frame
frame4_atomOptional fourth atom in the frame

◆ addVirtualSites()

std::vector< int > stormm::topology::AtomGraphStage::addVirtualSites ( const std::vector< int > & parameter_indices,
const std::vector< int > & frame1_atoms,
const std::vector< int > & frame2_atoms,
const std::vector< int > & frame3_atoms,
const std::vector< int > & frame4_atoms )

Add a series of virtual sites to a topology. A frame parameter set of each specified index must be added in advance, and any frame atoms referenced by the sites must also exist in advance. If other atoms or other virtual sites are added later, the atom indexing of existing virtual sites will update along with the growing topology. The indices of the added virtual sites in the topology are returned so that future additions, whether of virtual sites or other atoms, can be adjusted accordingly.

Parameters
parameter_indicesIndices of the virtual site parameter set to apply (this must match the number of frame atoms provided)
frame1_atomsThe parent atoms of each virtual site
frame2_atomsSecond atoms in each virtual site frame (all site types have at least two atoms)
frame3_atomsOptional third atoms in the virtual site frames
frame4_atomsOptional fourth atoms in the virtual site frames

◆ exportTopology()

AtomGraph stormm::topology::AtomGraphStage::exportTopology ( int minimum_solute_size = default_minimum_solute_atoms,
const std::vector< char4 > & solute_included_residue_names = {},
double hmass_repartition_factor = 1.0 ) const

Export the topology based on the internal data arrays. A final list of checks and fixes will be applied to ensure consistent parameterization.

Parameters
minimum_solute_sizeThe minimum size of a solute molecule. Anything larger will be automatically included in the solute.
solute_included_residues_namesNames of additional residues that should be included in the solute
hmass_repartition_factorHydrogen mass repartitioning factor. The mass of hydrogen atoms in the resulting topology will be scaled by this amount. The default of 1.0 implies no repartitioning. Simulations seeking to use a 4fs time step typically repartition with a factor of 3.0. The extra mass added to hydrogen atoms will be borrowed from their parent heavy atoms.

◆ setBond()

void stormm::topology::AtomGraphStage::setBond ( int atom_i,
int atom_j,
int parameter_index )

Set a bond between two atoms. The bond term lists, bond parameter lists, and number of bonds will be updated as appropriate.

Parameters
atom_iThe first atom in the bond
atom_jThe second atom in the bond
parameter_indexIndex of the parameter set to use for the bond

◆ setBonds()

void stormm::topology::AtomGraphStage::setBonds ( const std::vector< int > & atom_i,
const std::vector< int > & atom_j,
const std::vector< int > & parameter_index )

Set multiple bonds throughout the structure. The bond term lists, bond parameter lists, and number of bonds will be updated as appropriate.

Parameters
atom_iThe first atom in each of the bonds
atom_jThe second atom in each of the bonds
parameter_indexIndices of the parameter sets to use for each bond

◆ setElectrostatic14Screening()

void stormm::topology::AtomGraphStage::setElectrostatic14Screening ( double screening_factor_in)

Set the general screening factor on 1:4 electrostatic interactions.

Parameters
screening_factor_inThe screening factor to apply

◆ setOutputFile()

void stormm::topology::AtomGraphStage::setOutputFile ( const std::string & output_file_in)

Set the output file.

Parameters
output_file_inName of the output file to associate with the topology, the file to write

◆ setTitle()

void stormm::topology::AtomGraphStage::setTitle ( const std::string & title_in)

Set the title of the topology.

Parameters
title_inThe title to apply (will be checked for length, max 80 characters)

◆ setVanDerWaals14Screening()

void stormm::topology::AtomGraphStage::setVanDerWaals14Screening ( double screening_factor_in)

Set the general screening factor on 1:4 van-der Waals interactions.

Parameters
screening_factor_inThe screening factor to apply

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