Skip to content

Documentation for ChemicalSamplers.py

ChemicalSamplers

Sampling classes for ChemicalMixtureCreator

CRPMS2Sampler

CRPMS2Sampler(
    n_draws=1000,
    min_mz=MIN_MZ_MS2,
    min_proportion=0.1,
    max_proportion=0.8,
    alpha=1,
    base="uniform",
)

Bases: MS2Sampler

A sampler that generates MS2 peaks following the Chinese Restaurant Process (CRP), i.e. an MS2 peak that has been selected in one spectra has a higher likelihood to appear again elsewhere.

Create a CRP-based MS2 sampler. Args: n_draws: the number of draws from the CRP process min_mz: the minimum m/z value to consider min_proportion: the minimum proportion to consider max_proportion: the maximum proportion to consider alpha: CRP parameter base: base distribution for the CRP process

Source code in vimms/ChemicalSamplers.py
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
def __init__(
    self,
    n_draws=1000,
    min_mz=MIN_MZ_MS2,
    min_proportion=0.1,
    max_proportion=0.8,
    alpha=1,
    base="uniform",
):
    """
    Create a CRP-based MS2 sampler.
    Args:
        n_draws: the number of draws from the CRP process
        min_mz: the minimum m/z value to consider
        min_proportion: the minimum proportion to consider
        max_proportion: the maximum proportion to consider
        alpha: CRP parameter
        base: base distribution for the CRP process
    """
    self.n_draws = n_draws
    self.min_mz = min_mz
    self.min_proportion = min_proportion
    self.max_proportion = max_proportion
    self.alpha = alpha
    assert self.alpha > 0
    self.base = base
    assert self.base == "uniform"

sample

sample(chemical)

Sample MS2 spectra using chemical as the parent Args: chemical: the parent chemical

Returns: a tuple of (mz_list, intensity_list, parent_proportion)

Source code in vimms/ChemicalSamplers.py
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
def sample(self, chemical):
    """
    Sample MS2 spectra using chemical as the parent
    Args:
        chemical: the parent chemical

    Returns: a tuple of (mz_list, intensity_list, parent_proportion)

    """

    max_mz = chemical.mass
    unique_vals = [self._base_sample(max_mz)]
    counts = [1]
    for i in range(self.n_draws - 1):
        temp = counts + [self.alpha]
        s = sum(temp)
        probs = [t / s for t in temp]
        choice = np.random.choice(len(temp), p=probs)
        if choice == len(unique_vals):
            # new value
            unique_vals.append(self._base_sample(max_mz))
            counts.append(1)
        else:
            counts[choice] += 1

    mz_list = unique_vals
    s = sum(counts)
    intensity_list = [c / s for c in counts]
    parent_proportion = (
        np.random.rand() * (self.max_proportion - self.min_proportion) + self.min_proportion
    )

    return mz_list, intensity_list, parent_proportion

ChromatogramSampler

Bases: ABC

Base class for chromatogram sampler.

ConstantChromatogramSampler

Bases: ChromatogramSampler

A sampler to return constant chromatograms -- direct infusion

sample

sample(formula, rt, intensity)

Sample a constant chromatogram (present everywhere) Args: formula: formula, unused rt: RT, unused intensity: intensity, unused

Returns: a [vimms.Chromatograms.ConstantChromatogram] object.

Source code in vimms/ChemicalSamplers.py
459
460
461
462
463
464
465
466
467
468
469
470
def sample(self, formula, rt, intensity):
    """
    Sample a constant chromatogram (present everywhere)
    Args:
        formula: formula, unused
        rt: RT, unused
        intensity: intensity, unused

    Returns: a [vimms.Chromatograms.ConstantChromatogram] object.

    """
    return ConstantChromatogram()

DatabaseFormulaSampler

DatabaseFormulaSampler(database, min_mz=MIN_MZ, max_mz=MAX_MZ)

Bases: FormulaSampler

A sampler to draw formula from a database

Initiliases database formula sampler

Parameters:

Name Type Description Default
database

a list of Formula objects containing chemical formulae from e.g. HMDB

required
min_mz

the minimum m/z value of formulae to sample from

MIN_MZ
max_mz

the maximum m/z value of formulae to sample from

MAX_MZ
Source code in vimms/ChemicalSamplers.py
61
62
63
64
65
66
67
68
69
70
71
72
def __init__(self, database, min_mz=MIN_MZ, max_mz=MAX_MZ):
    """
    Initiliases database formula sampler

    Args:
        database: a list of Formula objects containing chemical
                  formulae from e.g. HMDB
        min_mz: the minimum m/z value of formulae to sample from
        max_mz: the maximum m/z value of formulae to sample from
    """
    super().__init__(min_mz=min_mz, max_mz=max_mz)
    self.database = database

