center_of_mass¶
- mdcraft.algorithm.molecule.center_of_mass(group: AtomGroup = None, grouping: str = None, *, masses: ndarray[float] | list[ndarray[float]] = None, positions: ndarray[float] | list[ndarray[float]] = None, images: ndarray[int] | list[ndarray[int]] = None, dimensions: ndarray[float] = None, n_groups: int = None, raw: bool = False) ndarray[float] | tuple[ndarray[float], ndarray[float], ndarray[float]][source]¶
Computes the centers of mass \(\mathbf{R}_\mathrm{com}\) for a collection of atoms.
For a group of \(N\) atoms with masses \(m_i\) and positions \(\mathbf{r}_i\), the center of mass is defined as
\[\mathbf{R}_\mathrm{com}=\frac{\sum_{i=1}^Nm_i \mathbf{r}_i}{\sum_{i=1}^Nm_i}\]Note
This function supports a wide variety of inputs, depending on how the atom information is provided and what should be calculated.
When an
MDAnalysis.core.groups.AtomGroupobject is provided in group, the atom masses and positions are retrieved from it and do not need to be provided in masses and positions, respectively. If theAtomGroupabides by the standard topological heirarchy, you can specify the desired grouping in grouping and the appropriate centers of mass will be calculated. Otherwise, if and only if theAtomGroupcontains equisized or identical groups corresponding to the desired grouping (i.e., theAtomGrouphas atoms that are or can be treated as nonbonded entities or topological groups with the same number of but not necessarily identical constituents), you can provide the total number of groups in n_groups and the atom masses and positions will be distributed accordingly.If the trajectory is not unwrapped, the number of periodic boundary crossings (and optionally, the system dimensions if they are not embedded in the
AtomGroup) can be provided in images (and dimensions).If the
AtomGroupdoes not have the correct structural information and the residues or segments do not contain the same number of atoms, the atom masses and positions can each be provided directly as anumpy.ndarrayor list in masses and positions, respectively. To calculate the overall center of mass, the array-like object holding the masses should be one-dimensional, while that containing the positions should be two-dimensional. To calculate centers of mass for multiple groups, the array-like object holding the masses should be two-dimensional (indices: group, atom), while that containing the positions should be three-dimensional (indices: group, atom, axis). When a list is used, the inner arrays do not have to be homogeneously shaped, thus allowing you to calculate the centers of mass for residues or segments with different numbers of atoms.You may also provide only one of the atom masses or positions, in which case the missing information will be retrieved from the
AtomGroupprovided in group. This is generally not recommended since the shapes of the provided and retrieved arrays may be incompatible.- Parameters:
- groupMDAnalysis.AtomGroup, optional
Collection of atoms to compute the centers of mass for. If not provided, the atom masses and posititions must be provided directly in masses and positions, respectively.
- groupingstr, optional
Determines which center of mass is calculated if atom masses and positions are retrieved from group.
Valid values:
None: Center of mass of all atoms in group."residues": Centers of mass for each residue or molecule in group."segments": Centers of mass for each segment or chain in group.
- massesarray-like, keyword-only, optional
Atom masses.
Shape:
The general ungrouped shape is \((N,)\).
For equisized or identical groups, the
numpy.ndarrayobject should have shape\((N,)\) for the overall center of mass,
\((N_\mathrm{res},\,N/N_\mathrm{res})\) for the residue centers of mass, where \(N_\mathrm{res}\) is the number of residues, or
\((N_\mathrm{seg},\,N/N_\mathrm{seg}\) for the segment centers of mass, where \(N_\mathrm{seg}\) is the number of segments.
For groups with different numbers of atoms, the list should contain inner array-like objects holding the masses of the atoms in each group.
Reference unit: \(\mathrm{g/mol}\).
- positionsarray-like, keyword-only, optional
Atom positions.
Shape:
The general ungrouped shape is \((N,\,3)\).
For equisized or identical groups, the
numpy.ndarrayobject should have shape\((N,\,3)\) for the overall center of mass,
\((N_\mathrm{res},\,N/N_\mathrm{res},\,3)\) for the residue centers of mass, or
\((N_\mathrm{seg},\,N/N_\mathrm{seg},\,3)\) for the segment centers of mass.
For groups with different numbers of atoms, the list should contain inner array-like objects holding the positions of the atoms in each group.
Reference unit: \(\mathrm{Å}\).
- imagesarray-like, keyword-only, optional
Image flags for the atoms. Must be provided to get correct results if the trajectory is wrapped.
Shape: Same as positions.
- dimensionsnumpy.ndarray, keyword-only, optional
System dimensions. Must be provided if images is provided and group is not provided or does not contain the system dimensions.
Shape: \((3,)\).
Reference unit: \(\mathrm{Å}\).
- n_groupsint, keyword-only, optional
Number of residues or segments. Must be provided if group has an irregular topological heirarchy or the masses and positions arrays have the general ungrouped shapes.
- rawbool, keyword-only, default:
False Determines whether atom masses and positions are returned.
- Returns:
- comnumpy.ndarray
Centers of mass.
Shape:
\((3,)\) for
grouping=None.\((N_\mathrm{res},\,3)\) for
grouping="residues".\((N_\mathrm{seg},\,3)\) for
grouping="segments".
- massesnumpy.ndarray
Atom masses. Only returned if group was provided and contains equisized or identical groups, and
raw=True.Shape:
\((N,)\) for
grouping=None.\((N_\mathrm{res},\,N/N_\mathrm{res})\) for
grouping="residues".\((N_\mathrm{seg},\,N/N_\mathrm{seg})\) for
grouping="segments".
Reference unit: \(\mathrm{g/mol}\).
- positionsnumpy.ndarray
Unwrapped atom positions. Only returned if group was provided and contains equisized or identical groups, and
raw=True.Shape:
\((N,\,3)\) for
grouping=None.\((N_\mathrm{res},\,N/N_\mathrm{res},\,3)\) for
grouping="residues".\((N_\mathrm{seg},\,N/N_\mathrm{seg},\,3)\) for
grouping="segments".
Reference unit: \(\mathrm{Å}\).
Examples
For an
MDAnalysis.core.groups.AtomGroupobject with all necessary topological information and an unwrapped trajectory, the overall, per-residue, and per-segment centers of mass can be calculated as follows:>>> universe = mda.Universe("topology.pdb", "trajectory.nc") >>> group = universe.select_atoms("all") >>> com_ovr = center_of_mass(group) >>> com_res = center_of_mass(group, "residues") >>> com_seg = center_of_mass(group, "segments")
If the
AtomGroupdoes not contain residue or segment information and the per-residue or per-segment centers of mass, respectively, are desired, the number of residues or segments can be provided in n_groups:>>> com_res = center_of_mass(group, "residues", n_groups=2) # 2 residues
If the trajectory is wrapped, the number of periodic boundary crossings must be provided in images. Additionally, if the system dimensions are not embedded in the
AtomGroup, they must also be provided in dimensions:>>> images = np.array(((1, 0, 0), (0, 0, 0), (0, -1, 0), ... (0, 1, 0), (0, 0, 1), (-1, 1, 0))) >>> dimensions = np.array((5.0, 8.0, 10.0)) >>> com_ovr = center_of_mass(group, images=images, dimensions=dimensions)
Alternatively, if the atom masses and positions are directly available as
numpy.ndarrayobjects, the overall center of mass can be calculated as follows:>>> masses = np.array((12.01, 1.01, 1.01, 12.01, 1.01, 1.01)) >>> positions = np.array(((0.0, -0.07579, 0.0), ... (0.86681, 0.60144, 0.0), ... (-0.86681, 0.60144, 0.0), ... (0.0, -0.07579, 1.0), ... (0.86681, 0.60144, 1.0), ... (-0.86681, 0.60144, 1.0))) >>> com_ovr = center_of_mass(masses=masses, positions=positions)
If the per-residue center of mass is desired, the number of residues can be provided in n_groups:
>>> com_res = center_of_mass(masses=masses, positions=positions, n_groups=2)
or the arrays containing the atom masses and positions can be reshaped to the appropriate shapes:
>>> n_groups = 2 >>> com_res = center_of_mass(masses=masses.reshape((n_groups, -1)), ... positions=positions.reshape((n_groups, -1, 3)))
Like before, if the trajectory is wrapped, the number of periodic boundary crossings and system dimensions must be provided in images and dimensions, respectively:
>>> images = np.array(((0, 0, 0), (0, 0, 0), (0, 0, 0), ... (0, 0, 1), (0, 0, 1), (0, 0, 1))) >>> dimensions = np.array((12.0, 12.0, 12.0)) >>> com_ovr = center_of_mass(masses=masses, positions=positions, ... images=images, dimensions=dimensions)
Finally, if the per-residue or per-segment center of mass is desired but the groups contain different numbers of atoms, the the atom masses and positions can be provided as lists of arrays holding the masses and positions of the atoms in each group:
>>> masses = [(12.01, 1.01, 1.01), (22.99,), (12.01, 1.01)] >>> positions = [((0.0, -0.07579, 0.0), ... (0.86681, 0.60144, 0.0), ... (-0.86681, 0.60144, 0.0)), ... ((0.0, 0.0, 0.0),), ... ((0.0, -0.07579, 1.0), ... (0.86681, 0.60144, 1.0))] >>> com_res = center_of_mass(masses=masses, positions=positions)
It is still possible to pass in the number of periodic boundary crossings if the trajectory is wrapped, but attention must be paid to the shape of the array:
>>> images = [((0, 0, 0), (0, 0, 0), (0, 0, 0)), ... ((1, 1, 1),), ... ((0, 1, 0), (0, 1, 0))] >>> dimensions = np.array((10.0, 10.0, 10.0)) >>> com_res = center_of_mass(masses=masses, positions=positions, ... images=images, dimensions=dimensions)