cmp_contact_count_matrix
- mdtools.structure.cmp_contact_count_matrix(cm, natms_per_refcmp=1, natms_per_selcmp=1, dtype=int)[source]
Take an
Atomcontact matrix and sum the contacts of allAtomsbelonging to the same compound.A compound is usually a chemically meaningful subgroup of an
AtomGroup. This can e.g. be aSegment,Residue,fragmentor a singleAtom. Refer to the MDAnalysis’ user guide for an explanation of these terms. Note that in any case, onlyAtomsbelonging to the inputAtomGroupare taken into account, even if the compound might comprise additionalAtomsthat are not contained in the inputAtomGroup.- Parameters:
cm (
array_like) – (Boolean) contact matrix of shape(m, n)as e.g. generated bymdtools.structure.contact_matrix(), wheremis the number of referenceAtomsandnis the number of selectionAtoms.natms_per_refcmp (
intorarray_like, optional) – Number ofAtomsper reference compound. Can be a single integer or an array of integers. If natms_per_refcmp is a single integer, all reference compounds are assumed to contain the same number ofAtoms. In this case, natms_per_refcmp must be an integer divisor ofcm.shape[0]. If natms_per_refcmp is an array of integers, it must contain the number of referenceAtomsfor each single reference compound. In this case,sum(natms_per_refcmp)must be equal tocm.shape[0].natms_per_selcmp (
intorarray_like, optional) – Same for selection compounds (natms_per_selcmp is checked againstcm.shape[1]).dtype (
dtype, optional) – Data type of the output array.
- Returns:
cmp_ccm (
numpy.ndarray) – A new (smaller) contact count matrix indicating how many contacts exist between theAtomsof a given reference and selection compound. If both natms_per_refcmp and natms_per_selcmp are1, the input contact matrix is returned instead with its data type converted to dtype.
See also
mdtools.structure.contact_matrix()Construct a boolean contact matrix for two MDAnalysis
AtomGroupsmdtools.structure.natms_per_cmp()Get the number of
Atomsof each compound in an MDAnalysisAtomGroupmdtools.structure.cmp_contact_matrix()Convert an
Atomcontact matrix to a compound contact matrixmdtools.structure.contact_hists()Bin the number of contacts between reference and selection compounds into histograms.
numpy.ufunc.reduceat()Perform a (local)
reduce()with specified slices over a single axis
Notes
Atomsbelonging to the same compound must form a contiguous set in the input contact matrix, otherwise the result will be wrong.If both natms_per_refcmp and natms_per_selcmp are
1,cmp_contact_count_matrix(cm, dtype=dtype)is equivalent tonp.asarray(cm, dtype=dtype).Examples
>>> cm = np.array([[ True, True, False, False], ... [ True, False, False, True], ... [False, False, True, True]]) >>> mdt.strc.cmp_contact_count_matrix(cm, natms_per_refcmp=3) array([[2, 1, 1, 2]]) >>> mdt.strc.cmp_contact_count_matrix(cm, natms_per_selcmp=2) array([[2, 0], [1, 1], [0, 2]]) >>> mdt.strc.cmp_contact_count_matrix( ... cm, natms_per_refcmp=[1, 2], natms_per_selcmp=[2, 2] ... ) array([[2, 0], [1, 3]]) >>> mdt.strc.cmp_contact_count_matrix( ... cm, ... natms_per_refcmp=[1, 2], ... natms_per_selcmp=[2, 2], ... dtype=np.uint32, ... ) array([[2, 0], [1, 3]], dtype=uint32)
If both natms_per_refcmp and natms_per_selcmp are
1,cmp_contact_count_matrix(cm, dtype=dtype)is equivalent tonp.asarray(cm, dtype=dtype):>>> mdt.strc.cmp_contact_count_matrix(cm, dtype=bool) is cm True >>> mdt.strc.cmp_contact_count_matrix(cm, dtype=int) is cm False
Edge cases:
>>> cm = np.array([], dtype=bool).reshape(0, 4) >>> mdt.strc.cmp_contact_count_matrix(cm, natms_per_refcmp=[]) array([], shape=(0, 4), dtype=int64) >>> mdt.strc.cmp_contact_count_matrix(cm, natms_per_refcmp=1) array([], shape=(0, 4), dtype=int64)
>>> cm = np.array([], dtype=bool).reshape(6, 0) >>> mdt.strc.cmp_contact_count_matrix( ... cm, natms_per_refcmp=3, natms_per_selcmp=[], dtype=np.uint32 ... ) array([], shape=(2, 0), dtype=uint32)