sample

sample(n_formulas)

Samples n_formulas from the specified database

Parameters:

Name Type Description Default
n_formulas

the number of formula to draw

required

Returns: a list of Formula objects

Source code in vimms/ChemicalSamplers.py
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
def sample(self, n_formulas):
    """
    Samples n_formulas from the specified database

    Args:
        n_formulas: the number of formula to draw

    Returns: a list of Formula objects

    """
    # filter database formulae to be within mz_range
    offset = 20  # to ensure that we have room for at least M+H
    formulas = list(set([(x.chemical_formula, x.name) for x in self.database]))
    sub_formulas = list(
        filter(
            lambda x: Formula(x[0]).mass >= self.min_mz
            and Formula(x[0]).mass <= self.max_mz - offset,
            formulas,
        )
    )
    logger.debug("{} unique formulas in filtered database".format(len(sub_formulas)))
    chosen_formula_positions = np.random.choice(
        len(sub_formulas), size=n_formulas, replace=False
    )
    logger.debug("Sampled formulas")
    return [
        (Formula(sub_formulas[f][0]), sub_formulas[f][1]) for f in chosen_formula_positions
    ]

DefaultScanTimeSampler

DefaultScanTimeSampler(scan_time_dict=None)

Bases: ScanTimeSampler

A scan time sampler that returns some fixed values that represent the average scan times for MS1 and MS2 scans.

Initialises a default scan time sampler object.

Parameters:

Name Type Description Default
scan_time_dict

A dictionary of scan times for each MS-level. It should look like this: {1: 0.4, 2: 0.2}. If not specified, then the default value is used. Note that this default is obtained from our Orbitrap instrument and would certainly differ from yours!

None
Source code in vimms/ChemicalSamplers.py
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
def __init__(self, scan_time_dict=None):
    """
    Initialises a default scan time sampler object.

    Args:
        scan_time_dict: A dictionary of scan times for each MS-level.
                        It should look like this: {1: 0.4, 2: 0.2}.
                        If not specified, then the default value is used.
                        Note that this default is obtained from our Orbitrap instrument and
                        would certainly differ from yours!
    """

    self.scan_time_dict = (
        scan_time_dict if scan_time_dict is not None else DEFAULT_SCAN_TIME_DICT
    )

sample

sample(current_level, next_level, current_rt)

Sample a scan duration given the MS levels of current and next scans. Args: current_level: the MS level of the current scan next_level: the MS level of the next scan current_rt: not used

Returns: a sampled scan duration value

Source code in vimms/ChemicalSamplers.py
945
946
947
948
949
950
951
952
953
954
955
956
def sample(self, current_level, next_level, current_rt):
    """
    Sample a scan duration given the MS levels of current and next scans.
    Args:
        current_level: the MS level of the current scan
        next_level: the MS level of the next scan
        current_rt: not used

    Returns: a sampled scan duration value

    """
    return self.scan_time_dict[current_level]

EvenMZFormulaSampler

EvenMZFormulaSampler()

Bases: FormulaSampler

A sampler that picks mz values evenly spaced, starting from where it left off. Useful for test cases

Create an even m/z formula sampler

Source code in vimms/ChemicalSamplers.py
165
166
167
168
169
170
def __init__(self):
    """
    Create an even m/z formula sampler
    """
    self.n_sampled = 0
    self.step = 100

sample

sample(n_formulas)

Sample up to n_formulas from this sampler Args: n_formulas: the number of formula to return

Returns: the list of formulae having evenly spaced m/z values

Source code in vimms/ChemicalSamplers.py
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
def sample(self, n_formulas):
    """
    Sample up to n_formulas from this sampler
    Args:
        n_formulas: the number of formula to return

    Returns: the list of formulae having evenly spaced m/z values

    """
    mz_list = []
    for i in range(n_formulas):
        new_mz = (self.n_sampled + 1) * self.step
        mz_list.append(new_mz)
        self.n_sampled += 1
    return [(DummyFormula(m), None) for m in mz_list]

ExactMatchMS2Sampler

ExactMatchMS2Sampler(
    mgf_file, min_proportion=0.1, max_proportion=0.8, id_field="SPECTRUMID"
)

Bases: MGFMS2Sampler

Exact match MS2 sampler allows us to have particular formulas and we have a particular spectrum for each exact formula...

TODO: not sure if this class is actually completed and fully tested.

Source code in vimms/ChemicalSamplers.py
794
795
796
797
798
799
800
def __init__(self, mgf_file, min_proportion=0.1, max_proportion=0.8, id_field="SPECTRUMID"):
    super().__init__(
        mgf_file,
        min_proportion=min_proportion,
        max_proportion=max_proportion,
        id_field=id_field,
    )

