4.7.1. Polymer analysis — MDAnalysis.analysis.polymer
- Author:
Richard J. Gowers
- Year:
2015, 2018
- Copyright:
Lesser GNU Public License v2.1+
This module contains various commonly used tools in analysing polymers.
- class MDAnalysis.analysis.polymer.PersistenceLength(atomgroups, **kwargs)[source]
Calculate the persistence length for polymer chains
The persistence length is the length at which two points on the polymer chain become decorrelated. This is determined by first measuring the autocorrelation (\(C(n)\)) of two bond vectors (\(\mathbf{a}_i, \mathbf{a}_{i + n}\)) separated by \(n\) bonds
\[C(n) = \langle \cos\theta_{i, i+n} \rangle = \langle \mathbf{a_i} \cdot \mathbf{a_{i+n}} \rangle\]An exponential decay is then fitted to this, which yields the persistence length
\[C(n) \approx \exp\left( - \frac{n \bar{l_B}}{l_P} \right)\]where \(\bar{l_B}\) is the average bond length, and \(l_P\) is the persistence length which is fitted
- Parameters:
atomgroups (iterable) – List of AtomGroups. Each should represent a single polymer chain, ordered in the correct order.
verbose (bool, optional) – Show detailed progress of the calculation if set to
True
.
- results.bond_autocorrelation
the measured bond autocorrelation
- Type:
- lb
Alias to the
results.lb
.Deprecated since version 2.0.0: Will be removed in MDAnalysis 3.0.0. Please use
results.lb
instead.- Type:
- results.x
length of the decorrelation predicted by lp
Added in version 2.0.0.
- Type:
- lp
Alias to the
results.lp
.Deprecated since version 2.0.0: Will be removed in MDAnalysis 3.0.0. Please use
results.lp
instead.- Type:
- results.fit
the modelled backbone decorrelation predicted by lp
Added in version 2.0.0.
- Type:
- fit
Alias to the
results.fit
.Deprecated since version 2.0.0: Will be removed in MDAnalysis 3.0.0. Please use
results.fit
instead.- Type:
See also
sort_backbone()
for producing the sorted AtomGroup required for input.
Example
from MDAnalysis.tests.datafiles import TRZ_psf, TRZ import MDAnalysis as mda from MDAnalysis.analysis import polymer u = mda.Universe(TRZ_psf, TRZ) # this system is a pure polymer melt of polyamide, # so we can select the chains by using the .fragments attribute chains = u.atoms.fragments # select only the backbone atoms for each chain backbones = [chain.select_atoms('not name O* H*') for chain in chains] # sort the chains, removing any non-backbone atoms sorted_backbones = [polymer.sort_backbone(bb) for bb in backbones] persistence_length = polymer.PersistenceLength(sorted_backbones) # Run the analysis, this will average over all polymer chains # and all timesteps in trajectory persistence_length = persistence_length.run() print(f'The persistence length is: {persistence_length.results.lp}') # always check the visualisation of this: persistence_length.plot()
Added in version 0.13.0.
Changed in version 0.20.0: The run method now automatically performs the exponential fit
Changed in version 1.0.0: Deprecated
PersistenceLength.perform_fit()
has now been removed.Changed in version 2.0.0: Former
results
are now stored asresults.bond_autocorrelation
.lb
,lp
,fit
are now stored in aMDAnalysis.analysis.base.Results
instance.- classmethod get_supported_backends()
Tuple with backends supported by the core library for a given class. User can pass either one of these values as
backend=...
torun()
method, or a custom object that hasapply
method (see documentation forrun()
):‘serial’: no parallelization
‘multiprocessing’: parallelization using multiprocessing.Pool
‘dask’: parallelization using dask.delayed.compute(). Requires installation of mdanalysis[dask]
If you want to add your own backend to an existing class, pass a
backends.BackendBase
subclass (see its documentation to learn how to implement it properly), and specifyunsupported_backend=True
.- Returns:
names of built-in backends that can be used in
run(backend=...)()
- Return type:
Added in version 2.8.0.
- property parallelizable
Boolean mark showing that a given class can be parallelizable with split-apply-combine procedure. Namely, if we can safely distribute
_single_frame()
to multiple workers and then combine them with a proper_conclude()
call. If set toFalse
, no backends except forserial
are supported.Note
If you want to check parallelizability of the whole class, without explicitly creating an instance of the class, see
_analysis_algorithm_is_parallelizable
. Note that you setting it to other value will break things if the algorithm behind the analysis is not trivially parallelizable.- Returns:
if a given
AnalysisBase
subclass instance is parallelizable with split-apply-combine, or not- Return type:
Added in version 2.8.0.
- plot(ax=None)[source]
Visualize the results and fit
- Parameters:
ax (matplotlib.Axes, optional) – if provided, the graph is plotted on this axis
- Returns:
ax
- Return type:
the axis that the graph was plotted on
- run(start: int | None = None, stop: int | None = None, step: int | None = None, frames: Iterable | None = None, verbose: bool | None = None, n_workers: int | None = None, n_parts: int | None = None, backend: str | BackendBase | None = None, *, unsupported_backend: bool = False, progressbar_kwargs=None)
Perform the calculation
- Parameters:
start (int, optional) – start frame of analysis
stop (int, optional) – stop frame of analysis
step (int, optional) – number of frames to skip between each analysed frame
frames (array_like, optional) –
array of integers or booleans to slice trajectory;
frames
can only be used instead ofstart
,stop
, andstep
. Setting bothframes
and at least one ofstart
,stop
,step
to a non-default value will raise aValueError
.Added in version 2.2.0.
verbose (bool, optional) – Turn on verbosity
progressbar_kwargs (dict, optional) – ProgressBar keywords with custom parameters regarding progress bar position, etc; see
MDAnalysis.lib.log.ProgressBar
for full list. Available only forbackend='serial'
backend (Union[str, BackendBase], optional) –
By default, performs calculations in a serial fashion. Otherwise, user can choose a backend:
str
is matched to a builtin backend (one ofserial
,multiprocessing
anddask
), or aMDAnalysis.analysis.results.BackendBase
subclass.Added in version 2.8.0.
n_workers (int) –
positive integer with number of workers (processes, in case of built-in backends) to split the work between
Added in version 2.8.0.
n_parts (int, optional) –
number of parts to split computations across. Can be more than number of workers.
Added in version 2.8.0.
unsupported_backend (bool, optional) –
if you want to run your custom backend on a parallelizable class that has not been tested by developers, by default False
Added in version 2.8.0.
Changed in version 2.2.0: Added ability to analyze arbitrary frames by passing a list of frame indices in the frames keyword argument.
Changed in version 2.5.0: Add progressbar_kwargs parameter, allowing to modify description, position etc of tqdm progressbars
Changed in version 2.8.0: Introduced
backend
,n_workers
,n_parts
andunsupported_backend
keywords, and refactored the method logic to support parallelizable execution.
- MDAnalysis.analysis.polymer.fit_exponential_decay(x, y)[source]
Fit a function to an exponential decay
\[y = \exp\left(- \frac{x}{a}\right)\]- Parameters:
x (array_like) – The two arrays of data
y (array_like) – The two arrays of data
- Returns:
a – The coefficient a for this decay
- Return type:
Notes
This function assumes that data starts at 1.0 and decays to 0.0
- MDAnalysis.analysis.polymer.sort_backbone(backbone)[source]
Rearrange a linear AtomGroup into backbone order
Requires that the backbone has bond information, and that only backbone atoms are provided (ie no side chains or hydrogens).
- Parameters:
backbone (AtomGroup) – the backbone atoms, not necessarily in order
- Returns:
sorted_backbone – backbone in order, so sorted_backbone[i] is bonded to sorted_backbone[i - 1] and sorted_backbone[i + 1]
- Return type:
Added in version 0.20.0.