corems.mass_spectrum.input.boosterHDF5

  1from io import BytesIO
  2
  3import h5py
  4from s3path import S3Path
  5
  6from corems.encapsulation.constant import Labels
  7from corems.encapsulation.factory.parameters import default_parameters
  8from corems.mass_spectrum.factory.MassSpectrumClasses import MassSpecProfile
  9from corems.mass_spectrum.input.baseClass import MassListBaseClass
 10
 11
 12class ReadHDF_BoosterMassSpectrum(MassListBaseClass):
 13    """The ReadHDF_BoosterMassSpectrum class parses the mass spectrum data from an HDF file and generate a mass spectrum object.
 14
 15    Parameters
 16    ----------
 17    file_location : str
 18        The path to the HDF file.
 19    isCentroid : bool, optional
 20        Specifies whether the mass spectrum is centroided or not. Default is False.
 21
 22    Attributes
 23    ----------
 24    polarity : int
 25        The polarity of the mass spectrum.
 26    h5pydata : h5py.File
 27        The HDF file object.
 28    scans : list
 29        The list of scan names in the HDF file.
 30
 31    Methods
 32    -------
 33    * get_data_profile(mz, abundance, auto_process). Returns a MassSpecProfile object from the given m/z and abundance arrays.
 34    * get_attr_data(scan, attr_srt). Returns the attribute value for the given scan and attribute name.
 35    * get_polarity(file_location). Returns the polarity of the mass spectrum.
 36    * get_mass_spectrum(auto_process). Returns the mass spectrum as a MassSpecProfile object.
 37    * get_output_parameters(). Returns the default output parameters for the mass spectrum.
 38    """
 39
 40    def __init__(self, file_location, isCentroid=False):
 41        self.polarity = self.get_polarity(file_location)
 42        super().__init__(file_location, isCentroid=False)
 43
 44    def get_data_profile(self, mz, abundance, auto_process) -> MassSpecProfile:
 45        """
 46        Returns a MassSpecProfile object from the given m/z and abundance arrays.
 47
 48        Parameters
 49        ----------
 50        mz : array_like
 51            The m/z values.
 52        abundance : array_like
 53            The abundance values.
 54        auto_process : bool
 55            Specifies whether to automatically process the mass spectrum.
 56
 57        Returns
 58        -------
 59        MassSpecProfile
 60            The MassSpecProfile object.
 61
 62        """
 63        data_dict = {Labels.mz: mz, Labels.abundance: abundance}
 64        output_parameters = self.get_output_parameters()
 65        return MassSpecProfile(data_dict, output_parameters, auto_process=auto_process)
 66
 67    def get_attr_data(self, scan, attr_srt):
 68        """
 69        Returns the attribute value for the given scan and attribute name.
 70
 71        Parameters
 72        ----------
 73        scan : int
 74            The scan index.
 75        attr_srt : str
 76            The attribute name.
 77
 78        Returns
 79        -------
 80        object
 81            The attribute value.
 82
 83        """
 84        return self.h5pydata[self.scans[scan]].attrs[attr_srt]
 85
 86    def get_polarity(self, file_location: str | S3Path) -> int:
 87        """
 88        Returns the polarity of the mass spectrum.
 89
 90        Parameters
 91        ----------
 92        file_location : str
 93            The path to the HDF file.
 94
 95        Returns
 96        -------
 97        int
 98            The polarity of the mass spectrum.
 99
100        """
101        if isinstance(file_location, S3Path):
102            data = BytesIO(file_location.open("rb").read())
103        else:
104            data = file_location
105
106        self.h5pydata = h5py.File(data, "r")
107        self.scans = list(self.h5pydata.keys())
108
109        polarity = self.get_attr_data(0, "r_h_polarity")
110
111        if polarity == "negative scan":
112            return -1
113        else:
114            return +1
115
116    def get_mass_spectrum(self, auto_process=True) -> MassSpecProfile:
117        """
118        Returns the mass spectrum as a MassSpecProfile object.
119
120        Parameters
121        ----------
122        auto_process : bool, optional
123            Specifies whether to automatically process the mass spectrum. Default is True.
124
125        Returns
126        -------
127        MassSpecProfile
128            The MassSpecProfile object.
129
130        """
131        if len(self.scans) == 1:
132            booster_data = self.h5pydata[self.scans[0]]
133
134            if self.isCentroid:
135                raise NotImplementedError
136            else:
137                mz = booster_data[0]
138                abun = booster_data[1]
139                return self.get_data_profile(mz, abun, auto_process)
140
141    def get_output_parameters(self) -> dict:
142        """
143        Returns the default output parameters for the mass spectrum.
144
145        Returns
146        -------
147        dict
148            The default output parameters.
149
150        """
151        d_params = default_parameters(self.file_location)
152        d_params["polarity"] = self.polarity
153        d_params["filename_path"] = self.file_location
154        d_params["mobility_scan"] = 0
155        d_params["mobility_rt"] = 0
156        d_params["scan_number"] = 0
157        d_params["rt"] = self.get_attr_data(0, "r_h_start_time")
158        d_params["label"] = Labels.booster_profile
159        d_params["Aterm"] = self.get_attr_data(0, "r_cparams")[0]
160        d_params["Bterm"] = self.get_attr_data(0, "r_cparams")[1]
161        return d_params
class ReadHDF_BoosterMassSpectrum(corems.mass_spectrum.input.baseClass.MassListBaseClass):
 13class ReadHDF_BoosterMassSpectrum(MassListBaseClass):
 14    """The ReadHDF_BoosterMassSpectrum class parses the mass spectrum data from an HDF file and generate a mass spectrum object.
 15
 16    Parameters
 17    ----------
 18    file_location : str
 19        The path to the HDF file.
 20    isCentroid : bool, optional
 21        Specifies whether the mass spectrum is centroided or not. Default is False.
 22
 23    Attributes
 24    ----------
 25    polarity : int
 26        The polarity of the mass spectrum.
 27    h5pydata : h5py.File
 28        The HDF file object.
 29    scans : list
 30        The list of scan names in the HDF file.
 31
 32    Methods
 33    -------
 34    * get_data_profile(mz, abundance, auto_process). Returns a MassSpecProfile object from the given m/z and abundance arrays.
 35    * get_attr_data(scan, attr_srt). Returns the attribute value for the given scan and attribute name.
 36    * get_polarity(file_location). Returns the polarity of the mass spectrum.
 37    * get_mass_spectrum(auto_process). Returns the mass spectrum as a MassSpecProfile object.
 38    * get_output_parameters(). Returns the default output parameters for the mass spectrum.
 39    """
 40
 41    def __init__(self, file_location, isCentroid=False):
 42        self.polarity = self.get_polarity(file_location)
 43        super().__init__(file_location, isCentroid=False)
 44
 45    def get_data_profile(self, mz, abundance, auto_process) -> MassSpecProfile:
 46        """
 47        Returns a MassSpecProfile object from the given m/z and abundance arrays.
 48
 49        Parameters
 50        ----------
 51        mz : array_like
 52            The m/z values.
 53        abundance : array_like
 54            The abundance values.
 55        auto_process : bool
 56            Specifies whether to automatically process the mass spectrum.
 57
 58        Returns
 59        -------
 60        MassSpecProfile
 61            The MassSpecProfile object.
 62
 63        """
 64        data_dict = {Labels.mz: mz, Labels.abundance: abundance}
 65        output_parameters = self.get_output_parameters()
 66        return MassSpecProfile(data_dict, output_parameters, auto_process=auto_process)
 67
 68    def get_attr_data(self, scan, attr_srt):
 69        """
 70        Returns the attribute value for the given scan and attribute name.
 71
 72        Parameters
 73        ----------
 74        scan : int
 75            The scan index.
 76        attr_srt : str
 77            The attribute name.
 78
 79        Returns
 80        -------
 81        object
 82            The attribute value.
 83
 84        """
 85        return self.h5pydata[self.scans[scan]].attrs[attr_srt]
 86
 87    def get_polarity(self, file_location: str | S3Path) -> int:
 88        """
 89        Returns the polarity of the mass spectrum.
 90
 91        Parameters
 92        ----------
 93        file_location : str
 94            The path to the HDF file.
 95
 96        Returns
 97        -------
 98        int
 99            The polarity of the mass spectrum.
100
101        """
102        if isinstance(file_location, S3Path):
103            data = BytesIO(file_location.open("rb").read())
104        else:
105            data = file_location
106
107        self.h5pydata = h5py.File(data, "r")
108        self.scans = list(self.h5pydata.keys())
109
110        polarity = self.get_attr_data(0, "r_h_polarity")
111
112        if polarity == "negative scan":
113            return -1
114        else:
115            return +1
116
117    def get_mass_spectrum(self, auto_process=True) -> MassSpecProfile:
118        """
119        Returns the mass spectrum as a MassSpecProfile object.
120
121        Parameters
122        ----------
123        auto_process : bool, optional
124            Specifies whether to automatically process the mass spectrum. Default is True.
125
126        Returns
127        -------
128        MassSpecProfile
129            The MassSpecProfile object.
130
131        """
132        if len(self.scans) == 1:
133            booster_data = self.h5pydata[self.scans[0]]
134
135            if self.isCentroid:
136                raise NotImplementedError
137            else:
138                mz = booster_data[0]
139                abun = booster_data[1]
140                return self.get_data_profile(mz, abun, auto_process)
141
142    def get_output_parameters(self) -> dict:
143        """
144        Returns the default output parameters for the mass spectrum.
145
146        Returns
147        -------
148        dict
149            The default output parameters.
150
151        """
152        d_params = default_parameters(self.file_location)
153        d_params["polarity"] = self.polarity
154        d_params["filename_path"] = self.file_location
155        d_params["mobility_scan"] = 0
156        d_params["mobility_rt"] = 0
157        d_params["scan_number"] = 0
158        d_params["rt"] = self.get_attr_data(0, "r_h_start_time")
159        d_params["label"] = Labels.booster_profile
160        d_params["Aterm"] = self.get_attr_data(0, "r_cparams")[0]
161        d_params["Bterm"] = self.get_attr_data(0, "r_cparams")[1]
162        return d_params