sample

sample(chemical)

Sample MS2 spectra using chemical as the parent Args: chemical: the parent chemical

Returns: a tuple of (mz_list, intensity_list, parent_proportion)

Source code in vimms/ChemicalSamplers.py
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
def sample(self, chemical):
    """
    Sample MS2 spectra using chemical as the parent
    Args:
        chemical: the parent chemical

    Returns: a tuple of (mz_list, intensity_list, parent_proportion)

    """

    spectrum = self.spectra_dict[chemical.database_accession]
    mz_list, intensity_list = zip(*spectrum.peaks)
    parent_proportion = (
        np.random.rand() * (self.max_proportion - self.min_proportion) + self.min_proportion
    )
    return mz_list, intensity_list, parent_proportion

FixedMS2Sampler

FixedMS2Sampler(n_frags=2)

Bases: MS2Sampler

Generates n_frags fragments, where each is chemical - i*10 mz

Create a fixed MS2 sampler Args: n_frags: the number of fragment peaks to generate

Source code in vimms/ChemicalSamplers.py
600
601
602
603
604
605
606
def __init__(self, n_frags=2):
    """
    Create a fixed MS2 sampler
    Args:
        n_frags: the number of fragment peaks to generate
    """
    self.n_frags = n_frags

sample

sample(chemical)

Sample MS2 spectra using chemical as the parent Args: chemical: the parent chemical

Returns: a tuple of (mz_list, intensity_list, parent_proportion)

Source code in vimms/ChemicalSamplers.py
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
def sample(self, chemical):
    """
    Sample MS2 spectra using chemical as the parent
    Args:
        chemical: the parent chemical

    Returns: a tuple of (mz_list, intensity_list, parent_proportion)

    """
    initial_mz = chemical.mass
    mz_list = []
    intensity_list = []
    parent_proportion = 0.5
    for i in range(self.n_frags):
        mz_list.append(initial_mz - (i + 1) * 10)
        intensity_list.append(1)
    s = sum(intensity_list)
    intensity_list = [i / s for i in intensity_list]
    return mz_list, intensity_list, parent_proportion

FormulaSampler

FormulaSampler(min_mz=MIN_MZ, max_mz=MAX_MZ)

Bases: ABC

Base class for formula sampler

Create a Formula sampler Args: min_mz: the minimum m/z value of formulae to sample from max_mz: the maximum m/z value of formulae to sample from

Source code in vimms/ChemicalSamplers.py
41
42
43
44
45
46
47
48
49
def __init__(self, min_mz=MIN_MZ, max_mz=MAX_MZ):
    """
    Create a Formula sampler
    Args:
        min_mz: the minimum m/z value of formulae to sample from
        max_mz: the maximum m/z value of formulae to sample from
    """
    self.min_mz = min_mz
    self.max_mz = max_mz

GaussianChromatogramSampler

GaussianChromatogramSampler(sigma=10)

Bases: ChromatogramSampler

A sampler to return Gaussian-shaped chromatogram

Create a Gaussian-shaped chromatogram sampler Args: sigma: parameter for the Gaussian distribution to sample from

Source code in vimms/ChemicalSamplers.py
430
431
432
433
434
435
436
437
def __init__(self, sigma=10):
    """
    Create a Gaussian-shaped chromatogram sampler
    Args:
        sigma: parameter for the Gaussian distribution to sample from
    """
    assert sigma > 0
    self.sigma = sigma

sample

sample(formula, rt, intensity)

Sample a Gaussian-shaped chromatogram

Parameters:

Name Type Description Default
formula

the formula to condition on (can be ignored)

required
rt

RT to condition on (can be ignored)

required
intensity

intensity to condition on (can be ignored)

required

Returns: a [vimms.Chromatograms.FunctionalChromatogram] object.

Source code in vimms/ChemicalSamplers.py
439
440
441
442
443
444
445
446
447
448
449
450
451
def sample(self, formula, rt, intensity):
    """
    Sample a Gaussian-shaped chromatogram

    Args:
        formula: the formula to condition on (can be ignored)
        rt: RT to condition on (can be ignored)
        intensity: intensity to condition on (can be ignored)

    Returns: a [vimms.Chromatograms.FunctionalChromatogram] object.

    """
    return FunctionalChromatogram("normal", [0, self.sigma])

MGFMS2Sampler

MGFMS2Sampler(
    mgf_file,
    min_proportion=0.1,
    max_proportion=0.8,
    max_peaks=0,
    replace=False,
    id_field="SPECTRUMID",
)

Bases: MS2Sampler

A sampler that generates MS2 spectra from real ones defined in some MGF file.

