diff --git a/package/CHANGELOG b/package/CHANGELOG index 133cfe3ea92..b53ffc987bf 100644 --- a/package/CHANGELOG +++ b/package/CHANGELOG @@ -21,6 +21,8 @@ The rules for this file: * 2.8.0 Fixes + * Fix `MDAnalysis.analysis.align.AlignTraj` not accepting writer kwargs + (Issue #4564, PR #4565) * Fix #4259 via removing argument `parallelizable` of `NoJump` transformation. * Fix doctest errors of analysis/pca.py related to rounding issues (Issue #3925, PR #4377) diff --git a/package/MDAnalysis/analysis/align.py b/package/MDAnalysis/analysis/align.py index dd0e0100e64..2bc80042b69 100644 --- a/package/MDAnalysis/analysis/align.py +++ b/package/MDAnalysis/analysis/align.py @@ -681,6 +681,7 @@ class AlignTraj(AnalysisBase): def __init__(self, mobile, reference, select='all', filename=None, prefix='rmsfit_', weights=None, tol_mass=0.1, match_atoms=True, strict=False, force=True, in_memory=False, + writer_kwargs=None, **kwargs): """Parameters ---------- @@ -720,6 +721,8 @@ def __init__(self, mobile, reference, select='all', filename=None, verbose : bool (optional) Set logger to show more information and show detailed progress of the calculation if set to ``True``; the default is ``False``. + writer_kwargs : dict (optional) + kwarg dict to be passed to the constructed writer Attributes @@ -780,6 +783,9 @@ def __init__(self, mobile, reference, select='all', filename=None, :attr:`rmsd` results are now stored in a :class:`MDAnalysis.analysis.base.Results` instance. + .. versionchanged:: 2.8.0 + Added ``writer_kwargs`` kwarg dict to pass to the writer + """ select = rms.process_selection(select) self.ref_atoms = reference.select_atoms(*select['reference']) @@ -816,10 +822,12 @@ def __init__(self, mobile, reference, select='all', filename=None, self.ref_atoms, self.mobile_atoms, tol_mass=tol_mass, strict=strict, match_atoms=match_atoms) + if writer_kwargs is None: + writer_kwargs = {} # with self.filename == None (in_memory), the NullWriter is chosen # (which just ignores input) and so only the in_memory trajectory is # retained - self._writer = mda.Writer(self.filename, natoms) + self._writer = mda.Writer(self.filename, natoms, **writer_kwargs) self._weights = get_weights(self.ref_atoms, weights) @@ -1661,4 +1669,4 @@ def get_atoms_byres(g, match_mask=np.logical_not(mismatch_mask)): logger.error(errmsg) raise SelectionError(errmsg) - return ag1, ag2 + return ag1, ag2 \ No newline at end of file diff --git a/testsuite/MDAnalysisTests/analysis/test_align.py b/testsuite/MDAnalysisTests/analysis/test_align.py index 78ba028a114..63dfd25db80 100644 --- a/testsuite/MDAnalysisTests/analysis/test_align.py +++ b/testsuite/MDAnalysisTests/analysis/test_align.py @@ -353,6 +353,17 @@ def test_AlignTraj_in_memory(self, universe, reference, tmpdir): self._assert_rmsd(reference, universe, 0, 6.929083044751061) self._assert_rmsd(reference, universe, -1, 0.0) + def test_AlignTraj_writer_kwargs(self, universe, reference, tmpdir): + # Issue 4564 + writer_kwargs = dict(precision=2) + with tmpdir.as_cwd(): + aligner = align.AlignTraj(universe, reference, + select='protein and name CA', + filename='aligned_traj.xtc', + writer_kwargs=writer_kwargs, + in_memory=False).run() + assert_equal(aligner._writer.precision, 2) + def _assert_rmsd(self, reference, fitted, frame, desired, weights=None): fitted.trajectory[frame] rmsd = rms.rmsd(reference.atoms.positions, fitted.atoms.positions,