com¶
-
mdtools.structure.
com
(ag, pbc=False, compound='group', make_whole=False, debug=False)[source]¶ Calculate the center of mass of (compounds of) a MDAnalysis
AtomGroup
.- Parameters
ag (MDAnalysis.core.groups.AtomGroup) – The MDAnalysis
AtomGroup
for which to calculate the center of mass.pbc (bool, optional) – If
True
and compound is'group'
, move allAtoms
to the primary unit cell before calculating the center of mass. IfTrue
and compound is not'group'
, the center of mass of each compound will be calculated without moving anyAtoms
to keep the compounds intact (if they were intact before). Instead, the resulting position vectors will be moved to the primary unit cell after calculating the center of mass.compound ({‘group’, ‘segments’, ‘residues’, ‘molecules’, ‘fragments’}, optional) – The compounds of ag for which to calculate the center of mass. If
'group'
, the center of mass of allAtoms
in ag will be returned as a single position vector. Else, the centers of mass of eachSegment
,Residue
, molecule, orfragment
contained in ag will be returned as an array of position vectors, i.e. a 2d array. Refer to the MDAnalysis’ user guide for an explanation of these terms. Note that in any case, even if compound is e.g.'residues'
, only theAtoms
belonging to ag are taken into account, even if the compound might comprise additionalAtoms
that are not contained in ag.make_whole (bool, optional) – If
True
, first allAtoms
in ag are wrapped back into the primary unit cell, to make sure that the algorithm is working properly. Then compounds whose bonds are broken across the box edges are made whole again. This means that making compounds whole in an unwrapped trajectory while keeping the trajectory unwrapped is not possible with this function.debug (bool, optional) – If
True
, run in debug mode.
- Returns
center – Center of mass positions of the compounds in ag. If compound was set to
'group'
, the output will be a single position vector. Else, the output will be a 2d array of shape(n, 3)
wheren
is the number of compounds in ag.- Return type
See also
MDAnalysis.core.groups.AtomGroup.center_of_mass()
Center of mass of (compounds of) the group
MDAnalysis.core.groups.AtomGroup.center()
Weighted center of (compounds of) the group
Notes
This function uses the
center_of_mass()
method of the inputAtomGroup
to calculate the center of mass.If make_whole is
True
, allAtoms
in ag are wrapped back into the primary unit cell usingmdtools.box.wrap()
before callingcenter_of_mass()
with the option unwrap set toTrue
. This is done to make sure that the unwrap algorithm (better called “make whole” algorithm) ofcenter_of_mass()
is working properly. This means that making compounds whole in an unwrapped trajectory while keeping the trajectory unwrapped is not possible with this function.Todo
Check if it is really necessary to wrap all
Atoms
back into the primary unit cell before callingcenter_of_mass()
with unwrap set toTrue
. The currently done back-wrapping is a serious problem, because it implies an inplace change of theAtom
coordinates.If the no wrapping is necessary, we can mark this function as deprecated and instead recommend the direct use of
center_of_mass()
.