Create an MGFMS2Sampler object. Args: mgf_file: input MGF file. min_proportion: the minimum proportion to consider max_proportion: the maximum proportion to consider max_peaks: the maximum number of peaks replace: whether to sample with replacement or not id_field: the ID field in the MGF file

Source code in vimms/ChemicalSamplers.py
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
def __init__(
    self,
    mgf_file,
    min_proportion=0.1,
    max_proportion=0.8,
    max_peaks=0,
    replace=False,
    id_field="SPECTRUMID",
):
    """
    Create an MGFMS2Sampler object.
    Args:
        mgf_file: input MGF file.
        min_proportion: the minimum proportion to consider
        max_proportion: the maximum proportion to consider
        max_peaks: the maximum number of peaks
        replace: whether to sample with replacement or not
        id_field: the ID field in the MGF file
    """
    self.mgf_file = mgf_file
    self.min_proportion = min_proportion
    self.max_proportion = max_proportion
    self.replace = replace  # sample with replacement

    # load the mgf
    self.spectra_dict = load_mgf(self.mgf_file, id_field=id_field)

    # turn into a list where the last item is the number of times
    # this one has been sampled
    self.spectra_list = [[s.precursor_mz, s, 0] for s in self.spectra_dict.values()]

    # filter to remove those with more than  max_peaks (if max_peaks > 0)
    if max_peaks > 0:
        self.spectra_list = list(
            filter(lambda x: len(x[1].peaks) <= max_peaks, self.spectra_list)
        )

    # sort by precursor mz
    self.spectra_list.sort(key=lambda x: x[0])
    logger.debug("Loaded {} spectra from {}".format(len(self.spectra_list), self.mgf_file))

sample

sample(chemical)

Sample MS2 spectra using chemical as the parent Args: chemical: the parent chemical

Returns: a tuple of (mz_list, intensity_list, parent_proportion)

Source code in vimms/ChemicalSamplers.py
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
def sample(self, chemical):
    """
    Sample MS2 spectra using chemical as the parent
    Args:
        chemical: the parent chemical

    Returns: a tuple of (mz_list, intensity_list, parent_proportion)

    """

    formula_mz = chemical.mass
    sub_spec = list(filter(lambda x: x[0] < formula_mz, self.spectra_list))
    if len(sub_spec) == 0:
        # if there aren't any smaller than the mz, we just take any one
        sub_spec = self.spectra_list

    # sample one. If replace == True we take any, if not we only
    # take those that have not been sampled before
    found_permissable = False
    n_attempts = 0
    while not found_permissable:
        n_attempts += 1
        spec = np.random.choice(len(sub_spec))
        if self.replace is True or sub_spec[spec][2] == 0 or n_attempts > 100:
            found_permissable = True

    sub_spec[spec][2] += 1  # add one to the count
    spectrum = sub_spec[spec][1]
    mz_list, intensity_list = zip(*spectrum.peaks)
    s = sum(intensity_list)
    intensity_list = [i / s for i in intensity_list]
    parent_proportion = (
        np.random.rand() * (self.max_proportion - self.min_proportion) + self.min_proportion
    )

    return mz_list, intensity_list, parent_proportion

MS2Sampler

Bases: ABC

Base class for MS2 sampler

MZMLChromatogramSampler

MZMLChromatogramSampler(mzml_file_name, roi_params=None)

Bases: ChromatogramSampler

A sampler to return chromatograms extracted from an existing mzML file. Useful to mimic the characteristics of actual experimental data.

Create an MZMLChromatogramSampler object. Args: mzml_file_name: the input mzML file. roi_params: parameters for ROI building, as defined in vimms.Roi.RoiBuilderParams.

Source code in vimms/ChemicalSamplers.py
479
480
481
482
483
484
485
486
487
488
489
490
491
def __init__(self, mzml_file_name, roi_params=None):
    """
    Create an MZMLChromatogramSampler object.
    Args:
        mzml_file_name: the input mzML file.
        roi_params: parameters for ROI building, as defined in [vimms.Roi.RoiBuilderParams][].
    """
    self.mzml_file_name = mzml_file_name
    self.roi_params = roi_params
    if self.roi_params is None:
        self.roi_params = RoiBuilderParams()

    self.good_rois = self._extract_rois()

sample

sample(formula, rt, intensity)

Sample an empirical chromatogram extracted from the mzML file Args: formula: formula, unused rt: RT, unused intensity: intensity, unused

Returns: a [vimms.Chromatograms.EmpiricalChromatogram] object.