The ReadHDF_BoosterMassSpectrum class parses the mass spectrum data from an HDF file and generate a mass spectrum object.

Parameters
  • file_location (str): The path to the HDF file.
  • isCentroid (bool, optional): Specifies whether the mass spectrum is centroided or not. Default is False.
Attributes
  • polarity (int): The polarity of the mass spectrum.
  • h5pydata (h5py.File): The HDF file object.
  • scans (list): The list of scan names in the HDF file.
Methods
  • get_data_profile(mz, abundance, auto_process). Returns a MassSpecProfile object from the given m/z and abundance arrays.
  • get_attr_data(scan, attr_srt). Returns the attribute value for the given scan and attribute name.
  • get_polarity(file_location). Returns the polarity of the mass spectrum.
  • get_mass_spectrum(auto_process). Returns the mass spectrum as a MassSpecProfile object.
  • get_output_parameters(). Returns the default output parameters for the mass spectrum.
ReadHDF_BoosterMassSpectrum(file_location, isCentroid=False)
41    def __init__(self, file_location, isCentroid=False):
42        self.polarity = self.get_polarity(file_location)
43        super().__init__(file_location, isCentroid=False)
polarity
def get_data_profile( self, mz, abundance, auto_process) -> corems.mass_spectrum.factory.MassSpectrumClasses.MassSpecProfile:
45    def get_data_profile(self, mz, abundance, auto_process) -> MassSpecProfile:
46        """
47        Returns a MassSpecProfile object from the given m/z and abundance arrays.
48
49        Parameters
50        ----------
51        mz : array_like
52            The m/z values.
53        abundance : array_like
54            The abundance values.
55        auto_process : bool
56            Specifies whether to automatically process the mass spectrum.
57
58        Returns
59        -------
60        MassSpecProfile
61            The MassSpecProfile object.
62
63        """
64        data_dict = {Labels.mz: mz, Labels.abundance: abundance}
65        output_parameters = self.get_output_parameters()
66        return MassSpecProfile(data_dict, output_parameters, auto_process=auto_process)

