Spectrum objects
ms2ml.spectrum
Implements classes to represent spectra.
This module implements classes to represent spectra and their annotations (when they have any).
There are broadly two types of spectra: 1. General Spectra 2. Annotated Spectra
Classes
ms2ml.spectrum.Spectrum
dataclass
Class to store the spectrum information.
Examples:
>>> spectrum = Spectrum(
... mz=np.array([1000.0, 1500.0, 2000.0]),
... intensity=np.array([1.0, 2.0, 3.0]),
... precursor_mz=1000.0,
... precursor_charge=2,
... ms_level=2,
... extras={"EXTRAS": ["extra1", "extra2"]},
... )
>>> spectrum
Spectrum(mz=array([1000., 1500., 2000.]), ...)
Attributes
mz: NDArray[np.float64]
class-attribute
intensity: NDArray[np.float32]
class-attribute
ms_level: int
class-attribute
precursor_mz: float | None
class-attribute
precursor_charge: int | None = None
class-attribute
instrument: str | None = None
class-attribute
analyzer: str | None = None
class-attribute
collision_energy: float = field(default=float('nan'))
class-attribute
activation: str | None = None
class-attribute
extras: dict = field(default_factory=dict)
class-attribute
retention_time: RetentionTime | float = field(default_factory=lambda : RetentionTime(rt=float('nan'), units='seconds'))
class-attribute
config: Config = field(repr=False, default_factory=get_default_config)
class-attribute
Functions
__post_init__()
replace(*args: Any, **kwargs: Any) -> Spectrum
Replaces the attributes of the spectrum with the given values.
Arguments are passed to the class Spectrum constructor.
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
Spectrum, A new Spectrum object with the replaced attributes. |
filter_top(n: int) -> Spectrum
Filters the spectrum to the top n peaks.
PARAMETER | DESCRIPTION |
---|---|
n |
The number of peaks to keep.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
The spectrum with the filtered top n peaks. |
Examples:
filter_mz_range(min_mz: float, max_mz: float) -> Spectrum
Filters the spectrum to a given m/z range.
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
The spectrum with the filtered m/z range. |
Examples:
remove_precursor() -> Spectrum
Removes the precursor peak from the spectrum
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
Spectrum, A new Spectrum object with the precursor peak removed. |
intensity_cutoff(cutoff: float) -> Spectrum
Filters the spectrum to a given intensity cutoff.
PARAMETER | DESCRIPTION |
---|---|
cutoff |
The intensity cutoff. All peaks with less than that intensity will be deleted.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
A new Spectrum object with the filtered intensity cutoff. |
Examples:
normalize_intensity(method: str = 'max') -> Spectrum
Normalizes the spectrum intensities.
PARAMETER | DESCRIPTION |
---|---|
method |
The method to use for normalization. Can be one of "max", "sum", "rank", "log".
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
The normalized spectrum. |
Examples:
>>> spectrum = Spectrum._sample()
>>> spectrum.intensity
array([ 50., 200., 1., 2., 3.])
>>> Spectrum._sample().normalize_intensity("max").intensity
array([0.25 , 1. , 0.005, 0.01 , 0.015])
>>> Spectrum._sample().normalize_intensity("sum").intensity
array([0.1953125 , 0.78125 , 0.00390625, 0.0078125 , 0.01171875])
>>> Spectrum._sample().normalize_intensity("log").intensity
array([3.91202301, 5.29831737, 0. , 0.69314718, 1.09861229])
>>> Spectrum._sample().normalize_intensity("sqrt").intensity
array([ 7.07106781, 14.14213562, 1. , 1.41421356, 1.73205081])
encode_spec_bins() -> NDArray[np.float32]
bin_spectrum(start: float, end: float, binsize: float | None = None, n_bins: int | None = None, relative: bool = False, offset: float = 0, get_breaks: bool = False) -> NDArray[np.float32] | tuple[NDArray[np.float32], NDArray[np.float64]]
Bins the spectrum.
PARAMETER | DESCRIPTION |
---|---|
start |
The start of the binning range. If missing will use the lowest mz value.
TYPE:
|
end |
The end of the binning range. If missing will use the highest mz value.
TYPE:
|
binsize |
The size of the bins. Cannot be used in conjunction with n_bins.
TYPE:
|
n_bins |
The number of bins. Cannot be used in conjunction with binsize.
TYPE:
|
relative |
Whether to use binning relative to the precursor mass.
TYPE:
|
offset |
The offset to use for relative binning.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
NDArray[np.float32] | tuple[NDArray[np.float32], NDArray[np.float64]]
|
An array of binned intensities. |
base_peak(value) -> None
tic(value) -> None
sic(mzs: NDArray[np.float32], resolution: Callable = sum) -> NDArray[np.float32]
Returns the selected ion current for a given set of m/z values.
PARAMETER | DESCRIPTION |
---|---|
mzs |
The m/z values to calculate the SIC for.
TYPE:
|
resolution |
The function used to resolve ambiguities when multiple
peaks match. possible options are
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
NDArray[np.float32]
|
An array of SIC values. This array will have the same length as the |
NDArray[np.float32]
|
input mzs array. |
Examples:
annotate(peptide: str | Peptide) -> AnnotatedPeptideSpectrum
Annotates the spectrum with the given peptide.
PARAMETER | DESCRIPTION |
---|---|
peptide |
The peptide to annotate the spectrum with.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
AnnotatedPeptideSpectrum
|
An AnnotatedPeptideSpectrum object. |
Examples:
to_sus()
plot(ax = None, **kwargs) -> plt.Axes
stack(spectra: Iterable[Spectrum]) -> Spectrum
classmethod
Stacks multiple spectra into a single spectrum.
PARAMETER | DESCRIPTION |
---|---|
spectra |
An iterable of spectra to stack.
TYPE:
|
RETURNS | DESCRIPTION |
---|---|
Spectrum
|
A stacked spectrum. |
Examples:
>>> np.random.seed(0)
>>> spectrum1 = Spectrum._sample()
>>> spectrum2 = Spectrum._sample(random=True)
>>> spectrum3 = Spectrum._sample(random=True)
>>> spectrum4 = Spectrum._sample(random=False)
>>> spectrum4.mz = spectrum1.mz + 1e-4
>>> mzs, ints = Spectrum.stack([spectrum1, spectrum2, spectrum3, spectrum4])
>>> mzs
array([ 50. , 147.11333 , 423.65479934, 437.58721126,
544.883183 , 592.84461823, 645.89411307, 1000. ,
1500. , 2000. ])
>>> ints
array([[ 50. , 0. , 0. , 50. ],
[200. , 0. , 0. , 200. ],
[ 0. , 0. , 963.6627605 , 0. ],
[ 0. , 0. , 791.72503808, 0. ],
[ 0. , 0. , 891.77300078, 0. ],
[ 0. , 844.26574858, 0. , 0. ],
[ 0. , 0. , 383.44151883, 0. ],
[ 1. , 0. , 0. , 1. ],
[ 2. , 0. , 0. , 2. ],
[ 3. , 0. , 0. , 3. ]])
>>> ints.shape
(10, 4)
ms2ml.spectrum.AnnotatedPeptideSpectrum
dataclass
Bases: Spectrum
Class to store the spectrum information.
In combination with the peptide information, it can be used to annotate the spectrum.
Examples:
>>> config = Config()
>>> peptide = Peptide.from_sequence("PEPPINK/2", config)
>>> spectrum = AnnotatedPeptideSpectrum(
... mz=np.array([50.0, 147.11333, 1000.0, 1500.0, 2000.0]),
... intensity=np.array([50.0, 200.0, 1.0, 2.0, 3.0]),
... ms_level=2,
... extras={"EXTRAS": ["extra1", "extra2"]},
... precursor_peptide=peptide,
... precursor_mz=147.11333,
... )
>>> spectrum
AnnotatedPeptideSpectrum(mz=array([ 50. ... precursor_isotope=0)
>>> spectrum.fragment_intensities
{'y1^1': 200.0}
>>> spectrum["y1^1"]
200.0
>>> spectrum.fragments
{'y1^1': AnnotatedIon(mass=147.11334,
charge=1, position=1, ion_series='y', intensity=200.0, neutral_loss=None)}
Attributes
precursor_peptide: Peptide | None = None
class-attribute
precursor_isotope: int | None = 0
class-attribute
precursor_charge: int | None = None
class-attribute
charge
property
mass_error
property
fragment_labels: list[str]
property
Functions
__post_init__(*args, **kwargs)
fragment_intensities() -> dict[str, float]
Returs a dictionary with the fragment ion names as keys and the corresponding intensities as values.
Note
The current implementation only keeps the last peak that matches the theoretical mass. future implementations should either keep all peaks or only the highest peak, or add the peaks.
Examples:
fragments() -> dict[str, AnnotatedIon]
__getitem__(index) -> float
encode_fragments() -> NDArray[np.float32]
Encodes the fragment ions as a numpy array
The order of the ions will be defined in the config file.
Examples:
>>> spec = AnnotatedPeptideSpectrum._sample()
>>> spec.encode_fragments()
array([200., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
dtype=float32)
decode_fragments(peptide: Peptide, fragment_vector: NDArray[np.float32])
staticmethod
Examples:
>>> spec = AnnotatedPeptideSpectrum._sample()
>>> pep = spec.precursor_peptide
>>> frags = spec.encode_fragments()
>>> AnnotatedPeptideSpectrum.decode_fragments(pep, frags)
AnnotatedPeptideSpectrum(mz=array([...]),
intensity=array([...], dtype=float32),
ms_level=2, precursor_mz=397.724526907315,
precursor_charge=2,
instrument=None,
analyzer=None,
collision_energy=nan,
activation=None,
extras={},
retention_time=RetentionTime(rt=nan, units='seconds',
run=None),
precursor_peptide=Peptide([...], {...}), precursor_isotope=0)
to_sus()
Return a spectrum object as a spectrum_utils.spectrum.Spectrum object