Source code in vimms/ChemicalSamplers.py
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
def sample(self, formula, rt, intensity):
    """
    Sample an empirical chromatogram extracted from the mzML file
    Args:
        formula: formula, unused
        rt: RT, unused
        intensity: intensity, unused

    Returns: a [vimms.Chromatograms.EmpiricalChromatogram] object.

    """
    roi_idx = np.random.choice(len(self.good_rois))
    r = self.good_rois[roi_idx]
    chromatogram = EmpiricalChromatogram(
        np.array(r.rt_list),
        np.array(r.mz_list),
        np.array(r.intensity_list),
        single_point_length=0.9,
    )
    return chromatogram

MZMLFormulaSampler

MZMLFormulaSampler(
    mzml_file_name, min_mz=MIN_MZ, max_mz=MAX_MZ, source_polarity=POSITIVE
)

Bases: FormulaSampler

A sampler to generate m/z values from a histogram of m/z taken from a user supplied mzML file

Create an mzML formula sampler Args: mzml_file_name: the source mzML file min_mz: the minimum m/z to consider max_mz: the maximum m/z to consider source_polarity: either POSITIVE or NEGATIVE

Source code in vimms/ChemicalSamplers.py
195
196
197
198
199
200
201
202
203
204
205
206
207
def __init__(self, mzml_file_name, min_mz=MIN_MZ, max_mz=MAX_MZ, source_polarity=POSITIVE):
    """
    Create an mzML formula sampler
    Args:
        mzml_file_name: the source mzML file
        min_mz: the minimum m/z to consider
        max_mz: the maximum m/z to consider
        source_polarity: either POSITIVE or NEGATIVE
    """
    super().__init__(min_mz=min_mz, max_mz=max_mz)
    self.mzml_file_name = mzml_file_name
    self.source_polarity = source_polarity
    self._get_distributions()

sample

sample(n_formulas)

Sample up to n_formulas from the m/z values in the mzML file Args: n_formulas: the number of formula to sample

Returns: a list of Formula objects

Source code in vimms/ChemicalSamplers.py
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
def sample(self, n_formulas):
    """
    Sample up to n_formulas from the m/z values in the mzML file
    Args:
        n_formulas: the number of formula to sample

    Returns: a list of Formula objects

    """
    mz_list = []
    for i in range(n_formulas):
        mz_bin_idx = np.random.choice(len(self.mz_bins), p=self.mz_probs)
        mz_bin = self.mz_bins[mz_bin_idx]
        mz = np.random.rand() * (mz_bin[1] - mz_bin[0]) + mz_bin[0]
        mz_list.append(mz)
    return [(DummyFormula(m), None) for m in mz_list]

MZMLMS2Sampler

MZMLMS2Sampler(
    mzml_file,
    min_n_peaks=1,
    min_total_intensity=1000.0,
    min_proportion=0.1,
    max_proportion=0.8,
    with_replacement=False,
)

Bases: MS2Sampler

A sampler that sample MS2 spectra from an actual mzML file.

Create an MZMLMS2Sampler object Args: mzml_file: the source mzML file min_n_peaks: the minimum number of peaks to consider for each frag. spectra min_total_intensity: the minimum total intensity min_proportion: the minimum proportion to consider max_proportion: the maximum proportion to consider with_replacement: whether to sample with replacement or not

Source code in vimms/ChemicalSamplers.py
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
def __init__(
    self,
    mzml_file,
    min_n_peaks=1,
    min_total_intensity=1e3,
    min_proportion=0.1,
    max_proportion=0.8,
    with_replacement=False,
):
    """
    Create an MZMLMS2Sampler object
    Args:
        mzml_file: the source mzML file
        min_n_peaks: the minimum number of peaks to consider for each frag. spectra
        min_total_intensity: the minimum total intensity
        min_proportion: the minimum proportion to consider
        max_proportion: the maximum proportion to consider
        with_replacement: whether to sample with replacement or not
    """
    self.mzml_file_name = mzml_file
    self.mzml_object = MZMLFile(str(mzml_file))
    self.min_n_peaks = min_n_peaks
    self.min_total_intensity = min_total_intensity
    self.with_replacement = with_replacement

    self.min_proportion = min_proportion
    self.max_proportion = max_proportion

    # only keep MS2 scans that have a least min_n_peaks and
    # a total intesity of at least min_total_intesity
    self._filter_scans()

sample

sample(chemical)

Sample MS2 spectra using chemical as the parent Args: chemical: the parent chemical

Returns: a tuple of (mz_list, intensity_list, parent_proportion)