Returns a MassSpecProfile object from the given m/z and abundance arrays.

Parameters
  • mz (array_like): The m/z values.
  • abundance (array_like): The abundance values.
  • auto_process (bool): Specifies whether to automatically process the mass spectrum.
Returns
  • MassSpecProfile: The MassSpecProfile object.
def get_attr_data(self, scan, attr_srt):
68    def get_attr_data(self, scan, attr_srt):
69        """
70        Returns the attribute value for the given scan and attribute name.
71
72        Parameters
73        ----------
74        scan : int
75            The scan index.
76        attr_srt : str
77            The attribute name.
78
79        Returns
80        -------
81        object
82            The attribute value.
83
84        """
85        return self.h5pydata[self.scans[scan]].attrs[attr_srt]

Returns the attribute value for the given scan and attribute name.

Parameters
  • scan (int): The scan index.
  • attr_srt (str): The attribute name.
Returns
  • object: The attribute value.
def get_polarity(self, file_location: str | s3path.S3Path) -> int:
 87    def get_polarity(self, file_location: str | S3Path) -> int:
 88        """
 89        Returns the polarity of the mass spectrum.
 90
 91        Parameters
 92        ----------
 93        file_location : str
 94            The path to the HDF file.
 95
 96        Returns
 97        -------
 98        int
 99            The polarity of the mass spectrum.
100
101        """
102        if isinstance(file_location, S3Path):
103            data = BytesIO(file_location.open("rb").read())
104        else:
105            data = file_location
106
107        self.h5pydata = h5py.File(data, "r")
108        self.scans = list(self.h5pydata.keys())
109
110        polarity = self.get_attr_data(0, "r_h_polarity")
111
112        if polarity == "negative scan":
113            return -1
114        else:
115            return +1

