Pteros  2.0
Molecular modeling library for human beings!
pteros::Selection Class Reference

Selection class. More...

#include <selection.h>

Classes

class  iterator
 Random-access forward iterator for Selection. More...
 

Public Member Functions

Modification of existing selection
void append (const Selection &sel)
 Append another selection to this one.
 
void append (int ind)
 Append absolute index to selection.
 
void append (const std::vector< Selection > &sel_vec)
 Append several selections to this one.
 
void remove (const Selection &sel)
 Remove all atoms of sel from current selection.
 
void remove (int ind)
 Remove given absolute index from current selection.
 
void invert ()
 Inverts selection in place by selecting those atoms which were not selected.
 
void set_system (const System &sys)
 Sets new system for selection. More...
 
void modify (std::string str, int fr=0)
 Modifies selection string in existing selection. More...
 
void modify (int ind1, int ind2)
 Modifies selection using the range of indexes.
 
void modify (const std::vector< int > &ind)
 Modifies selection using vector of indexes Vector may be in any order and may contain duplicates.
 
void modify (std::vector< int >::iterator it1, std::vector< int >::iterator it2)
 Modifies selection using pair of iterators to index vector Vector may be in any order and may contain duplicates.
 
void modify (const std::function< void(const System &, int, std::vector< int > &)> &callback, int fr=0)
 Modifies selection using user-defined callback. More...
 
void modify (const System &sys, std::string str, int fr=0)
 Convenience method, which combines set_system and modify(str)
 
void modify (const System &sys, int ind1, int ind2)
 Convenience method, which combines set_system and modify(int,int)
 
void modify (const System &sys, const std::vector< int > &ind)
 Convenience method, which combines set_system and modify(vector<int>)
 
void modify (const System &sys, std::vector< int >::iterator it1, std::vector< int >::iterator it2)
 Convenience method, which combines set_system and modify(iter,iter)
 
void modify (const System &sys, const std::function< void(const System &, int, std::vector< int > &)> &callback, int fr=0)
 Convenience method, which combines set_system and modify(callback)
 
void apply ()
 Recomputes selection without re-parsing selection text. More...
 
void update ()
 Recomputes selection completely. More...
 
void clear ()
 Clears selection and frees memory. More...
 
Sub-selections
Selection select (std::string str)
 Sub-selections allow selecting atoms inside existing selection (narrowing or refining existing selection in other terms). More...
 
Selection operator() (std::string str)
 
Selection select (int ind1, int ind2)
 Local selection indexes are used!
 
Selection operator() (int ind1, int ind2)
 
Selection select (const std::vector< int > &ind)
 Local selection indexes are used!
 
Selection operator() (const std::vector< int > &ind)
 
Iterator access
std::vector< int >::const_iterator index_begin () const
 Get const iterator for begin of index.
 
std::vector< int >::const_iterator index_end () const
 Get const iterator for the end of index.
 
std::back_insert_iterator< std::vector< int > > index_back_inserter ()
 Get back_insert_iterator for index.
 
iterator begin ()
 Begin iterator.
 
iterator end ()
 End iterator.
 
Get and set functions
int get_frame () const
 Get current frame selection is pointing to.
 
void set_frame (int fr)
 Set current frame for selection. More...
 
Systemget_system () const
 Get pointer to the system, which is pointed by this selection.
 
std::string get_text () const
 Get selection text If selection is not textual returns generated index-based text "index i j k...".
 
std::vector< int > get_index () const
 Get vector of all indexes in selection.
 
std::vector< char > get_chain (bool unique=false) const
 Get vector of all chains in selection.
 
void set_chain (const std::vector< char > &data)
 Set chains from supplied vector. More...
 
void set_chain (char data)
 Sets chain of all selected atoms to the same given value.
 
std::vector< int > get_resid (bool unique=false) const
 Get vector of all resid's in selection. More...
 
void set_resid (const std::vector< int > &data)
 Set resid's in selection from supplied vector. More...
 
void set_resid (int data)
 Sets resid of all selected atoms to the same given value.
 
std::vector< std::string > get_name (bool unique=false) const
 Get vector of all atom names in selection.
 
void set_name (const std::vector< std::string > &data)
 Set atom names in selection from supplied vector. More...
 
void set_name (std::string data)
 Sets atom names of all selected atoms to the same given value.
 
std::vector< std::string > get_resname (bool unique=false) const
 Get vector of all resnames in selection.
 