Source code in vimms/ChemicalSamplers.py
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
def sample(self, chemical):
    """
    Sample MS2 spectra using chemical as the parent
    Args:
        chemical: the parent chemical

    Returns: a tuple of (mz_list, intensity_list, parent_proportion)

    """

    assert len(self.ms2_scans) > 0, (
        "MS2 sampler ran out of scans. "
        "Consider an alternative, or "
        "setting with_replacement to True"
    )
    # pick a scan and removoe
    scan_idx = np.random.choice(len(self.ms2_scans), 1)[0]
    scan = self.ms2_scans[scan_idx]
    if not self.with_replacement:
        del self.ms2_scans[scan_idx]

    parent_proportion = (
        np.random.rand() * (self.max_proportion - self.min_proportion) + self.min_proportion
    )

    mz_list, intensity_list = zip(*scan.peaks)

    return mz_list, intensity_list, parent_proportion

MZMLRTandIntensitySampler

MZMLRTandIntensitySampler(
    mzml_file_name,
    n_intensity_bins=10,
    min_rt=0,
    max_rt=1600,
    min_log_intensity=np.log(10000.0),
    max_log_intensity=np.log(10000000.0),
    roi_params=None,
)

Bases: RTAndIntensitySampler

A sampler to sample RT and intensity values from an existing mzML file. Useful to mimic the characteristics of actual experimental data.

Create an instance of MZMLRTandIntensitySampler. Args: mzml_file_name: the source mzML filename n_intensity_bins: number of bins for intensities min_rt: the minimum RT to consider max_rt: the maximum RT to consider min_log_intensity: the minimum intensity (in log) to consider max_log_intensity: the maximum intensity (in log) to consider roi_params: parameters for ROI building, as defined in vimms.Roi.RoiBuilderParams.

Source code in vimms/ChemicalSamplers.py
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
def __init__(
    self,
    mzml_file_name,
    n_intensity_bins=10,
    min_rt=0,
    max_rt=1600,
    min_log_intensity=np.log(1e4),
    max_log_intensity=np.log(1e7),
    roi_params=None,
):
    """
    Create an instance of MZMLRTandIntensitySampler.
    Args:
        mzml_file_name: the source mzML filename
        n_intensity_bins: number of bins for intensities
        min_rt: the minimum RT to consider
        max_rt: the maximum RT to consider
        min_log_intensity: the minimum intensity (in log) to consider
        max_log_intensity: the maximum intensity (in log) to consider
        roi_params: parameters for ROI building, as defined in [vimms.Roi.RoiBuilderParams][].
    """
    self.min_rt = min_rt
    self.max_rt = max_rt
    self.min_log_intensity = min_log_intensity
    self.max_log_intensity = max_log_intensity
    self.mzml_file_name = mzml_file_name
    self.roi_params = roi_params
    self.n_intensity_bins = n_intensity_bins
    if self.roi_params is None:
        self.roi_params = RoiBuilderParams()
    self._get_distributions()

sample

sample(formula)

Sample RT and intensity value from this sampler Args: formula: the chemical formula, unused for now.

Returns: a tuple of (RT, intensity) values.

Source code in vimms/ChemicalSamplers.py
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
def sample(self, formula):
    """
    Sample RT and intensity value from this sampler
    Args:
        formula: the chemical formula, unused for now.

    Returns: a tuple of (RT, intensity) values.

    """
    rt_bin_idx = np.random.choice(len(self.rt_bins), p=self.rt_probs)
    rt_bin = self.rt_bins[rt_bin_idx]
    rt = np.random.rand() * (rt_bin[1] - rt_bin[0]) + rt_bin[0]

    intensity_bin_idx = np.random.choice(len(self.intensity_bins), p=self.intensity_probs)
    intensity_bin = self.intensity_bins[intensity_bin_idx]
    log_intensity = np.random.rand() * (intensity_bin[1] - intensity_bin[0]) + intensity_bin[0]
    return rt, np.exp(log_intensity)

MzMLScanTimeSampler

MzMLScanTimeSampler(mzml_file, num_bins=1)

Bases: ScanTimeSampler

A scan time sampler that obtains its values from an existing MZML file.

Initialises a MZML scan time sampler object.

Parameters:

Name Type Description Default
num_bins

the number of bins to sample scan durations from

1
Source code in vimms/ChemicalSamplers.py
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
def __init__(self, mzml_file, num_bins=1):
    """
    Initialises a MZML scan time sampler object.

    Args:
        num_bins: the number of bins to sample scan durations from
    """

    self.mzml_file = str(mzml_file)
    self.num_bins = num_bins
    self.time_dict, self.bin_edges = self._extract_timing(self.mzml_file)
    self.is_frag_file = self._is_frag_file(self.time_dict)

    if self.is_frag_file and len(self.time_dict[(1, 1)]) == 0:
        # this could sometimes happen if there's not enough MS1 scan
        # followed by another MS1 scan
        default = DEFAULT_SCAN_TIME_DICT[1]
        logger.warning(
            "Not enough MS1 scans to compute (1, 1) scan duration. "
            "The default of %f will be used" % default
        )
        self.time_dict[(1, 1)] = [default]

