Skip to content

Refactor RMSD analyses into MDAnalysis AnaysisBase classes#90

Open
hannahbaumann wants to merge 10 commits intomainfrom
rmsd_refactor_analysisbase
Open

Refactor RMSD analyses into MDAnalysis AnaysisBase classes#90
hannahbaumann wants to merge 10 commits intomainfrom
rmsd_refactor_analysisbase

Conversation

@hannahbaumann
Copy link
Copy Markdown
Contributor

@hannahbaumann hannahbaumann commented Feb 20, 2026

Refactor the rmsd analysis using this example: https://docs.mdanalysis.org/2.7.0/documentation_pages/analysis/base.html

  • Longer term we may not want to have the gather_rms_data function but access the individual analysis classes directly in the openfe Protocol
  • Should the make_Universe function go into utils?

@codecov
Copy link
Copy Markdown

codecov bot commented Feb 20, 2026

Codecov Report

✅ All modified and coverable lines are covered by tests.
✅ Project coverage is 96.37%. Comparing base (12122bd) to head (90ed39e).

Additional details and impacted files
@@            Coverage Diff             @@
##             main      #90      +/-   ##
==========================================
+ Coverage   96.09%   96.37%   +0.28%     
==========================================
  Files           6        6              
  Lines         333      359      +26     
==========================================
+ Hits          320      346      +26     
  Misses         13       13              

☔ View full report in Codecov by Sentry.
📢 Have feedback on the report? Share it here.

🚀 New features to boost your workflow:
  • ❄️ Test Analytics: Detect flaky tests, report on failures, and find test suite problems.

@hannahbaumann
Copy link
Copy Markdown
Contributor Author

@talagayev and @jthorton : This is a first go at the refactor into the individual RMSD classes, based on the MDAnalysis AnalysisBase. Please let me know what you think!

@hannahbaumann hannahbaumann changed the title Refactor RMSD analyses into MDAnalysis AnaysisBase classes Refactor RMSD analyses into MDAnalysis AnaysisBase classes Feb 20, 2026
Comment thread src/openfe_analysis/rmsd.py Outdated
Comment on lines +177 to +203
class LigandRMSD(AnalysisBase):
"""
1D RMSD time series for a ligand AtomGroup.
"""

def __init__(self, atomgroup, **kwargs):
super(LigandRMSD, self).__init__(atomgroup.universe.trajectory, **kwargs)

self._ag = atomgroup

def _prepare(self):
self.results.rmsd = []
self._reference = self._ag.positions
self._weights = self._ag.masses / np.mean(self._ag.masses)

def _single_frame(self):
rmsd = rms.rmsd(
self._ag.positions,
self._reference,
self._weights,
center=False,
superposition=False,
)
self.results.rmsd.append(rmsd)

def _conclude(self):
self.results.rmsd = np.asarray(self.results.rmsd)
Copy link
Copy Markdown
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

2 initial thoughts:

  • Can we make a more general RMSD class that can be reused for protein and ligand analysis, basically it should work on any atom group and does this not already exist in MDAnalysis?
  • Do we not want to switch to use the symmetry RMSD (spyrmsd I think its called)?

Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@jthorton : There's a difference between the RMSD class in MDAnalysis in what we've been doing, and I'm not sure yet what we want. In MDAnalysis, by default they are doing a superposition (rotational and translational), while we didn't do that. Now I'm not sure, what are we actually interested in, esp. for the ligand RMSD. Should this be the RMSD of the internal conformation of the ligand, or the RMSD of the ligand in the binding pocket/ligand pose? I think in the solvent we would definitely be more interested in the internal conformation of the ligand. In the binding pocket, I'm not sure. We also calculate the ligand COM displacement as a measure of stability in the pocket, but I'm not sure what we would want for the RMSD. What do you think?

Copy link
Copy Markdown
Collaborator

@talagayev talagayev left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@hannahbaumann I like it, looks very good to me :) also the structure is like for the MDAnalysis AnalysisBase classes, which is I think @IAlibay wanted to have and also adressing it that it takes any atom group adressing @jthorton comment is good.

Overall Looks good to me :)

Comment thread src/openfe_analysis/rmsd.py Outdated
prot_rmsd = RMSDAnalysis(prot).run(step=skip)
output["protein_RMSD"].append(prot_rmsd.results.rmsd)
# # Using the MDAnalysis RMSD class instead
# gs = ["protein and name CA", "protein"]
Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Here is an example of how we could do the same analysis using the MDAnalysis RMSD class.

Comment thread src/openfe_analysis/rmsd.py Outdated
self, atomgroup, reference=None, mass_weighted=False, superposition=False, **kwargs
):
super(RMSDAnalysis, self).__init__(atomgroup.universe.trajectory, **kwargs)

Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

_analysis_algorithm_is_parallelizable

@hannahbaumann
Copy link
Copy Markdown
Contributor Author

Add tests for mass weighting and superposition. Potentially import test data from MDAnalysis.

prot_rmsd2d = Protein2DRMSD(prot).run(step=skip)
output["protein_2D_RMSD"].append(prot_rmsd2d.results.rmsd2d)
# # Using the MDAnalysis DistanceMatrix class
# prot_rmsd2d = diffusionmap.DistanceMatrix(u, select="protein and name CA")
Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This MDA code is much slower, on the test data 10s vs. 0.4s.

output["protein_RMSD"].append(prot_rmsd.results.rmsd)
# # Using the MDAnalysis RMSD class instead
# gs = ["protein and name CA"]
# prot_rmsd = rms.RMSD(
Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The two RMSD classes are approximately equal in timing (on the test data)

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

4 participants