natms_per_cmp
- mdtools.structure.natms_per_cmp(ag, cmp, return_array=False, return_cmp_ix=False, check_contiguous=False)[source]
Get the number of
Atoms
of each compound in an MDAnalysisAtomGroup
.- Parameters:
ag (
MDAnalysis.core.groups.AtomGroup
) – The MDAnalysisAtomGroup
for which to get the number ofAtoms
per compound.cmp (
{'group', 'segments', 'residues', 'molecules', 'fragments', 'atoms'}
) – The compounds of ag for which to get the number ofAtoms
. If'atoms'
, the output will simply be1
or an array of ones (depending on the value of return_array). Else, the number ofAtoms
in the entire group or of eachSegment
,Residue
, molecule orfragment
in ag is returned. Refer to the MDAnalysis’ user guide for an explanation of these terms. Note that in any case, even if cmp 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.return_array (
bool
, optional) – IfTrue
, natms_per_cmp (the output of this function) will always be an array with a length equal to the number of compounds in ag. IfFalse
, natms_per_cmp will be an integer if all compounds in ag have the same number ofAtoms
. However, if this is not the case, natms_per_cmp will be an array even if return_array isFalse
.return_cmp_ix (
bool
, optional) – IfTrue
, additionally return the unique indices of the compounds as assigned by MDAnalysis. If cmp is e.g.'residues'
, this isnp.unique(ag.resindices)
.check_contiguous (
bool
, optional) – IfTrue
, check ifAtoms
belonging to the same compound form a contiguous set in the inputAtomGroup
. This is e.g. important when constructing compound contact matrices fromAtom
contact matrices withmdtools.structure.cmp_contact_matrix()
.
- Returns:
natms_per_cmp (
int
ornumpy.ndarray
) – Number ofAtoms
of each compound contained in ag. If return_array isFalse
and all compounds in ag have the same number ofAtoms
, a single integer is returned. If return_array isTrue
or if not all compounds in ag have the same number ofAtoms
, an array is returned containing the number ofAtoms
of each single compound in ag.cmp_ix (
array_like
) – Unique compound indices as assigned by MDAnalysis. Only returned if return_cmp_ix isTrue
.
See also
mdtools.structure.cmp_attr()
Get attributes of an MDAnalysis
AtomGroup
compound-wise.mdtools.structure.cmp_contact_matrix()
Convert an
Atom
contact matrix to a compound contact matrixmdtools.structure.cmp_contact_count_matrix()
Take an
Atom
contact matrix and sum the contacts of allAtoms
belonging to the same compoundmdtools.structure.contact_hists()
Bin the number of contacts between reference and selection compounds into histograms
Examples
Create an MDAnalysis Universe from scratch for the following examples.
>>> import MDAnalysis as mda >>> # Number of segments. >>> n_seg = 1 >>> # Number of residues per segment. >>> n_res = [n+2 for n in range(n_seg)] >>> # Number of atoms per residue. >>> n_atms = [n+2 for n_r in n_res for n in range(n_r)] >>> u = mda.Universe.empty( ... n_atoms=sum(n_atms), ... n_residues=sum(n_res), ... n_segments=n_seg, ... atom_resindex=[ ... ix for ix, n_a in enumerate(n_atms) for _ in range(n_a) ... ], ... residue_segindex=[ ... ix for ix, n_r in enumerate(n_res) for _ in range(n_r) ... ], ... trajectory=True, ... velocities=True, ... forces=True, ... ) >>> bonds = [ ... (sum(n_atms[:i]), j) ... for i in range(len(n_atms)) ... for j in range(sum(n_atms[:i])+1, sum(n_atms[:i+1])) ... ] >>> u.add_TopologyAttr("bonds", bonds) >>> ag = u.atoms
>>> mdt.strc.natms_per_cmp(ag, cmp="group") 5 >>> mdt.strc.natms_per_cmp(ag, cmp="group", return_array=True) array([5]) >>> mdt.strc.natms_per_cmp(ag, cmp="segments") 5 >>> mdt.strc.natms_per_cmp(ag, cmp="residues") array([2, 3]) >>> mdt.strc.natms_per_cmp(ag, cmp="residues", return_cmp_ix=True) (array([2, 3]), array([0, 1])) >>> mdt.strc.natms_per_cmp(ag, cmp="fragments") array([2, 3]) >>> mdt.strc.natms_per_cmp(ag, cmp="atoms") 1