LOOS 4.1.0
The Lightweight Object Oriented Structural analysis library/toolkit
|
Class for handling groups of Atoms (pAtoms, actually) More...
#include <AtomicGroup.hpp>
Classes | |
class | Iterator |
Our own simple iterator for stepping over all managed atoms. More... | |
Public Types | |
typedef std::vector< pAtom >::iterator | iterator |
typedef std::vector< pAtom >::const_iterator | const_iterator |
typedef pAtom | value_type |
Public Member Functions | |
AtomicGroup (const int n) | |
Creates a new AtomicGroup with n un-initialized atoms. | |
AtomicGroup (const AtomicGroup &g) | |
Copy constructor (atoms and box shared) | |
AtomicGroup | copy (void) const |
Creates a deep copy of this group. | |
virtual AtomicGroup * | clone (void) const |
Creates a lightweight clone of this group (for polymorphism) | |
uint | length (void) const |
uint | size (void) const |
bool | empty (void) const |
pAtom | getAtom (const int i) const |
Get the ith atom from this group. | |
pAtom & | operator[] (const int i) |
Same as getAtom(i) | |
const pAtom & | operator[] (const int i) const |
AtomicGroup & | append (pAtom pa) |
Append the atom onto the group. | |
AtomicGroup & | append (std::vector< pAtom > pas) |
Append a vector of atoms. | |
AtomicGroup & | append (const AtomicGroup &grp) |
Append an entire AtomicGroup onto this one (concatenation) | |
AtomicGroup & | remove (pAtom pa) |
Delete a single atom. | |
AtomicGroup & | remove (std::vector< pAtom > pas) |
Deletes a set of atoms. | |
AtomicGroup & | remove (const AtomicGroup &grp) |
Deletes all atoms in the passed grp that are also in the current group. | |
AtomicGroup & | operator+= (const AtomicGroup &rhs) |
AtomicGroup & | operator+= (const pAtom &rhs) |
AtomicGroup | operator+ (const AtomicGroup &rhs) |
AtomicGroup | operator+ (const pAtom &rhs) |
bool | operator== (AtomicGroup &rhs) |
Equality test for two groups. | |
bool | operator!= (AtomicGroup &rhs) |
Inequality test for two groups. | |
bool | operator== (const AtomicGroup &rhs) const |
Equality test for const groups. | |
bool | operator!= (const AtomicGroup &rhs) const |
Inequality test for two groups. | |
AtomicGroup | subset (const int offset, const int len=0) |
subset() and excise() args are patterned after perl's substr... | |
AtomicGroup | excise (const int offset, const int len=0) |
excise returns the excised atoms as a group... | |
template<class EqualsOp > | |
bool | contains (const pAtom &p, const EqualsOp &op) const |
Determines if a pAtom is contained in this group using the EqualsOp atom-equality policy. | |
bool | contains (const pAtom &p) const |
Determines if a pAtom is contained in this group using the AtomEquals policy (ie the default comparison policy) | |
template<class EqualsOp > | |
bool | contains (const AtomicGroup &g, const EqualsOp &op) const |
Determines if the passed group is a subset of the current group using the EqualsOp atom-equality policy. | |
bool | contains (const AtomicGroup &g) const |
Determines if a group is a subset of the current group using the default AtomEquals policy. | |
template<class EqualsOp > | |
bool | containsAny (const AtomicGroup &g, const EqualsOp &op) const |
Determines if a group contains any atom. | |
bool | containsAny (const AtomicGroup &g) const |
Determines if a group contains any atom using the default AtomEquals policy. | |
template<class EqualsOp > | |
AtomicGroup | intersect (const AtomicGroup &g, const EqualsOp &op) |
Computes the intersection of two groups using the EqualsOp atom-equality policy. | |
AtomicGroup | intersect (const AtomicGroup &g) |
Intersection of two groups. | |
template<class EqualsOp > | |
AtomicGroup | merge (const AtomicGroup &g, const EqualsOp &op) |
Union of two groups using the specified atom-equality policy. | |
AtomicGroup | merge (const AtomicGroup &g) |
Union of two groups using the default AtomEquals atom-equality policy. | |
AtomicGroup | select (const AtomSelector &sel) const |
Return a group consisting of atoms for which sel predicate returns true... | |
std::vector< AtomicGroup > | splitByUniqueSegid (void) const |
Returns a vector of AtomicGroups split from the current group based on segid. | |
std::vector< AtomicGroup > | splitByMolecule (void) const |
Returns a vector of AtomicGroups split based on bond connectivity. | |
std::vector< AtomicGroup > | splitByMolecule (const std::string &selection) const |
std::vector< AtomicGroup > | splitByResidue (void) const |
Returns a vector of AtomicGroups, each comprising a single residue. | |
std::map< std::string, AtomicGroup > | splitByName (void) const |
Returns a vector of AtomicGroups, each containing atoms with the same name. | |
AtomicGroup | centrifyByMolecule () const |
Replace a group with the center of masses of contained molecules. | |
AtomicGroup | centrifyByResidue () const |
Replace a group with the cente of masses of contained residues (see centrifyByMolecule()) | |
pAtom | findById (const int id) const |
Find a contained atom by its atomid. | |
AtomicGroup | groupFromID (const std::vector< int > &id_list) const |
Create a new group from a vector of atomids. | |
AtomicGroup | groupFromID (const std::pair< int, int > &id_pair) const |
Create a new group from a pair of atomids. | |
AtomicGroup | getResidue (pAtom res) |
std::string | asString () const |
void | renumber (const int start=1, const int stride=1) |
Renumber the atomid's of the contained atoms... | |
int | minId (void) const |
int | maxId (void) const |
int | minResid (void) const |
int | maxResid (void) const |
int | numberOfResidues (void) const |
int | numberOfSegids (void) const |
bool | allHaveProperty (const Atom::bits &property) const |
True if all atoms in the group have the passed property(ies) | |
bool | anyHaveProperty (const Atom::bits &property) const |
True if any atom in the group have the passed property(ies) | |
bool | hasBonds (void) const |
Does any atom in the group have bond information??? | |
bool | hasCoords (void) const |
Do all the atoms in the group have coordinates? | |
void | clearBonds (void) |
Remove any bonding information present in contained atoms. | |
void | pruneBonds () |
Attempt to prune connectivity (only retain bonds to atoms within this AtomicGroup) | |
void | resetAtomIndices () |
Reset the atom indices (used for interfacing with trajectories) | |
uint | deduceAtomicNumberFromMass (const double tol=0.1) |
Deduce atomic number from mass (if present), returning number of atoms assigned. | |
double | deduceMassFromAtomicNumber () |
Deduce mass from atomic number (if present), returning the number of atoms assigned. | |
bool | sorted (void) const |
Is the array of atoms already sorted??? | |
void | sort (void) |
Sort based on atomid. | |
bool | isPeriodic (void) const |
Test whether or not periodic boundary conditions are set. | |
GCoord | periodicBox (void) const |
Fetch the periodic boundary conditions. | |
void | periodicBox (const GCoord &c) |
Set the periodic boundary conditions. | |
void | periodicBox (const greal x, const greal y, const greal z) |
Set the periodic boundary conditions. | |
const greal | ocf (uint offset) |
compute OCF for all atom-pairs in AG of distance offset from one another | |
loos::SharedPeriodicBox | sharedPeriodicBox () const |
Provide access to the underlying shared periodic box... | |
void | removePeriodicBox () |
Remove periodicity. | |
void | reimage () |
void | reimageByAtom () |
Reimage atoms individually into the primary cell. | |
void | mergeImage (pAtom &p) |
Takes a group that's split across a periodic boundary and reimages it so it's all together. | |
void | mergeImage () |
Takes a group that's split across a periodic boundary and reimages it so it's all together, using the first atom in the AtomicGroup as the reference. | |
AtomicGroup | within (const double dist, AtomicGroup &grp) const |
Find atoms in the current group that are within dist angstroms of any atom in grp. | |
AtomicGroup | within (const double dist, AtomicGroup &grp, const GCoord &box) const |
Find atoms in grp that are within dist angstroms of atoms in the current group, considering periodicity. | |
bool | contactWith (const double dist, const AtomicGroup &grp, const uint min=1) const |
Returns true if any atom of current group is within dist angstroms of grp. | |
bool | contactWith (const double dist, const AtomicGroup &grp, const GCoord &box, const uint min=1) const |
Returns true if any atom of current group is within dist angstroms of grp. | |
std::vector< std::pair< int, int > > | getBondsIDs () const |
return a list of atom ID pairs that correspond to all unique bonds. | |
std::vector< AtomicGroup > | getBondsAGs () const |
return a list of atom index pairs corresponding to all unique bonds. | |
void | findBonds (const double dist, const GCoord &box) |
Distance-based search for bonds. | |
void | findBonds (const double dist) |
void | findBonds (const GCoord &box) |
void | findBonds () |
template<class T > | |
T | apply (T func) |
Apply a functor or a function to each atom in the group. | |
iterator | begin (void) |
iterator | end (void) |
const_iterator | begin (void) const |
const_iterator | end (void) const |
std::vector< GCoord > | boundingBox (void) const |
Bounding box for the group... | |
GCoord | centerAtOrigin (void) |
Translates the group so that the centroid is at the origin. | |
GCoord | centroid (void) const |
Centroid of atoms (ignores mass, operates in group coordinates) | |
greal | radius (const bool use_atom_as_reference=false) const |
Maximum radius from centroid of all atoms (not gyration) | |
GCoord | centerOfMass (void) const |
Center of mass of the group (in group coordinates) | |
GCoord | centerOfElectrons (void) const |
Analogous to center of mass. | |
GCoord | dipoleMoment (void) const |
Dipole moment, relative to group's centroid. | |
greal | totalCharge (void) const |
greal | totalMass (void) const |
greal | radiusOfGyration (void) const |
greal | sphericalVariance (const pAtom) const |
Spherical variance of group with respect to target atom. | |
greal | sphericalVariance (const GCoord) const |
greal | stacking (const AtomicGroup &, const GCoord &box, const double threshold) const |
Estimate stacking, as between two nucleobases. | |
greal | rmsd (const AtomicGroup &) |
Compute the RMSD between two groups. | |
greal | kineticEnergy () |
Compute kinetic energy of group. | |
std::vector< GCoord > | getTransformedCoords (const XForm &) const |
Returns a vector of coordinates transformed by the passed XForm. | |
std::vector< GCoord > | differenceVectors (const AtomicGroup &other) |
Compute difference vectors between two AtomicGroups. | |
void | translate (const GCoord &v) |
Translate an atomic group by vector v. | |
void | rotate (const GCoord &axis, const greal angle_in_degrees) |
Rotate group's coordinates (right-handed, about centroid) | |
void | applyTransform (const XForm &) |
Apply the given transform to the group's coordinates... | |
void | copyCoordinatesWithIndex (const std::vector< GCoord > &coords) |
Copy coordinates from a vector of GCoords using the atom index as an index into the vector. | |
void | copyVelocitiesWithIndex (const std::vector< GCoord > &velocities) |
Copy velocities from a vector of GCoords using the atom index as an index into the vector. | |
void | copyCoordinatesFrom (const AtomicGroup &g, const uint offset=0, const uint length=0) |
Copy coordinates from g into current group. | |
std::vector< uint > | atomOrderMapFrom (const AtomicGroup &g) |
Map the order of atoms in AtomicGroup g into the current group. | |
void | copyMappedCoordinatesFrom (const AtomicGroup &g, const std::vector< uint > &order) |
Given a mapping of atom order, copy the coordinates into the current group. | |
void | copyMappedCoordinatesFrom (const AtomicGroup &g) |
Copy the coordinates from the group mapping the atom order. | |
void | perturbCoords (const greal) |
Each atom is moved in a random direction by a vector of the passed size. | |
std::vector< GCoord > | principalAxes (void) const |
Compute the principal axes of a group. | |
double | principalAxesOrder (void) const |
Computes order parameter based on principalAxes. | |
std::vector< GCoord > | momentsOfInertia (void) const |
Computes the moments of inertia for a group. | |
GMatrix | superposition (const AtomicGroup &) |
Calculates the transformation matrix for superposition of groups. | |
GMatrix | alignOnto (const AtomicGroup &) |
Superimposes the current group onto the passed group. | |
void | orientAlong (const GCoord &) |
Orient the principal axis of this group along the supplied vector. | |
void | setCoords (double *seq, int m, int n) |
void | getCoords (double **outseq, int *m, int *n) |
std::vector< double > | coordsAsVector () const |
double | packingScore (const AtomicGroup &other, const GCoord &box, bool norm) const |
double | logisticContact (const AtomicGroup &group, double radius, int sigma, const GCoord &box) const |
double | logisticContact2D (const AtomicGroup &group, double radius, int sigma, const GCoord &box) const |
double | hardContact (const AtomicGroup &group, double radius, const GCoord &box) const |
double | hardContact2D (const AtomicGroup &group, double radius, const GCoord &box) const |
std::vector< double > | scattering (const double qmin, const double max, const uint numValues, loos::FormFactorSet &formFactors) |
Static Public Attributes | |
static const double | superposition_zero_singular_value = 1e-10 |
Protected Member Functions | |
void | setGroupConnectivity () |
Protected Attributes | |
std::vector< pAtom > | atoms |
loos::SharedPeriodicBox | box |
Friends | |
std::ostream & | operator<< (std::ostream &os, const AtomicGroup &grp) |
Output the group in pseudo-XML format... | |
Class for handling groups of Atoms (pAtoms, actually)
This class contains a collection of shared pointers to Atoms (i.e. pAtoms). Copying an AtomicGroup is a light-copy. You can, however, perform a deep copy by using the AtomicGroup::copy() method. Note that atomid's are assumed to be unique for any given AtomicGroup.
Valid operators are '+' and '+=' and can combine either AtomicGroup objects or pAtom objects.
AtomicGroups also support periodic boundary conditions via the periodicBox() method. If a box has been set, then isPeriodic() will return true. The periodic box is shared between the parent group and all derived groups. AtomicGroup copies have non-shared periodic boxes...
|
inline |
Creates a new AtomicGroup with n un-initialized atoms.
The atoms will all have ascending atomid's beginning with 1, but otherwise no other properties will be set.
GMatrix loos::AtomicGroup::alignOnto | ( | const AtomicGroup & | grp | ) |
Superimposes the current group onto the passed group.
Calls superposition to calculate the transformation matrix to superimpose the current group onto the passed one, then applies the transformation to the current group's coordinates.
|
inline |
Apply a functor or a function to each atom in the group.
apply() let's you apply a functor or a function pointer to each atom in the group. The functor is passed a pAtom. The functor object is also returned (in case it retained state). For example, the following code snippet shows how to calculate the centroid of a group using apply and a functor...
std::vector< uint > loos::AtomicGroup::atomOrderMapFrom | ( | const AtomicGroup & | g | ) |
Map the order of atoms in AtomicGroup g into the current group.
Note that the order is only checked within a residue. The residues must appear in the same order between the two groups. This addresses edge issues such as when psfgen reorders the atoms within a residue. The map is an index into the AtomicGroup g that puts g into the same order as the current group.
std::vector< GCoord > loos::AtomicGroup::boundingBox | ( | void | ) | const |
Bounding box for the group...
Returns a vector containing 2 GCoords, one containing (minx, miny, minz) and the other (maxx, maxy, maxz)
GCoord loos::AtomicGroup::centerAtOrigin | ( | void | ) |
Translates the group so that the centroid is at the origin.
Returns the old centroid of the group
AtomicGroup loos::AtomicGroup::centrifyByMolecule | ( | ) | const |
Replace a group with the center of masses of contained molecules.
The AtomicGroup is split into molecules. A new group is constructed where each atom is the center of mass of one molecule. The atom metadata is taken from the first atom of the associated molecule, but with the atom name "CEN".
|
virtual |
Creates a lightweight clone of this group (for polymorphism)
Despite the name, this is meant for polymorphic use. It is not a deep copy. If you don't understand what any of this means, then you almost certainly want to be using the copy() method instead.
Reimplemented in loos::Amber, loos::CHARMM, loos::MDTraj, loos::MMCIF, loos::PDB, loos::PSF, and loos::TinkerXYZ.
|
inline |
Returns true if any atom of current group is within dist angstroms of grp.
min is the minimum number of pair-wise contacts required to be considered in contact
|
inline |
Returns true if any atom of current group is within dist angstroms of grp.
min is the minimum number of pair-wise contacts required to be considered in contact
|
inline |
Determines if a pAtom is contained in this group using the EqualsOp atom-equality policy.
The problem with determining containment/intersection/etc is how to define when two atoms are equal... This is done by specifying a comparison functor (the EqualsOp) as a policy. There are two comparison policies currently in LOOS: AtomEquals and AtomCoordsEquals. The default behavior is to use AtomEquals which only compares a subset of the available Atom metadata. You can specify the more restrictive policy (or an user-defined policy) like:
Or as another example, comparing only residue numbers...
AtomicGroup loos::AtomicGroup::copy | ( | void | ) | const |
Creates a deep copy of this group.
This creates a non-polymorphic deep copy of an AtomicGroup. The additional catch is that it may end up involving extra data-movement as the copy is constructed and then copied back out to wherever you're putting it.
void loos::AtomicGroup::copyCoordinatesFrom | ( | const AtomicGroup & | g, |
const uint | offset = 0, | ||
const uint | length = 0 ) |
Copy coordinates from g into current group.
The offset is relative to the start of the current group and length is the number of coordinates to copy. If length + offset is greater than the number of atoms in the current group, then the excess coordinates will not be copied.
If length is 0, then all coordinates in g will be copied.
It is assumed that the atoms in g are in the appropriate order relative to the current group for the copy to make sense.
void loos::AtomicGroup::copyMappedCoordinatesFrom | ( | const AtomicGroup & | g | ) |
Copy the coordinates from the group mapping the atom order.
See AtomicGroup::atomOrderFrom(const AtomicGroup& g) for more information
If you know that the atoms are in the same order in both groups, then AtomicGroup::copyCoordinatesFrom() will be faster...
void loos::AtomicGroup::copyMappedCoordinatesFrom | ( | const AtomicGroup & | g, |
const std::vector< uint > & | order ) |
Given a mapping of atom order, copy the coordinates into the current group.
See AtomicGroup::atomOrderMapFrom(const AtomicGroup& g) for more information
If you know that the atoms are in the same order in both groups, then AtomicGroup::copyCoordinatesFrom() will be faster...
void loos::AtomicGroup::copyVelocitiesWithIndex | ( | const std::vector< GCoord > & | velocities | ) |
Copy velocities from a vector of GCoords using the atom index as an index into the vector.
This can be used to update a group's velocities if they come from a separate trajectory...
uint loos::AtomicGroup::deduceAtomicNumberFromMass | ( | const double | tol = 0.1 | ) |
Deduce atomic number from mass (if present), returning number of atoms assigned.
If an atom has a mass, then this is used to look up it's atomic number. Note that LOOS only has the first 96 elements in its tables. If a mass is not found in the LOOS table, then the atomic number is not modified (or set), otherwise any existing atomic number is overwritten.
std::vector< GCoord > loos::AtomicGroup::differenceVectors | ( | const AtomicGroup & | other | ) |
Compute difference vectors between two AtomicGroups.
Does not align the coordinates first
|
inline |
Distance-based search for bonds.
Searches for bonds within an AtomicGroup based on distance. does NOT clear the existing bond list prior to building new bonds. The default distance cutoff is 1.65. If a box (GCoord) is passed, then periodicity is taken into consideration.
pAtom loos::AtomicGroup::findById | ( | const int | id | ) | const |
Find a contained atom by its atomid.
The default behavior is to assume that the atoms in the AtomicGroup are not in order of increasing atomid and to therefore use a linear search. If the atoms are sorted (AtomicGroup::sort()), then the more efficient binary search will be used.
void loos::AtomicGroup::getCoords | ( | double ** | outseq, |
int * | m, | ||
int * | n ) |
This function is meant for Numpy/swig use. It will store the current model's coordinates into a newly allocated array (using malloc). The caller is expected to manage the memory.
AtomicGroup loos::AtomicGroup::getResidue | ( | pAtom | res | ) |
Given an Atom, return a group of all the atoms contained by its containing residue
Returns a vector of coordinates transformed by the passed XForm.
Does not alter the group's coordinates...
AtomicGroup loos::AtomicGroup::groupFromID | ( | const std::vector< int > & | id_list | ) | const |
Create a new group from a vector of atomids.
Note: when calling this, you'll want to make sure you use the outermost group (eg the psf or pdb you used to create things, rather than using a subselection, unless you're sure the subsection contains these atoms as well. The main use of this routine is to create a group of atoms bound to another atom.
Any missing atoms are ignored... This is in contrast with the previous behavior where missing atoms would throw an exception
double loos::AtomicGroup::hardContact | ( | const AtomicGroup & | group, |
double | radius, | ||
const GCoord & | box ) const |
Compute contact value of another AG with respect to the centroid of a given AG, using a hard step function S = 1; iff dist <= radius; else 0
double loos::AtomicGroup::hardContact2D | ( | const AtomicGroup & | group, |
double | radius, | ||
const GCoord & | box ) const |
Compute contact value of another AG with respect to the centroid of a given AG, using a hard step function S = 1; iff dist <= radius; else 0
|
inline |
Computes the intersection of two groups using the EqualsOp atom-equality policy.
See AtomicGroup::contains(const pAtom&, const EqualsOp&) for more details
greal loos::AtomicGroup::kineticEnergy | ( | ) |
Compute kinetic energy of group.
Assumes mass and velocity have been set. Output units are kcal/mol
double loos::AtomicGroup::logisticContact | ( | const AtomicGroup & | group, |
double | radius, | ||
int | sigma, | ||
const GCoord & | box ) const |
Compute the degree of contact between the centroid of this AG and the centroid of another AG, using a smooth logistic function S = 1/(1 + dist/radius)**sigma
double loos::AtomicGroup::logisticContact2D | ( | const AtomicGroup & | group, |
double | radius, | ||
int | sigma, | ||
const GCoord & | box ) const |
Compute the number of contacts between the centroid of this AG and the centroid of another AG, using a smooth logistic function S = 1/(1 + dist/radius)**sigma
|
inline |
Union of two groups using the specified atom-equality policy.
Note that the periodic box of the current group is unchanged by this operation
void loos::AtomicGroup::mergeImage | ( | ) |
Takes a group that's split across a periodic boundary and reimages it so it's all together, using the first atom in the AtomicGroup as the reference.
Does the same as the other mergeImage, only using the first atom in the AtomicGroup as the reference atom.
void loos::AtomicGroup::mergeImage | ( | pAtom & | p | ) |
Takes a group that's split across a periodic boundary and reimages it so it's all together.
Works by translating the system so one atom is in the center of the box, reimaging by atom (so now the group is all in the middle of the box), and then translating back.
If you don't want to give it a reference atom, call the version that takes no argument; it uses the first atom in the AtomicGroup.
std::vector< GCoord > loos::AtomicGroup::momentsOfInertia | ( | void | ) | const |
Computes the moments of inertia for a group.
Calculates the principal moments and principal axes (from the moment of inertia). This is distinct from the principalAxes() function which calculates the distribution of points about the centroid.
bool loos::AtomicGroup::operator== | ( | AtomicGroup & | rhs | ) |
Equality test for two groups.
The test for equality is based on whether or not the contained atom pointers are the same. This operator will also force both sides of the equation to be sorted.
bool loos::AtomicGroup::operator== | ( | const AtomicGroup & | rhs | ) | const |
Equality test for const groups.
Similar to the non-const version, but it will sort copies of the atom lists if they are not already sorted...
void loos::AtomicGroup::orientAlong | ( | const GCoord & | vec | ) |
Orient the principal axis of this group along the supplied vector.
The supplied vector does not need to be normalized.
double loos::AtomicGroup::packingScore | ( | const AtomicGroup & | other, |
const GCoord & | box, | ||
bool | norm = false ) const |
The packing score is the sum of 1/r^6 over all pairs of atoms, respecting periodicity. Quantity first defined in Grossfield, A., et al, Proc. Nat. Acad. Sci. USA, 2006, 103, 4888-4893
std::vector< GCoord > loos::AtomicGroup::principalAxes | ( | void | ) | const |
Compute the principal axes of a group.
Calculates the eigendecomposition of AA' where A is column-wise concatenation of coordinates from all atoms in the group. The mean coordinate is automatically subtracted from A... Returns a vector of GCoord's in order of decreasing magnitude of the corresponding eigenvalue. The eigenvalues are returned as a GCoord after the eigenvectors.
Example
Notes
double loos::AtomicGroup::principalAxesOrder | ( | void | ) | const |
Computes order parameter based on principalAxes.
Computes the molecular order parameter by taking 2nd and 3rd principal axes and transforming using a 2nd order Legendre polynomial. This produces a whole-chain quantity with values comparable to deuterium order parameters.
If the AtomicGroup is not planar, this returns the average of the order parameter using the 2nd and 3rd axes. If the molecule is planar, indicated by a near-zero 3rd moment, it simply returns the 2nd moment's order parameter.
This implementation assumes you've already done any needed transformations to the coordinates (eg reimaging by molecule)
void loos::AtomicGroup::pruneBonds | ( | ) |
Attempt to prune connectivity (only retain bonds to atoms within this AtomicGroup)
The connectivity list is searched for each atom and if a bond is not found in the current group, then it is removed from the bond-list for that atom. Note that this means that any AtomicGroup sharing the atom in question will also now have the modified bond list. It's therefore recommended that this function be called on a copy (AtomicGroup::copy()). Also note that FindById() does not implicitly sort the atoms for more efficient searching. You may want to call AtomicGroup::sort() prior to AtomicGroup::pruneBonds() if the exact atom order does not matter.
greal loos::AtomicGroup::radius | ( | const bool | use_atom_as_reference = false | ) | const |
Maximum radius from centroid of all atoms (not gyration)
If optional argument is true, uses coordinates of atom 0 instead of centroid. Argument is false by default.
void loos::AtomicGroup::reimage | ( | ) |
Translate the entire group so that the centroid is in the primary cell
void loos::AtomicGroup::resetAtomIndices | ( | ) |
Reset the atom indices (used for interfacing with trajectories)
The Atom index is the original ordering of atoms from whatever file format the model came from. This is used as an index into each frame of the trajectory for corresponding atom properties (such as coordinates). If an AtomicGroup is a subset, then it may be necessary to reset the indices when working with a subsetted trajectory as well. This function will reset the atom indices to be sequential, beginning with 0.
greal loos::AtomicGroup::rmsd | ( | const AtomicGroup & | v | ) |
Compute the RMSD between two groups.
Assumes a 1:1 correspondence between ith atoms. Does NOT transform the coordinates in any way.
std::vector< double > loos::AtomicGroup::scattering | ( | const double | qmin, |
const double | max, | ||
const uint | numValues, | ||
loos::FormFactorSet & | formFactors ) |
Computes X-ray scattering as a function of q, using I(q) = \sum_(atom pair) F_i(q) F_j(q) sin (q d_ij)/ (q d_ij)
This approximates scattering off of individual atoms. If you use this with explicit solvent, you will get truncation artifacts from the periodic box (although the code computes all distances using periodicity).
Form factors are from Szaloki, X-ray Spectrometry (1996), V25, 21-28
void loos::AtomicGroup::setCoords | ( | double * | seq, |
int | m, | ||
int | n ) |
This function is meant for Numpy/swig use in setting the model's coordinates. The passed array is row-major.
|
inline |
Is the array of atoms already sorted???
While we make some effort to ensure that alterations to the AtomicGroup will invalidate the sorted status, it's a good idea to explicitly sort if you want to make sure that the group is in fact sorted.
greal loos::AtomicGroup::sphericalVariance | ( | const pAtom | target | ) | const |
Spherical variance of group with respect to target atom.
spherical variance as a measure of how much atom "target" is inside this atomic group Mezei, J Mol Graph Modeling, 2003, 21, 463-472
|
inline |
Takes selection string as argument to be applied to each group after splitting. Returns a vector of AtomicGroups split based on bond connectivity;
std::vector< AtomicGroup > loos::AtomicGroup::splitByResidue | ( | void | ) | const |
Returns a vector of AtomicGroups, each comprising a single residue.
Splits an AtomicGroup into individual residues. The residue boundary is marked by either a change in the resid or in the segid.
std::vector< AtomicGroup > loos::AtomicGroup::splitByUniqueSegid | ( | void | ) | const |
Returns a vector of AtomicGroups split from the current group based on segid.
The groups that are returned will be in the same order that the segids appear in the source AtomicGroup
greal loos::AtomicGroup::stacking | ( | const AtomicGroup & | other, |
const GCoord & | box, | ||
const double | threshold = 5.0 ) const |
Estimate stacking, as between two nucleobases.
Algorithm: n1=normal to self; n2=normal to other, dx = difference between centroids stacking = (n1*n2)^2 *[(n1 + n2)/2 * dx]/|dx| * 1/1 + (dx/threshold)^6
AtomicGroup loos::AtomicGroup::subset | ( | const int | offset, |
const int | len = 0 ) |
GMatrix loos::AtomicGroup::superposition | ( | const AtomicGroup & | grp | ) |
Calculates the transformation matrix for superposition of groups.
Uses the Kabsch alignment method (via SVD) to calculate the transformation matrix that superimposes the current group onto the passed group. Returns the matrix.
If too few atoms are given for aligning, the correlation matrix may become singular and return fewer than three eigenpairs. If this is detected, superposition() will throw a NumericalError. The threshold for a zero-eigenvalue (really, a zero singular value) is set in AtomicGroup::superposition_zero_singular_value