void set_resname (const std::vector< std::string > &data)
 Set resnames in selection from supplied vector. More...
 
void set_resname (std::string data)
 Sets resnames of all selected atoms to the same given value.
 
std::vector< int > get_resindex (bool unique=false) const
 Get vector of all resindexes in selection. More...
 
Eigen::MatrixXf get_xyz (bool make_row_major_matrix=false) const
 Get coordinates of all atoms in this selection for the current frame By default columnt of the matrix contain atom coordinates. More...
 
void set_xyz (MatrixXf_const_ref coord)
 Set coordinates of this selection for current frame.
 
std::vector< float > get_mass () const
 Get masses of all atoms in selection.
 
void set_mass (const std::vector< float > &m)
 Set atom masses in selection to the values from supplied vector. More...
 
void set_mass (float data)
 Sets masses of all selected atoms to the same given value.
 
std::vector< float > get_beta () const
 Get beta.
 
void set_beta (const std::vector< float > &data)
 Set beta in selection to the values from supplied vector. More...
 
void set_beta (float data)
 Sets beta of all selected atoms to the same given value.
 
std::vector< float > get_occupancy () const
 Get occupancy.
 
void set_occupancy (const std::vector< float > &data)
 Set occupancy in selection to the values from supplied vector. More...
 
void set_occupancy (float data)
 Sets occupancy of all selected atoms to the same given value.
 
std::vector< float > get_charge () const
 Get charge.
 
void set_charge (const std::vector< float > &data)
 Set charge in selection to the values from supplied vector. More...
 
void set_charge (float data)
 Sets charge of all selected atoms to the same given value.
 
float get_total_charge () const
 Computes total charge of selection.
 
std::vector< std::string > get_tag (bool unique=false) const
 Get tags.
 
void set_tag (const std::vector< std::string > &data)
 Set tags in selection to the values from supplied vector. More...
 
void set_tag (std::string data)
 Set tags of all selected atoms to the same given value.
 
Eigen::MatrixXf get_vel (bool make_row_major_matrix=false) const
 Get velocities of all atoms in this selection for the current frame By default columnt of the matrix contain atom coordinates. More...
 
void set_vel (MatrixXf_const_ref data)
 Set velocities of this selection for current frame.
 
Eigen::MatrixXf get_force (bool make_row_major_matrix=false) const
 Get forces of all atoms in this selection for the current frame By default columnt of the matrix contain atom coordinates. More...
 
void set_force (MatrixXf_const_ref coord)
 Set forces of this selection for current frame.
 
Computing properties of selection
bool is_large ()
 Checks if selection is larger than 1/2 of the box size in any dimension. More...
 
Eigen::Vector3f center (bool mass_weighted=false, Array3i_const_ref pbc=noPBC, int pbc_atom=-1) const
 Get the center of selection. More...
 
Eigen::Vector3f center (const std::vector< float > &weights, Array3i_const_ref pbc=noPBC, int pbc_atom=-1) const
 
void minmax (Vector3f_ref min, Vector3f_ref max) const
 Get minimal and maximal coordinates in selection.
 
float powersasa (float probe_r=0.14, std::vector< float > *area_per_atom=nullptr, float *total_volume=nullptr, std::vector< float > *volume_per_atom=nullptr) const
 Get the SASA using powersasa algorithm. Returns area and computes volume and per-atom values if asked.
 
