EndToEndVector¶
- class mdcraft.analysis.polymer.EndToEndVector(groups: AtomGroup | tuple[AtomGroup], groupings: str | tuple[str] = 'atoms', n_chains: int | tuple[int] = None, n_monomers: int | tuple[int] = None, *, n_blocks: int = 1, dt: float | Quantity | Quantity = None, fft: bool = True, unwrap: bool = False, verbose: bool = True, **kwargs)[source]¶
Bases:
_PolymerAnalysisBaseA serial implementation to calculate the end-to-end vector autocorrelation function (ACF) \(C_\mathrm{ee}(t)\) and estimate the orientational relaxation time \(\tau_\mathrm{r}\) of a polymer.
The end-to-end vector ACF is defined as
\[C_\mathrm{ee}(t)=\frac{\langle\mathbf{R}_\mathrm{ee}(t) \cdot\mathbf{R}_\mathrm{ee}(0)\rangle} {\langle\mathbf{R}_\mathrm{ee}^2\rangle}\]where \(\mathbf{R}_\mathrm{ee}=\mathbf{r}_N-\mathbf{r}_1\) is the end-to-end vector.
The orientational relaxation time can then be estimated by fitting a stretched exponential function
\[C_\mathrm{ee}=\exp{\left[-(t/\tau)^\beta\right]}\]to the end-to-end vector ACF and evaluating
\[\tau_\mathrm{r}=\int_0^\infty C_\mathrm{ee}\,dt =\frac{\tau}{\beta}\Gamma\left(\frac{1}{\beta}\right)\]- Parameters:
- groupsMDAnalysis.AtomGroup or array-like
Groups of polymers to be analyzed.
Note
All polymers in each group must have the same chain length.
- groupingsstr or array-like, default:
"atoms" Determines whether the centers of mass are used in lieu of individual atom positions. If groupings is a str, the same value is used for all groups.
Note
In a standard trajectory file, segments (or chains) contain residues (or molecules), and residues contain atoms. This heirarchy must be adhered to for this analysis module to function correctly. If your trajectory file does not contain the correct residue or segment information, provide the number of chains and chain lengths in n_chains and n_monomers, respectively.
Valid values:
"atoms": Atom positions (generally or for coarse-grained simulations)."residues": Residues’ centers of mass (for atomistic simulations).
- n_chainsint or array-like, optional
Number of chains \(M\) in each polymer group. Must be provided if the trajectory does not adhere to the standard container heirarchy. If an int is provided, the same value is used for all groups.
Shape: \((N_\mathrm{groups},)\).
- n_monomersint or array-like, optional
Number of monomers \(N\) in each chain in each polymer group. Must be provided if the trajectory does not adhere to the standard container heirarchy. If an int is provided, the same value is used for all groups.
Shape: \((N_\mathrm{groups},)\).
- n_blocksint, keyword-only, default:
1 Number of blocks to split the trajectory into.
- dtfloat, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional
Time between frames \(\Delta t\). While this is normally determined from the trajectory, the trajectory may not have the correct information if the data is in reduced units. For example, if the reduced timestep is \(0.01\) and trajectory data was outputted every \(10,000\) timesteps, then \(\Delta t=100\).
Reference unit: \(\mathrm{ps}\).
- fftbool, keyword-only, default:
True Determines whether fast Fourier transforms (FFT) are used to evaluate the ACFs.
- unwrapbool, keyword-only, default:
False Determines whether atom positions are unwrapped.
- verbosebool, keyword-only, default:
True Determines whether detailed progress is shown.
- **kwargs
Additional keyword arguments to pass to
MDAnalysis.analysis.base.AnalysisBase.
- Attributes:
- universeMDAnalysis.Universe
MDAnalysis.core.universe.Universeobject containing all information describing the system.- results.unitsdict
Reference units for the results. For example, to get the reference units for results.times, call
results.units["times"].- results.timesnumpy.ndarray
Changes in time \(t-t_0\).
Shape: \((N_\mathrm{frames},)\).
Reference unit: \(\mathrm{ps}\).
- results.acfnumpy.ndarray
End-to-end vector ACFs \(C_\mathrm{ee}(t)\).
Shape: \((N_\mathrm{groups},\,N_\mathrm{blocks},\,N_\mathrm{frames})\).
- results.relaxation_timesnumpy.ndarray
Average orientational relaxation times \(\tau_\mathrm{r}\).
Shape: \((N_\mathrm{groups},\,N_\mathrm{blocks})\).
Reference units: \(\mathrm{ps}\).
Methods
Calculates the orientational relaxation times.
Tuple with backends supported by the core library for a given class.
Performs the calculation.
Saves results to a binary or archive file in NumPy format.
- 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 hasapplymethod (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.BackendBasesubclass (see its documentation to learn how to implement it properly), and specifyunsupported_backend=True.- Returns:
- tuple
names of built-in backends that can be used in
run(backend=...)()
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 forserialare 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:
- bool
if a given
AnalysisBasesubclass instance is parallelizable with split-apply-combine, or not
Added in version 2.8.0: ..
- run(start: int = None, stop: int = None, step: int = None, frames: slice | ndarray[int] = None, verbose: bool = None, **kwargs) SerialAnalysisBase | ParallelAnalysisBase¶
Performs the calculation.
See also
For parallel-specific keyword arguments, see
ParallelAnalysisBase.run().- Parameters:
- startint, optional
Starting frame for analysis.
- stopint, optional
Ending frame for analysis.
- stepint, optional
Number of frames to skip between each analyzed frame.
- framesslice or array-like, optional
Index or logical array of the desired trajectory frames.
- verbosebool, optional
Determines whether detailed progress is shown.
- **kwargs
Additional keyword arguments to pass to
MDAnalysis.lib.log.ProgressBar.
- Returns:
- selfSerialAnalysisBase or ParallelAnalysisBase
Analysis object with results.
- save(file: str | TextIO, archive: bool = True, compress: bool = True, **kwargs) None¶
Saves results to a binary or archive file in NumPy format.
- Parameters:
- filestr or file
Filename or file-like object where the data will be saved. If file is a str, the
.npyor.npzextension will be appended automatically if not already present.- archivebool, default:
True Determines whether the results are saved to a single archive file. If True, the data is stored in a
.npzfile. Otherwise, the data is saved to multiple.npyfiles.- compressbool, default:
True Determines whether the
.npzfile is compressed. Has no effect whenarchive=False.- **kwargs
Additional keyword arguments to pass to
numpy.save(),numpy.savez(), ornumpy.savez_compressed(), depending on the values of archive and compress.