Skip to content

Commit 2e2ad9f

Browse files
authored
Fix weights deprecations (#1356)
* deprecate mass_weighted in psa * update CHANGELOG
1 parent e12a647 commit 2e2ad9f

File tree

2 files changed

+62
-31
lines changed

2 files changed

+62
-31
lines changed

package/CHANGELOG

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -51,6 +51,8 @@ Fixes
5151
Changes
5252
* Enable various pylint warnings to increase python 3 compatibility
5353
* Change Mathjax cdn (Issue #1313)
54+
* Change mass_weight to weights for PSA analysis
55+
* Move mass_weights deprecation to version 0.18
5456
* Docs moved to http://docs.mdanalysis.org (Issue #1315)
5557
and made responsive (Alabaster theme with readable-sphinx CSS)
5658
(Issue #378)

package/MDAnalysis/analysis/psa.py

Lines changed: 60 additions & 31 deletions
Original file line numberDiff line numberDiff line change
@@ -804,7 +804,7 @@ def __init__(self, universe, reference, ref_select='name CA',
804804

805805
def fit_to_reference(self, filename=None, prefix='', postfix='_fit',
806806
rmsdfile=None, targetdir=os.path.curdir,
807-
mass_weighted=False, tol_mass=0.1):
807+
mass_weighted=None, weights=None, tol_mass=0.1):
808808
"""Align each trajectory frame to the reference structure
809809
810810
Parameters
@@ -817,8 +817,10 @@ def fit_to_reference(self, filename=None, prefix='', postfix='_fit',
817817
prefix for auto-generating the new output filename
818818
rmsdfile : str (optional)
819819
file name for writing the RMSD time series [``None``]
820-
mass_weighted : bool (optional)
821-
do a mass-weighted RMSD fit, default is ``False``
820+
mass_weighted : bool (deprecated)
821+
do a mass-weighted RMSD fit
822+
weights : str/array_like (optional)
823+
choose weights. If 'str' uses masses as weights
822824
tol_mass : float (optional)
823825
Reject match if the atomic masses for matched atoms differ by more
824826
than `tol_mass` [0.1]
@@ -831,7 +833,18 @@ def fit_to_reference(self, filename=None, prefix='', postfix='_fit',
831833
Notes
832834
-----
833835
Uses :class:`MDAnalysis.analysis.align.AlignTraj` for the fitting.
836+
837+
838+
.. deprecated:: 0.16.1
839+
Instead of ``mass_weighted=True`` use new ``weights='mass'``;
840+
refactored to fit with AnalysisBase API
834841
"""
842+
if mass_weighted is not None:
843+
warnings.warn("mass weighted is deprecated argument. Please use "
844+
" 'weights=\"mass\" instead. Will be removed in 0.17.0",
845+
category=DeprecationWarning)
846+
if mass_weighted:
847+
weights = 'mass'
835848
head, tail = os.path.split(self.trj_name)
836849
oldname, ext = os.path.splitext(tail)
837850
filename = filename or oldname
@@ -842,7 +855,7 @@ def fit_to_reference(self, filename=None, prefix='', postfix='_fit',
842855
select=self.ref_select,
843856
filename=self.newtrj_name,
844857
prefix=prefix,
845-
mass_weighted=mass_weighted,
858+
weights=weights,
846859
tol_mass=tol_mass).run()
847860
if rmsdfile is not None:
848861
aligntrj.save(rmsdfile)
@@ -903,7 +916,7 @@ def to_path(self, fitted=False, select=None, flat=False):
903916

904917

905918
def run(self, align=False, filename=None, postfix='_fit', rmsdfile=None,
906-
targetdir=os.path.curdir, mass_weighted=False, tol_mass=0.1,
919+
targetdir=os.path.curdir, mass_weighted=None, weights=None, tol_mass=0.1,
907920
flat=False):
908921
r"""Generate a path from a trajectory and reference structure.
909922
@@ -937,8 +950,10 @@ def run(self, align=False, filename=None, postfix='_fit', rmsdfile=None,
937950
prefix for auto-generating the new output filename
938951
rmsdfile : str (optional)
939952
file name for writing the RMSD time series [``None``]
940-
mass_weighted : bool (optional)
941-
do a mass-weighted RMSD fit
953+
mass_weighted : bool (deprecated)
954+
do a mass-weighted RMSD fit
955+
weights : str/array_like (optional)
956+
choose weights. If 'str' uses masses as weights
942957
tol_mass : float (optional)
943958
Reject match if the atomic masses for matched atoms differ by more
944959
than *tol_mass* [0.1]
@@ -951,12 +966,23 @@ def run(self, align=False, filename=None, postfix='_fit', rmsdfile=None,
951966
-------
952967
topology_trajectory : tuple
953968
A tuple of the topology name and new trajectory name.
969+
970+
971+
.. deprecated:: 0.16.1
972+
Instead of ``mass_weighted=True`` use new ``weights='mass'``;
973+
refactored to fit with AnalysisBase API
954974
"""
975+
if mass_weighted is not None:
976+
warnings.warn("mass weighted is deprecated argument. Please use "
977+
" 'weights=\"mass\" instead. Will be removed in 0.17.0",
978+
category=DeprecationWarning)
979+
if mass_weighted:
980+
weights = 'mass'
955981
if align:
956-
self.u_fitted = self.fit_to_reference( \
957-
filename=filename, postfix=postfix, \
958-
rmsdfile=rmsdfile, targetdir=targetdir, \
959-
mass_weighted=False, tol_mass=0.1)
982+
self.u_fitted = self.fit_to_reference(
983+
filename=filename, postfix=postfix,
984+
rmsdfile=rmsdfile, targetdir=targetdir,
985+
weights=weights, tol_mass=0.1)
960986
self.path = self.to_path(fitted=align, flat=flat)
961987
return self.top_name, self.newtrj_name
962988

@@ -1263,8 +1289,6 @@ def __init__(self, universes, reference=None, ref_select='name CA',
12631289
can also each be a list of selection strings (to generate an
12641290
AtomGroup with defined atom order as described under
12651291
:ref:`ordered-selections-label`).
1266-
mass_weighted : bool
1267-
do a mass-weighted RMSD fit [``False``]
12681292
tol_mass : float
12691293
Reject match if the atomic masses for matched atoms differ by more
12701294
than *tol_mass* [0.1]
@@ -1350,7 +1374,8 @@ def __init__(self, universes, reference=None, ref_select='name CA',
13501374
self._psa_pairs = None # (distance vector order) list of all PSAPairs
13511375

13521376

1353-
def generate_paths(self, **kwargs):
1377+
def generate_paths(self, align=False, filename='fitted', infix='', mass_weighted=None, weights=None,
1378+
tol_mass=False, ref_frame=None, flat=False, save=True, store=True):
13541379
"""Generate paths, aligning each to reference structure if necessary.
13551380
13561381
Parameters
@@ -1364,8 +1389,10 @@ def generate_paths(self, **kwargs):
13641389
infix : str
13651390
additional tag string that is inserted into the output filename of
13661391
the fitted trajectory files ['']
1367-
mass_weighted : bool
1368-
do a mass-weighted RMSD fit
1392+
mass_weighted : bool (deprecated)
1393+
do a mass-weighted RMSD fit
1394+
weights : str/array_like (optional)
1395+
choose weights. If 'str' uses masses as weights
13691396
tol_mass : float
13701397
Reject match if the atomic masses for matched atoms differ by more
13711398
than *tol_mass*
@@ -1393,28 +1420,31 @@ def generate_paths(self, **kwargs):
13931420
persistence and can be accessed as the attribute
13941421
:attr:`PSAnalysis.paths`.
13951422
1423+
1424+
.. deprecated:: 0.16.1
1425+
Instead of ``mass_weighted=True`` use new ``weights='mass'``;
1426+
refactored to fit with AnalysisBase API
13961427
"""
1397-
align = kwargs.pop('align', False)
1398-
filename = kwargs.pop('filename', 'fitted')
1399-
infix = kwargs.pop('infix', '')
1400-
mass_weighted = kwargs.pop('mass_weighted', False)
1401-
tol_mass = kwargs.pop('tol_mass', False)
1402-
ref_frame = kwargs.pop('ref_frame', self.ref_frame)
1403-
flat = kwargs.pop('flat', False)
1404-
save = kwargs.pop('save', True)
1405-
store = kwargs.pop('store', False)
1428+
if mass_weighted is not None:
1429+
warnings.warn("mass weighted is deprecated argument. Please use "
1430+
" 'weights=\"mass\" instead. Will be removed in 0.17.0",
1431+
category=DeprecationWarning)
1432+
if mass_weighted:
1433+
weights = 'mass'
1434+
if ref_frame is None:
1435+
ref_frame = self.ref_frame
14061436

14071437
paths = []
14081438
fit_trj_names = []
14091439
for i, u in enumerate(self.universes):
1410-
p = Path(u, self.u_reference, ref_select=self.ref_select, \
1440+
p = Path(u, self.u_reference, ref_select=self.ref_select,
14111441
path_select=self.path_select, ref_frame=ref_frame)
14121442
trj_dir = self.targetdir + self.datadirs['fitted_trajs']
14131443
postfix = '{0}{1}{2:03n}'.format(infix, '_psa', i+1)
1414-
top_name, fit_trj_name = p.run(align=align, filename=filename, \
1415-
postfix=postfix, \
1416-
targetdir=trj_dir, \
1417-
mass_weighted=mass_weighted, \
1444+
top_name, fit_trj_name = p.run(align=align, filename=filename,
1445+
postfix=postfix,
1446+
targetdir=trj_dir,
1447+
weights=weights,
14181448
tol_mass=tol_mass, flat=flat)
14191449
paths.append(p.path)
14201450
fit_trj_names.append(fit_trj_name)
@@ -1426,7 +1456,6 @@ def generate_paths(self, **kwargs):
14261456
with open(self._fit_trjs_pkl, 'wb') as output:
14271457
cPickle.dump(self.fit_trj_names, output)
14281458
if store:
1429-
filename = kwargs.pop('filename', None)
14301459
self.save_paths(filename=filename)
14311460

14321461

0 commit comments

Comments
 (0)