float sasa (float probe_r=0.14, std::vector< float > *area_per_atom=nullptr, int n_sphere_points=960) const
 Get SASA using Shrake and Rupley algorithm (slower than powersasa and can't compute volumes)
 
Eigen::MatrixXf average_structure (int b=0, int e=-1, bool make_row_major_matrix=false) const
 Computes average structure over the range of frames.
 
Eigen::MatrixXf atom_traj (int ind, int b=0, int e=-1, bool make_row_major_matrix=false) const
 Extracts X,Y,Z for given atom index for specified range of frames (gets trajectory of given atom). More...
 
void inertia (Vector3f_ref moments, Matrix3f_ref axes, Array3i_const_ref pbc=noPBC, int pbc_atom=-1) const
 Computes the central momens of inertia and principal axes of inertia. More...
 
float gyration (Array3i_const_ref pbc=noPBC, int pbc_atom=-1) const
 Computes radius of gyration for selection. More...
 
Eigen::Vector3f dipole (bool as_charged=false, Array3i_const_ref pbc=fullPBC, int pbc_atom=-1) const
 
float distance (int i, int j, Array3i_const_ref pbc=fullPBC) const
 Get distance between two atoms (periodic in given dimensions if needed). More...
 
float angle (int i, int j, int k, Array3i_const_ref pbc=fullPBC) const
 Get angle in degrees between three atoms (periodic in given dimensions if needed). More...
 
float dihedral (int i, int j, int k, int l, Array3i_const_ref pbc=fullPBC) const
 Get dihedral angle in degrees between three atoms (periodic in given dimensions if needed). More...
 
Geometry transformation functions
void translate (Vector3f_const_ref v)
 Translate selection by given vector.
 
void translate_to (Vector3f_const_ref p, bool mass_weighted=false, Array3i_const_ref pbc=noPBC, int pbc_atom=-1)
 Translate center of masses of selection to given point.
 
void rotate (Vector3f_const_ref pivot, Vector3f_const_ref axis, float angle)
 Rotate selection around given axis relative to given pivot. More...
 
void mirror (int dim, float around=0)
 Mirror selection along one of coordinates around given point.
 
void mirror (Vector3f_const_ref normal, Vector3f_const_ref point=Eigen::Vector3f::Zero())
 Mirror selection along arbitrary vector relative to given point.
 
void wrap (Array3i_const_ref pbc=fullPBC)
 Wraps whole selection to the periodic box.
 
void unwrap (Array3i_const_ref pbc=fullPBC, int pbc_atom=-1)
 Unwraps selection to make it whole if possible (without jumps over periodic box boundary). More...
 
int unwrap_bonds (float d, Array3i_const_ref pbc=fullPBC, int pbc_atom=-1)
 Unwraps selection to make it whole (without jumps over periodic box boundary). More...
 
Eigen::Affine3f principal_transform (Array3i_const_ref pbc=noPBC, int pbc_atom=-1) const
 Get transform for orienting selection by principal axes. More...
 
void principal_orient (Array3i_const_ref pbc=noPBC, int pbc_atom=-1)
 Orient molecule by its principal axes. More...
 
int num_residues () const
 Returns total number of residues in selection. Partial residues are also included.
 
Energy functions
Eigen::Vector2f non_bond_energy (float cutoff=0, bool pbc=true) const
 Self-energy of selection computed within given interaction cut-off. More...
 
File IO
void write (std::string fname, int b=-1, int e=-1) const
 Write structure or trajectory for selection. More...
 
void write (const std::unique_ptr< FileHandler > &handler, FileContent what, int b=-1, int e=-1) const
 
Splitting selection into parts
void split_by_connectivity (float d, std::vector< Selection > &res, bool periodic=true)
 Split current selection into several selections according to the interatomic distances. More...
 
void split_by_residue (std::vector< Selection > &res)
 Split selection by residue index.
 
void split_by_molecule (std::vector< Selection > &res)
 Split selection by molecule definition from topology.
 
void split_by_chain (std::vector< Selection > &chains)
 Split selection by chain.
 
void split_by_contiguous_index (std::vector< Selection > &parts)
 Split selection into contiguous ranges of indexes.
 
void split_by_contiguous_residue (std::vector< Selection > &parts)
 Split selection into contiguous ranges of resindexes.
 
template<class F >
void split (std::vector< Selection > &parts, F callback)
 Split selection by the values returned by custom callback function. More...
 
void each_residue (std::vector< Selection > &sel) const
 Selects each residue, which is referenced by selection. More...
 
Secondary structure methods

These methods take protein residues referenced by selection and compute DSSP on them.

void dssp (std::string fname) const
 Determines secondary structure with DSSP algorithm and writes detailed report to file.
 
void dssp (std::ostream &os) const
 Determines secondary structure with DSSP algorithm and writes detailed report to stream.
 
std::string dssp () const
 Determines secondary structure with DSSP algorithm and return it as a code string. More...
 
Inline accessors

Used to access the properties of particular atom in selection.

The ind passed to these functions is is the selection index, not the global system index. I.e. passing 10 will extract the property of atom with global index index[10]
All these function except index could be used as lvalue, which means that it is possible to assign value to them. For example x(10) = 3.14

float & x (int ind)
 Extracts X for current frame.
 
float x (int ind) const
 
float & x (int ind, int fr)
 Extracts X for given frame frame fr.
 
float x (int ind, int fr) const
 
float & y (int ind)
 Extracts Y for current frame.
 
float y (int ind) const
 
float & y (int ind, int fr)
 Extracts Y for given frame frame fr.
 
float y (int ind, int fr) const
 
float & z (int ind)
 Extracts Z for current frame.
 
float z (int ind) const
 
float & z (int ind, int fr)
 Extracts Z for given frame frame fr.
 
float z (int ind, int fr) const
 
Eigen::Vector3f & xyz (int ind)
 Extracts X,Y and Z for current frame.
 
const Eigen::Vector3f & xyz (int ind) const
 
Eigen::Vector3f & xyz (int ind, int fr)
 Extracts X,Y and Z for given frame fr.
 
const Eigen::Vector3f & xyz (int ind, int fr) const
 
Eigen::Vector3f * xyz_ptr (int ind) const
 Returns pointer to the coordinates of given atom for current frame. More...
 
Eigen::Vector3f * xyz_ptr (int ind, int fr) const
 
Eigen::Vector3f & vel (int ind)
 Extracts velocity for current frame.
 
const Eigen::Vector3f & vel (int ind) const
 
Eigen::Vector3f & vel (int ind, int fr)
 Extracts velocity for given frame fr.
 
const Eigen::Vector3f & vel (int ind, int fr) const
 
Eigen::Vector3f & force (int ind)
 Extracts force for current frame.
 
const Eigen::Vector3f & force (int ind) const
 
Eigen::Vector3f & force (int ind, int fr)
 Extracts force for given frame fr.
 
const Eigen::Vector3f & force (int ind, int fr) const
 
int index (int ind) const
 Extracts type. More...
 
Atomatom (int ind)
 Extracts whole atom.
 
const Atomatom (int ind) const
 
float vdw (int ind) const
 Extracts resindex. More...
 
std::string element_name (int ind) const
 Extracts element name based on element_number. Read only.
 
PeriodicBoxbox ()
 Returns periodic box of the frame pointed by selection The same as: More...
 
const PeriodicBoxbox () const
 
float & time ()
 Returns time stamp of the frame pointed by selection The same as: More...
 
const float & time () const
 

Constructors and operators

 Selection ()
 Default constructor for empty selection.
 
 Selection (const System &sys)
 Associates selection with the system. More...
 
 Selection (const System &sys, std::string str, int fr=0)
 Textual selection constructor. More...
 
 Selection (const System &sys, int ind1, int ind2)
 Constructor, which creates selection from the interval of indexes. More...
 
 Selection (const System &sys, const std::vector< int > &ind)
 Constructor from the vector of indexes Vector may be in any order and may contain duplicates.
 
 Selection (const System &sys, std::vector< int >::iterator it1, std::vector< int >::iterator it2)
 Constructor from the pair of iterators to int sequence Sequence may be in any order and may contain duplicates.
 
 Selection (const System &sys, const std::function< void(const System &, int, std::vector< int > &)> &callback, int fr=0)
 Constructor which takes user-defined callback function. More...
 
 Selection (const Selection &other)
 Copy constructor.
 
virtual ~Selection ()
 Destructor.
 
Selectionoperator= (const Selection &other)
 Assignment operator.
 
bool operator== (const Selection &other) const
 Equality operator Selection are compared by their indexes.
 
bool operator!= (const Selection &other) const
 Inequality operator Selection are compared by their indexes.
 
AtomHandler operator[] (int ind)
 Indexing operator. More...
 
AtomHandler operator[] (const std::pair< int, int > &ind_fr)
 Indexing operator which sets both index and frame. More...
 
Selection operator~ () const
 Creates new Selection, which is a logical negation of existing one. More...
 

Fitting and RMSD functions

float rmsd (int fr) const
 RMSD between current and another frame.
 
float rmsd (int fr1, int fr2) const
 RMSD between two frames.
 
void fit_trajectory (int ref_frame=0, int b=0, int e=-1)
 Fit specified frames in the trajectory to reference frame.
 
Eigen::Affine3f fit_transform (int fr1, int fr2) const
 Returns fit transformation between frames fr1 and fr2.
 
void fit (int fr1, int fr2)
 Fits frame fr1 to fr2.
 
void apply_transform (const Eigen::Affine3f &t)
 Apply transformation.
 

Utility functions

int size () const
 Get the size of selection.
 
bool text_based () const
 Returns true if selection was created from text string and false if it was constructed 'by hand' by appending indexes or other selections.
 
bool coord_dependent () const
 Returns true if selection is coordinate-dependent and is able to recompute itself on the change of frame.
 
void flatten ()
 "Flattens" selection by removing coordinate dependence and making it not text-based. More...
 
std::string to_gromacs_ndx (std::string name) const
 Returns a string formatted as Gromacs ndx file containing the current selection with given name. More...
 
int find_index (int global_index) const
 Finds local index in selection of provided global index Returns -1 if not found.
 
std::vector< std::vector< int > > get_internal_bonds (float d, bool periodic=true) const
 Get all bonds within this selection returned as local selection indexes in form 1->[2,3,..]. More...
 

Detailed Description

Selection class.

Selections are key objects in Pteros. Technically speaking the selection is just an array, which contains indexes of selected atoms in particular system. Selection does not hold the copies of the atoms or their coordinates, it just points to them serving like an alias for certain subset of atoms. Thus selections may overlap arbitrarily. Selections are used to perform various operations on the group of selected atoms. The changes become immediately visible to all other selections, which point to some of changed atoms. Each selection is bound to particular System. There are neither 'parentless' selection nor the selections, which combine the atoms from different systems. Selections are created using the syntax, which is very similar to those used in VMD but with many useful extensions.

Constructor & Destructor Documentation

◆ Selection() [1/4]

Selection::Selection ( const System sys)
explicit

Associates selection with the system.

Parameters
sys,Selectioncontent should be set later by calling modify().

◆ Selection() [2/4]

Selection::Selection ( const System sys,
std::string  str,
int  fr = 0 
)

Textual selection constructor.

Parameters
sysSystem pointed by this selection
strSelection string

◆ Selection() [3/4]

Selection::Selection ( const System sys,
int  ind1,
int  ind2 
)

Constructor, which creates selection from the interval of indexes.

It is much faster then parsing selection corresponding string, but is limited to contigous interval of indexes.

Parameters
sysSystem pointed by this selection
ind1First index in interval
ind2Last index in interval (inclusive!)

◆ Selection() [4/4]

pteros::Selection::Selection ( const System sys,
const std::function< void(const System &, int, std::vector< int > &)> &  callback,
int  fr = 0 
)

Constructor which takes user-defined callback function.

Callback takes the system as first argument, target frame number as the second and the vector to be filled by selected atom indexes. Vector may be filled in any order and may contain duplicates.

Warning
Resulting selection is neither coordinate-dependent nor text-based. It won't recompute itself on the frame change even if it involves coordinates.

Member Function Documentation

◆ operator[]() [1/2]

AtomHandler Selection::operator[] ( int  ind)

Indexing operator.

Returns an Atom_proxy object, which incapsulates atom and its coordinates for current frame Can only be used as r-value (can not be assigned to) however individual fields of the proxy object are writable

Selection sel(sys,"name CA");
sel[1] = sel[0]; // ERROR! Can't assign to proxy object returned by []!
sel[0].name() = "A"; // OK to write fields of proxy object

Main goal of [] operator is usage in range-based for loops:

Selection sel(sys,"name CA");
for(auto a: sel){
cout << a.name() << " " << a.x() << endl;
}

Otherwise it is slower than conventional syntax like sel.name(i)

◆ operator[]() [2/2]

AtomHandler Selection::operator[] ( const std::pair< int, int > &  ind_fr)

Indexing operator which sets both index and frame.

Takes a pair of values.

Selection sel(sys,"name CA");
int at = 123;
// Assign x coorditate of atom at for frame 0 from frame 1.
sel[{at,0}].x() = sel[{at,1}].x();

◆ operator~()

Selection Selection::operator~ ( ) const

Creates new Selection, which is a logical negation of existing one.

Parent selection is not modified

◆ set_system()

void Selection::set_system ( const System sys)

Sets new system for selection.

Warning
This clears selection index and leaves it empty!

◆ modify() [1/2]

void Selection::modify ( std::string  str,
int  fr = 0 
)

Modifies selection string in existing selection.

Parameters
strNew value of selection text. Selection is re-parsed immediately with this new value.

◆ modify() [2/2]

void pteros::Selection::modify ( const std::function< void(const System &, int, std::vector< int > &)> &  callback,
int  fr = 0 
)

Modifies selection using user-defined callback.

Callback takes the system as first argument, target frame number as the second and the vector to be filled by selected atom indexes. Vector may be filled in any order and may contain duplicates.

Warning
Resulting selection is neither coordinate-dependent nor text-based. It will not recompute itself on the frame change even if it involves atom coordinates.

◆ apply()

void Selection::apply ( )

Recomputes selection without re-parsing selection text.

Only makes sense for coordinate-dependent selections when the coordinates change. Called authomatically for coordinate-dependent selections by set_frame() If selection is not coordinate-dependent does nothing.

◆ update()

void Selection::update ( )

Recomputes selection completely.

May be used when new file is loaded into the system, or when atoms are created/deleted. Forces re-parsing of selection text. For non-textual selections does nothing.

◆ clear()

void Selection::clear ( )

Clears selection and frees memory.

Selection is still assigned to its system after clearing. Subsequent call of modify() may populate it again.

◆ select()

Selection Selection::select ( std::string  str)

Sub-selections allow selecting atoms inside existing selection (narrowing or refining existing selection in other terms).

Sub-selections could be very useful in the following situation. Suppose that we need to create separate selections for N,C,CA and CB atoms of particular protein residue. With "normal" selections the following code could be used:

Selection sel_N(sys,"protein and resid 1 and name N");
Selection sel_C(sys,"protein and resid 1 and name C");
Selection sel_CA(sys,"protein and resid 1 and name CA");
Selection sel_CB(sys,"protein and resid 1 and name CB");

The problem with this code is that we are looping over all atoms in the system four times, ones in each selection. This is very inefficient since we only need to find our residue with "protein and resid 1" (one loop over all atoms) and then we need to search inside this residue four times (looping over ~10 atoms only). This problem is not apparent for small systems but becomes very painful for the systems with millions of atoms. Subselections solve this problem:

Selection residue1(sys,"protein and resid 1");
auto sel_N = residue1.select("name N");
auto sel_C = residue1.select("name C");
auto sel_CA = residue1.select("name CA");
auto sel_CB = residue1.select("name CB");

Subselections inherit the system and frame from the parent. The search in sub-selections is performed over selected atoms of the parent only (the only exception from this rule are within selections which involve seacrh over all atoms by design).

◆ set_frame()

void Selection::set_frame ( int  fr)

Set current frame for selection.

Note
If selection is coordinate-dependent it is re-evaluated by calling apply()

◆ set_chain()

void pteros::Selection::set_chain ( const std::vector< char > &  data)

Set chains from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ get_resid()

std::vector<int> pteros::Selection::get_resid ( bool  unique = false) const

Get vector of all resid's in selection.

Warning
Same resid's could be present in different chains!

◆ set_resid()

void pteros::Selection::set_resid ( const std::vector< int > &  data)

Set resid's in selection from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ set_name()

void pteros::Selection::set_name ( const std::vector< std::string > &  data)

Set atom names in selection from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ set_resname()

void pteros::Selection::set_resname ( const std::vector< std::string > &  data)

Set resnames in selection from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ get_resindex()

std::vector<int> pteros::Selection::get_resindex ( bool  unique = false) const

Get vector of all resindexes in selection.

Resindexes are unique regardless the number of the chains.

◆ get_xyz()

MatrixXf Selection::get_xyz ( bool  make_row_major_matrix = false) const

Get coordinates of all atoms in this selection for the current frame By default columnt of the matrix contain atom coordinates.

Optional flag turns this to row-major format (mainly used by python binding)

◆ set_mass()

void pteros::Selection::set_mass ( const std::vector< float > &  m)

Set atom masses in selection to the values from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ set_beta()

void pteros::Selection::set_beta ( const std::vector< float > &  data)

Set beta in selection to the values from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ set_occupancy()

void pteros::Selection::set_occupancy ( const std::vector< float > &  data)

Set occupancy in selection to the values from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ set_charge()

void pteros::Selection::set_charge ( const std::vector< float > &  data)

Set charge in selection to the values from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ set_tag()

void pteros::Selection::set_tag ( const std::vector< std::string > &  data)

Set tags in selection to the values from supplied vector.

Note
Vector size must be the save as the size of selection.

◆ get_vel()

MatrixXf Selection::get_vel ( bool  make_row_major_matrix = false) const

Get velocities of all atoms in this selection for the current frame By default columnt of the matrix contain atom coordinates.

Optional flag turns this to row-major format (mainly used by python binding)

◆ get_force()

MatrixXf Selection::get_force ( bool  make_row_major_matrix = false) const

Get forces of all atoms in this selection for the current frame By default columnt of the matrix contain atom coordinates.

Optional flag turns this to row-major format (mainly used by python binding)

◆ is_large()

bool Selection::is_large ( )

Checks if selection is larger than 1/2 of the box size in any dimension.

Returns true if so.

◆ center()

Vector3f Selection::center ( bool  mass_weighted = false,
Array3i_const_ref  pbc = noPBC,
int  pbc_atom = -1 
) const

Get the center of selection.

Parameters
mass_weightedUse mass-weighting
periodicAccount for periodic boundary conditions.
Warning
If the size of selection is larger than 1/2 of the box size in any dimension you will get incorrect results if periodic is set to true. This is not checked automatically! In this case use one of unwrapping options first.

◆ atom_traj()

MatrixXf Selection::atom_traj ( int  ind,
int  b = 0,
int  e = -1,
bool  make_row_major_matrix = false 
) const

Extracts X,Y,Z for given atom index for specified range of frames (gets trajectory of given atom).

Result is returned as MatrixXf, where i-th column (or row) is an XYZ vector for frame i.

◆ inertia()

void Selection::inertia ( Vector3f_ref  moments,
Matrix3f_ref  axes,
Array3i_const_ref  pbc = noPBC,
int  pbc_atom = -1 
) const

Computes the central momens of inertia and principal axes of inertia.

Warning
If the size of selection is larger than 1/2 of the box size in any dimension you will get incorrect results if periodic is set to true. This is not checked automatically! In this case use one of unwrapping options first.

◆ gyration()

float Selection::gyration ( Array3i_const_ref  pbc = noPBC,
int  pbc_atom = -1 
) const

Computes radius of gyration for selection.

Warning
If the size of selection is larger than 1/2 of the box size in any dimension you will get incorrect results if periodic is set to true. This is not checked automatically! In this case use one of unwrapping options first.

◆ distance()

float Selection::distance ( int  i,
int  j,
Array3i_const_ref  pbc = fullPBC 
) const

Get distance between two atoms (periodic in given dimensions if needed).

Note
This function takes selection indexes, not absolute indexes.

◆ angle()

float Selection::angle ( int  i,
int  j,
int  k,
Array3i_const_ref  pbc = fullPBC 
) const

Get angle in degrees between three atoms (periodic in given dimensions if needed).

Note
This function takes selection indexes, not absolute indexes.

◆ dihedral()

float Selection::dihedral ( int  i,
int  j,
int  k,
int  l,
Array3i_const_ref  pbc = fullPBC 
) const

Get dihedral angle in degrees between three atoms (periodic in given dimensions if needed).

Note
This function takes selection indexes, not absolute indexes.

◆ rotate()

void Selection::rotate ( Vector3f_const_ref  pivot,
Vector3f_const_ref  axis,
float  angle 
)

Rotate selection around given axis relative to given pivot.

Parameters
pivotRotation around this pivot
axisRotate around this vector
angleRotation angle in radians

◆ unwrap()

void Selection::unwrap ( Array3i_const_ref  pbc = fullPBC,
int  pbc_atom = -1 
)

Unwraps selection to make it whole if possible (without jumps over periodic box boundary).

The periodic center of mass is used as an anchor point if pbc_atom<0.

Warning
If the size of selection is larger than 1/2 of the box size in any dimension unwrap() will not work as expected and will not make selection "compact"! This is not checked automatically.

◆ unwrap_bonds()

int Selection::unwrap_bonds ( float  d,
Array3i_const_ref  pbc = fullPBC,
int  pbc_atom = -1 
)

Unwraps selection to make it whole (without jumps over periodic box boundary).

based on preserving all bonds. This method works reliably in any case, but is much slower than unwrap()

Parameters
dMaximal bond length. If 0 bonds from topology are used (if present).
pbc_atomLocal index of the reference atom, which doesn't move.
Returns
Number of disconnected pieces after unwrapping. 1 means solid selection.

◆ principal_transform()

Eigen::Affine3f Selection::principal_transform ( Array3i_const_ref  pbc = noPBC,
int  pbc_atom = -1 
) const

Get transform for orienting selection by principal axes.

Warning
If the size of selection is larger than 1/2 of the box size in any dimension you will get funny results if
Parameters
is_periodicis set to true.

◆ principal_orient()

void Selection::principal_orient ( Array3i_const_ref  pbc = noPBC,
int  pbc_atom = -1 
)

Orient molecule by its principal axes.

The same as

Eigen::Affine3f tr = sel.principal_transform();
sel.apply_transform(tr);

◆ non_bond_energy()

Vector2f Selection::non_bond_energy ( float  cutoff = 0,
bool  pbc = true 
) const

Self-energy of selection computed within given interaction cut-off.

If cutoff is 0 the cutoff from topology is used.

◆ write()

void Selection::write ( std::string  fname,
int  b = -1,
int  e = -1 
) const

Write structure or trajectory for selection.

File type is determined by extension. Frames from b to e are written. If

Parameters
bis not set or -1 it means current frame If
eis not set or -1 it means the last frame

◆ flatten()

void Selection::flatten ( )

"Flattens" selection by removing coordinate dependence and making it not text-based.

Resulting selection is equivalent to plain set of indexes "index i1 i2 i3..." Useful to avoid recomputing selection on frame change when tracking given set of atoms

◆ to_gromacs_ndx()

string Selection::to_gromacs_ndx ( std::string  name) const

Returns a string formatted as Gromacs ndx file containing the current selection with given name.

Warning
Indexes in Gromacs ndx are starting from 1! Thus one is added to all pteros indexes!

◆ get_internal_bonds()

std::vector< std::vector< int > > Selection::get_internal_bonds ( float  d,
bool  periodic = true 
) const

Get all bonds within this selection returned as local selection indexes in form 1->[2,3,..].

If d>0 it is used as cutoff if d==0 the bonds from topology are used and periodicity is ignored

◆ split_by_connectivity()

void Selection::split_by_connectivity ( float  d,
std::vector< Selection > &  res,
bool  periodic = true 
)

Split current selection into several selections according to the interatomic distances.

Each resulting selection is a group of atoms connected by distances less than d. if d==0 the bonds from topology are used instead and periodic argument is ignored

◆ split()

template<class F >
void pteros::Selection::split ( std::vector< Selection > &  parts,
callback 
)
inline

Split selection by the values returned by custom callback function.

Signature of callback is: T cb(const Selection&,int) Callback is called for each atom in selection (local index is passed as its second argument)
Callback can return any type usable as a key of std::map. Selection is split into parts according to returned values - all identical values go to the same part.

◆ each_residue()

void Selection::each_residue ( std::vector< Selection > &  sel) const

Selects each residue, which is referenced by selection.

All selections for residues are placed into supplied vector. Handles multiple chains correctly.

◆ dssp()

string Selection::dssp ( ) const

Determines secondary structure with DSSP algorithm and return it as a code string.

Returns
Code string The code is the same as in DSSP: alphahelix: 'H' betabridge: 'B' strand: 'E' helix_3: 'G' helix_5: 'I' turn: 'T' bend: 'S' loop: ' '

◆ xyz_ptr()

Eigen::Vector3f* pteros::Selection::xyz_ptr ( int  ind) const
inline

Returns pointer to the coordinates of given atom for current frame.

Used internally in Grid_searcher.

◆ index()

int pteros::Selection::index ( int  ind) const
inline

Extracts type.

Extracts typename Extracts residue name Extracts chain Extracts atom name Extracts atom mass Extracts atom charge Extracts B-factor Extracts occupancy field Extracts residue number Extracts tag Extracts atom index in the system, which is pointed by selection

◆ vdw()

float Selection::vdw ( int  ind) const

Extracts resindex.

Extracts atomic number Computes VDW radius. Read only. If atomic number is set uses whole periodic table from VMD If no atomic number uses rough guess from Gromacs vdwradii.dat

◆ box()

PeriodicBox& pteros::Selection::box ( )
inline

Returns periodic box of the frame pointed by selection The same as:

sel.get_system()->box(sel.get_frame());

This is a convenience method. The same box is returned by all selection which point to the same frame.

◆ time()

float& pteros::Selection::time ( )
inline

Returns time stamp of the frame pointed by selection The same as:

sel.get_system()->Time(sel.get_frame());

This is a convenience method. The same time is returned by all selection which point to the same frame.


The documentation for this class was generated from the following files:
pteros::Selection::Selection
Selection()
Default constructor for empty selection.
Definition: selection.cpp:89