Returns the polarity of the mass spectrum.

Parameters
  • file_location (str): The path to the HDF file.
Returns
  • int: The polarity of the mass spectrum.
def get_mass_spectrum( self, auto_process=True) -> corems.mass_spectrum.factory.MassSpectrumClasses.MassSpecProfile:
117    def get_mass_spectrum(self, auto_process=True) -> MassSpecProfile:
118        """
119        Returns the mass spectrum as a MassSpecProfile object.
120
121        Parameters
122        ----------
123        auto_process : bool, optional
124            Specifies whether to automatically process the mass spectrum. Default is True.
125
126        Returns
127        -------
128        MassSpecProfile
129            The MassSpecProfile object.
130
131        """
132        if len(self.scans) == 1:
133            booster_data = self.h5pydata[self.scans[0]]
134
135            if self.isCentroid:
136                raise NotImplementedError
137            else:
138                mz = booster_data[0]
139                abun = booster_data[1]
140                return self.get_data_profile(mz, abun, auto_process)

Returns the mass spectrum as a MassSpecProfile object.

Parameters
  • auto_process (bool, optional): Specifies whether to automatically process the mass spectrum. Default is True.
Returns
  • MassSpecProfile: The MassSpecProfile object.
def get_output_parameters(self) -> dict:
142    def get_output_parameters(self) -> dict:
143        """
144        Returns the default output parameters for the mass spectrum.
145
146        Returns
147        -------
148        dict
149            The default output parameters.
150
151        """
152        d_params = default_parameters(self.file_location)
153        d_params["polarity"] = self.polarity
154        d_params["filename_path"] = self.file_location
155        d_params["mobility_scan"] = 0
156        d_params["mobility_rt"] = 0
157        d_params["scan_number"] = 0
158        d_params["rt"] = self.get_attr_data(0, "r_h_start_time")
159        d_params["label"] = Labels.booster_profile
160        d_params["Aterm"] = self.get_attr_data(0, "r_cparams")[0]
161        d_params["Bterm"] = self.get_attr_data(0, "r_cparams")[1]
162        return d_params

Returns the default output parameters for the mass spectrum.

Returns
  • dict: The default output parameters.