sample

sample(current_level, next_level, current_rt)

Sample a scan duration given the MS levels of current and next scans. Args: current_level: the MS level of the current scan next_level: the MS level of the next scan current_rt: the current retention time of the current scan

Returns: a sampled scan duration value

Source code in vimms/ChemicalSamplers.py
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
def sample(self, current_level, next_level, current_rt):
    """
    Sample a scan duration given the MS levels of current and next scans.
    Args:
        current_level: the MS level of the current scan
        next_level: the MS level of the next scan
        current_rt: the current retention time of the current scan

    Returns: a sampled scan duration value

    """

    # Determine the appropriate bin based on the current RT
    current_bin = self._find_bin(current_rt, self.bin_edges)

    # sample a scan duration value extracted from the mzML based
    # on the current and next level
    # note: the same value could be selected again by np.random.choice next time
    values = self.time_dict[current_bin][(current_level, next_level)]
    try:
        sampled = np.random.choice(values, replace=False, size=1)
        return sampled[0]
    except ValueError:  # no value to sample, just return the default
        default = DEFAULT_SCAN_TIME_DICT[current_level]
        return default

PickEverythingFormulaSampler

PickEverythingFormulaSampler(database, min_mz=MIN_MZ, max_mz=MAX_MZ)

Bases: DatabaseFormulaSampler

A sampler that returns everything in the database

Initiliases a Pick-Everything formula sampler

Parameters:

Name Type Description Default
database

a list of Formula objects containing chemical formulae from e.g. HMDB

required
min_mz

the minimum m/z value of formulae to sample from

MIN_MZ
max_mz

the maximum m/z value of formulae to sample from

MAX_MZ
Source code in vimms/ChemicalSamplers.py
130
131
132
133
134
135
136
137
138
139
140
141
def __init__(self, database, min_mz=MIN_MZ, max_mz=MAX_MZ):
    """
    Initiliases a Pick-Everything formula sampler

    Args:
        database: a list of Formula objects containing chemical
                  formulae from e.g. HMDB
        min_mz: the minimum m/z value of formulae to sample from
        max_mz: the maximum m/z value of formulae to sample from
    """
    super().__init__(min_mz=min_mz, max_mz=max_mz)
    self.database = database

sample

sample(n_formulas)

Just return everything from the database

Parameters:

Name Type Description Default
n_formulas

ignored?

required

Returns: all formulae from the database

Source code in vimms/ChemicalSamplers.py
143
144
145
146
147
148
149
150
151
152
153
154
155
156
def sample(self, n_formulas):
    """
    Just return everything from the database

    Args:
        n_formulas: ignored?

    Returns: all formulae from the database

    """
    formula_list = [(Formula(x.chemical_formula), x.name) for x in self.database]
    return list(
        filter(lambda x: x[0].mass >= self.min_mz and x[0].mass <= self.max_mz, formula_list)
    )

RTAndIntensitySampler

Bases: ABC

Base class for RT and intensity sampler. Usually used when initialising a formula object.

ScanTimeSampler

Bases: ABC

Base class for scan time sampler

UniformMS2Sampler

UniformMS2Sampler(
    poiss_peak_mean=10,
    min_mz=MIN_MZ_MS2,
    min_proportion=0.1,
    max_proportion=0.8,
)

Bases: MS2Sampler

A sampler that generates MS2 peaks uniformly between min_mz and the mass of the formula.

Initialises uniform MS2 sampler

Parameters:

Name Type Description Default
poiss_peak_mean

the mean of the Poisson distribution used to draw the number of peaks

10
min_mz

minimum m/z value

MIN_MZ_MS2
min_proportion

minimum proportion from the parent MS1 peak intensities

0.1
max_proportion

maximum proportion from the parent MS1 peak intensities

0.8
Source code in vimms/ChemicalSamplers.py
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
def __init__(
    self, poiss_peak_mean=10, min_mz=MIN_MZ_MS2, min_proportion=0.1, max_proportion=0.8
):
    """
    Initialises uniform MS2 sampler

    Args:
        poiss_peak_mean: the mean of the Poisson distribution used
                         to draw the number of peaks
        min_mz: minimum m/z value
        min_proportion: minimum proportion from the parent MS1
                        peak intensities
        max_proportion: maximum proportion from the parent MS1
                        peak intensities
    """
    self.poiss_peak_mean = poiss_peak_mean
    self.min_mz = min_mz
    # proportion of parent intensity shared by MS2
    self.min_proportion = min_proportion
    self.max_proportion = max_proportion

sample

sample(chemical)

