diff --git a/package/CHANGELOG b/package/CHANGELOG index c6d2f3f4b3b..b87116b2940 100644 --- a/package/CHANGELOG +++ b/package/CHANGELOG @@ -19,6 +19,7 @@ The rules for this file: * 2.8.0 Fixes + * Fix bug in PCA preventing use of `frames=...` syntax (PR #4423) * Fix `analysis/diffusionmap.py` iteration through trajectory to iteration over `self._sliced_trajectory`, hence supporting `DistanceMatrix.run(frames=...)` (PR #4433) diff --git a/package/MDAnalysis/analysis/pca.py b/package/MDAnalysis/analysis/pca.py index cad6946b2dc..4ed5486e4c6 100644 --- a/package/MDAnalysis/analysis/pca.py +++ b/package/MDAnalysis/analysis/pca.py @@ -142,7 +142,8 @@ class PCA(AnalysisBase): generates the principal components of the backbone of the atomgroup and then transforms those atomgroup coordinates by the direction of those variances. Please refer to the :ref:`PCA-tutorial` for more detailed - instructions. + instructions. When using mean selections, the first frame of the selected + trajectory slice is used as a reference. Parameters ---------- @@ -230,6 +231,12 @@ class PCA(AnalysisBase): ``mean`` input now accepts coordinate arrays instead of atomgroup. :attr:`p_components`, :attr:`variance` and :attr:`cumulated_variance` are now stored in a :class:`MDAnalysis.analysis.base.Results` instance. + .. versionchanged:: 2.8.0 + ``self.run()`` can now appropriately use ``frames`` parameter (bug + described by #4425 and fixed by #4423). Previously, behaviour was to + manually iterate through ``self._trajectory``, which would + incorrectly handle cases where the ``frame`` argument + was passed. """ def __init__(self, universe, select='all', align=False, mean=None, @@ -247,7 +254,7 @@ def __init__(self, universe, select='all', align=False, mean=None, def _prepare(self): # access start index - self._u.trajectory[self.start] + self._sliced_trajectory[0] # reference will be start index self._reference = self._u.select_atoms(self._select) self._atoms = self._u.select_atoms(self._select) @@ -275,7 +282,7 @@ def _prepare(self): self._ref_atom_positions -= self._ref_cog if self._calc_mean: - for ts in ProgressBar(self._u.trajectory[self.start:self.stop:self.step], + for ts in ProgressBar(self._sliced_trajectory, verbose=self._verbose, desc="Mean Calculation"): if self.align: mobile_cog = self._atoms.center_of_geometry() diff --git a/testsuite/MDAnalysisTests/analysis/test_pca.py b/testsuite/MDAnalysisTests/analysis/test_pca.py index e157bee994e..c72971c30ed 100644 --- a/testsuite/MDAnalysisTests/analysis/test_pca.py +++ b/testsuite/MDAnalysisTests/analysis/test_pca.py @@ -133,6 +133,12 @@ def test_no_frames(u): PCA(u, select=SELECTION).run(stop=1) +def test_can_run_frames(u): + atoms = u.select_atoms(SELECTION) + u.transfer_to_memory() + PCA(u, select=SELECTION, mean=None).run(frames=[0, 1]) + + def test_transform(pca, u): ag = u.select_atoms(SELECTION) pca_space = pca.transform(ag, n_components=1)