Samples n_peaks of MS2 peaks uniformly between min_mz and the exact mass of the formula. The intensity is also randomly sampled between between min_proportion and max_proportion of the parent formula intensity

Parameters:

Name Type Description Default
chemical

the chemical to compute max m/z value from

required

Returns: a tuple of (mz_list, intensity_list, parent_proportion)

Source code in vimms/ChemicalSamplers.py
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
def sample(self, chemical):
    """
    Samples n_peaks of MS2 peaks uniformly between min_mz and
    the exact mass of the formula. The intensity is also randomly sampled
    between between min_proportion and max_proportion of the parent
    formula intensity

    Args:
        chemical: the chemical to compute max m/z value from

    Returns: a tuple of (mz_list, intensity_list, parent_proportion)

    """
    n_peaks = np.random.poisson(self.poiss_peak_mean)
    max_mz = chemical.mass
    mz_list = uniform_list(n_peaks, self.min_mz, max_mz)
    intensity_list = uniform_list(n_peaks, 0, 1)

    s = sum(intensity_list)
    intensity_list = [i / s for i in intensity_list]
    parent_proportion = (
        np.random.rand() * (self.max_proportion - self.min_proportion) + self.min_proportion
    )

    return mz_list, intensity_list, parent_proportion

UniformMZFormulaSampler

UniformMZFormulaSampler(min_mz=MIN_MZ, max_mz=MAX_MZ)

Bases: FormulaSampler

A sampler to generate formula uniformly between min_mz to max_mz, so just mz rather then formulas. Resulting in UnknownChemical objects instead of known_chemical ones.

Source code in vimms/ChemicalSamplers.py
41
42
43
44
45
46
47
48
49
def __init__(self, min_mz=MIN_MZ, max_mz=MAX_MZ):
    """
    Create a Formula sampler
    Args:
        min_mz: the minimum m/z value of formulae to sample from
        max_mz: the maximum m/z value of formulae to sample from
    """
    self.min_mz = min_mz
    self.max_mz = max_mz

sample

sample(n_formulas)

Samples n_formulas uniformly between min_mz and max_mz

Parameters:

Name Type Description Default
n_formulas

the number of formula to draw

required

Returns: a list of Formula objects

Source code in vimms/ChemicalSamplers.py
111
112
113
114
115
116
117
118
119
120
121
122
def sample(self, n_formulas):
    """
    Samples n_formulas uniformly between min_mz and max_mz

    Args:
        n_formulas: the number of formula to draw

    Returns: a list of Formula objects

    """
    mz_list = np.random.rand(n_formulas) * (self.max_mz - self.min_mz) + self.min_mz
    return [(DummyFormula(m), None) for m in mz_list]

UniformRTAndIntensitySampler

UniformRTAndIntensitySampler(
    min_rt=0,
    max_rt=1600,
    min_log_intensity=np.log(10000.0),
    max_log_intensity=np.log(10000000.0),
)

Bases: RTAndIntensitySampler

A sampler to sample RT and log intensity uniformly. See class def for min and max log intensity. Returns actual intensity, but samples in log space.

Initialises uniform RT and intensity sampler

Parameters:

Name Type Description Default
min_rt

minimum RT

0
max_rt

maximum RT

1600
min_log_intensity

minimum log intensity

log(10000.0)
max_log_intensity

maximum log intensity

log(10000000.0)
Source code in vimms/ChemicalSamplers.py
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
def __init__(
    self, min_rt=0, max_rt=1600, min_log_intensity=np.log(1e4), max_log_intensity=np.log(1e7)
):
    """
    Initialises uniform RT and intensity sampler

    Args:
        min_rt: minimum RT
        max_rt: maximum RT
        min_log_intensity: minimum log intensity
        max_log_intensity: maximum log intensity
    """
    self.min_rt = min_rt
    self.max_rt = max_rt
    self.min_log_intensity = min_log_intensity
    self.max_log_intensity = max_log_intensity

sample

sample(formula)

Samples RT and log intensity uniformly between (min_rt, max_rt) and (min_log_intensity, max_log_intensity)

Parameters:

Name Type Description Default
formula

the formula to condition on (can be ignored)

required

Returns: a tuple of (RT, intensity)

Source code in vimms/ChemicalSamplers.py
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
def sample(self, formula):
    """
    Samples RT and log intensity uniformly between (min_rt, max_rt) and
    (min_log_intensity, max_log_intensity)

    Args:
        formula: the formula to condition on (can be ignored)

    Returns: a tuple of (RT, intensity)

    """
    rt = np.random.rand() * (self.max_rt - self.min_rt) + self.min_rt
    diff = self.max_log_intensity - self.min_log_intensity
    log_intensity = np.random.rand() * (diff) + self.min_log_intensity
    return rt, np.exp(log_intensity)