corems.mass_spectra.input.brukerSolarix
1__author__ = "Yuri E. Corilo" 2__date__ = "Oct 29, 2019" 3 4from threading import Thread 5from pathlib import Path 6from s3path import S3Path 7 8# import h5py 9 10# from corems.encapsulation.constant import Labels 11from corems.mass_spectra.factory.lc_class import LCMSBase 12 13# from corems.encapsulation.factory.parameters import default_parameters 14from corems.transient.input.brukerSolarix import ReadBrukerSolarix 15 16 17class ReadBruker_SolarixTransientMassSpectra(Thread): 18 """ 19 Class for reading Bruker Solarix Transient Mass Spectra. 20 21 Parameters 22 ---------- 23 d_directory_location : str, pathlib.Path, or s3path.S3Path 24 Path object from pathlib containing the file location. 25 analyzer : str, optional 26 Type of analyzer used in the mass spectrometer. Defaults to "ICR". 27 instrument_label : str, optional 28 Label for the instrument. Defaults to "15T". 29 auto_process : bool, optional 30 Flag indicating whether to automatically process the mass spectra. Defaults to True. 31 keep_profile : bool, optional 32 Flag indicating whether to keep the profile data in the mass spectra. Defaults to False. 33 """ 34 35 def __init__( 36 self, 37 d_directory_location: str | Path | S3Path, 38 analyzer="ICR", 39 instrument_label="15T", 40 auto_process=True, 41 keep_profile=False, 42 ): 43 Thread.__init__(self) 44 45 if isinstance(d_directory_location, str): 46 # if obj is a string it defaults to create a Path obj, pass the S3Path if needed 47 d_directory_location = Path(d_directory_location) 48 49 if not d_directory_location.exists(): 50 raise FileNotFoundError("File does not exist: " + str(d_directory_location)) 51 52 self.scan_attr = d_directory_location / "scan.xml" 53 54 if not self.scan_attr.exists(): 55 raise FileExistsError( 56 "%s does not seem to be a valid Solarix Mass Spectra Experiment,\ 57 maybe an Imaging experiment?\ 58 please ReadBruker_SolarixTransientImage class for Imaging dataset " 59 % d_directory_location 60 ) 61 62 self.lcms = LCMSBase(d_directory_location, analyzer, instrument_label) 63 64 self.auto_process = auto_process 65 self.keep_profile = keep_profile 66 67 def get_scan_attr(self) -> dict: 68 """ 69 Get the scan attributes from the scan.xml file. 70 71 Returns 72 ------- 73 dict 74 Dictionary containing the scan number as key and a tuple of retention time and TIC as value. 75 """ 76 from bs4 import BeautifulSoup 77 78 soup = BeautifulSoup(self.scan_attr.open(), "xml") 79 80 list_rt = [float(rt.text) for rt in soup.find_all("minutes")] 81 list_tic = [float(tic.text) for tic in soup.find_all("tic")] 82 list_scan = [int(scan.text) for scan in soup.find_all("count")] 83 84 dict_scan_rt_tic = dict(zip(list_scan, zip(list_rt, list_tic))) 85 86 return dict_scan_rt_tic 87 88 def import_mass_spectra(self) -> None: 89 """ 90 Import the mass spectra from the scan.xml file. 91 """ 92 dict_scan_rt_tic = self.get_scan_attr() 93 94 list_rt, list_tic = ( 95 list(), 96 list(), 97 ) 98 99 list_scans = sorted(list(dict_scan_rt_tic.keys())) 100 101 for scan_number in list_scans: 102 mass_spec = self.get_mass_spectrum(scan_number) 103 104 self.lcms.add_mass_spectrum(mass_spec) 105 106 list_rt.append(dict_scan_rt_tic.get(scan_number)[0]) 107 108 list_tic.append(dict_scan_rt_tic.get(scan_number)[1]) 109 110 self.lcms.retention_time = list_rt 111 self.lcms.tic = list_tic 112 self.lcms.scans_number = list_scans 113 114 def get_mass_spectrum(self, scan_number: int): 115 """ 116 Get the mass spectrum for a given scan number. 117 118 Parameters 119 ---------- 120 scan_number : int 121 Scan number. 122 123 """ 124 bruker_reader = ReadBrukerSolarix(self.lcms.file_location) 125 126 bruker_transient = bruker_reader.get_transient(scan_number) 127 128 mass_spec = bruker_transient.get_mass_spectrum( 129 plot_result=False, 130 auto_process=self.auto_process, 131 keep_profile=self.keep_profile, 132 ) 133 134 return mass_spec 135 136 def run(self): 137 """ 138 Run the import_mass_spectra method. 139 """ 140 self.import_mass_spectra() 141 142 def get_lcms_obj(self): 143 """ 144 Get the LCMSBase object. 145 146 Raises 147 ------ 148 Exception 149 If the LCMSBase object is empty. 150 """ 151 if self.lcms: 152 return self.lcms 153 else: 154 raise Exception("Returning an empty LCMSBase class.")
18class ReadBruker_SolarixTransientMassSpectra(Thread): 19 """ 20 Class for reading Bruker Solarix Transient Mass Spectra. 21 22 Parameters 23 ---------- 24 d_directory_location : str, pathlib.Path, or s3path.S3Path 25 Path object from pathlib containing the file location. 26 analyzer : str, optional 27 Type of analyzer used in the mass spectrometer. Defaults to "ICR". 28 instrument_label : str, optional 29 Label for the instrument. Defaults to "15T". 30 auto_process : bool, optional 31 Flag indicating whether to automatically process the mass spectra. Defaults to True. 32 keep_profile : bool, optional 33 Flag indicating whether to keep the profile data in the mass spectra. Defaults to False. 34 """ 35 36 def __init__( 37 self, 38 d_directory_location: str | Path | S3Path, 39 analyzer="ICR", 40 instrument_label="15T", 41 auto_process=True, 42 keep_profile=False, 43 ): 44 Thread.__init__(self) 45 46 if isinstance(d_directory_location, str): 47 # if obj is a string it defaults to create a Path obj, pass the S3Path if needed 48 d_directory_location = Path(d_directory_location) 49 50 if not d_directory_location.exists(): 51 raise FileNotFoundError("File does not exist: " + str(d_directory_location)) 52 53 self.scan_attr = d_directory_location / "scan.xml" 54 55 if not self.scan_attr.exists(): 56 raise FileExistsError( 57 "%s does not seem to be a valid Solarix Mass Spectra Experiment,\ 58 maybe an Imaging experiment?\ 59 please ReadBruker_SolarixTransientImage class for Imaging dataset " 60 % d_directory_location 61 ) 62 63 self.lcms = LCMSBase(d_directory_location, analyzer, instrument_label) 64 65 self.auto_process = auto_process 66 self.keep_profile = keep_profile 67 68 def get_scan_attr(self) -> dict: 69 """ 70 Get the scan attributes from the scan.xml file. 71 72 Returns 73 ------- 74 dict 75 Dictionary containing the scan number as key and a tuple of retention time and TIC as value. 76 """ 77 from bs4 import BeautifulSoup 78 79 soup = BeautifulSoup(self.scan_attr.open(), "xml") 80 81 list_rt = [float(rt.text) for rt in soup.find_all("minutes")] 82 list_tic = [float(tic.text) for tic in soup.find_all("tic")] 83 list_scan = [int(scan.text) for scan in soup.find_all("count")] 84 85 dict_scan_rt_tic = dict(zip(list_scan, zip(list_rt, list_tic))) 86 87 return dict_scan_rt_tic 88 89 def import_mass_spectra(self) -> None: 90 """ 91 Import the mass spectra from the scan.xml file. 92 """ 93 dict_scan_rt_tic = self.get_scan_attr() 94 95 list_rt, list_tic = ( 96 list(), 97 list(), 98 ) 99 100 list_scans = sorted(list(dict_scan_rt_tic.keys())) 101 102 for scan_number in list_scans: 103 mass_spec = self.get_mass_spectrum(scan_number) 104 105 self.lcms.add_mass_spectrum(mass_spec) 106 107 list_rt.append(dict_scan_rt_tic.get(scan_number)[0]) 108 109 list_tic.append(dict_scan_rt_tic.get(scan_number)[1]) 110 111 self.lcms.retention_time = list_rt 112 self.lcms.tic = list_tic 113 self.lcms.scans_number = list_scans 114 115 def get_mass_spectrum(self, scan_number: int): 116 """ 117 Get the mass spectrum for a given scan number. 118 119 Parameters 120 ---------- 121 scan_number : int 122 Scan number. 123 124 """ 125 bruker_reader = ReadBrukerSolarix(self.lcms.file_location) 126 127 bruker_transient = bruker_reader.get_transient(scan_number) 128 129 mass_spec = bruker_transient.get_mass_spectrum( 130 plot_result=False, 131 auto_process=self.auto_process, 132 keep_profile=self.keep_profile, 133 ) 134 135 return mass_spec 136 137 def run(self): 138 """ 139 Run the import_mass_spectra method. 140 """ 141 self.import_mass_spectra() 142 143 def get_lcms_obj(self): 144 """ 145 Get the LCMSBase object. 146 147 Raises 148 ------ 149 Exception 150 If the LCMSBase object is empty. 151 """ 152 if self.lcms: 153 return self.lcms 154 else: 155 raise Exception("Returning an empty LCMSBase class.")
Class for reading Bruker Solarix Transient Mass Spectra.
Parameters
- d_directory_location (str, pathlib.Path, or s3path.S3Path): Path object from pathlib containing the file location.
- analyzer (str, optional): Type of analyzer used in the mass spectrometer. Defaults to "ICR".
- instrument_label (str, optional): Label for the instrument. Defaults to "15T".
- auto_process (bool, optional): Flag indicating whether to automatically process the mass spectra. Defaults to True.
- keep_profile (bool, optional): Flag indicating whether to keep the profile data in the mass spectra. Defaults to False.
36 def __init__( 37 self, 38 d_directory_location: str | Path | S3Path, 39 analyzer="ICR", 40 instrument_label="15T", 41 auto_process=True, 42 keep_profile=False, 43 ): 44 Thread.__init__(self) 45 46 if isinstance(d_directory_location, str): 47 # if obj is a string it defaults to create a Path obj, pass the S3Path if needed 48 d_directory_location = Path(d_directory_location) 49 50 if not d_directory_location.exists(): 51 raise FileNotFoundError("File does not exist: " + str(d_directory_location)) 52 53 self.scan_attr = d_directory_location / "scan.xml" 54 55 if not self.scan_attr.exists(): 56 raise FileExistsError( 57 "%s does not seem to be a valid Solarix Mass Spectra Experiment,\ 58 maybe an Imaging experiment?\ 59 please ReadBruker_SolarixTransientImage class for Imaging dataset " 60 % d_directory_location 61 ) 62 63 self.lcms = LCMSBase(d_directory_location, analyzer, instrument_label) 64 65 self.auto_process = auto_process 66 self.keep_profile = keep_profile
This constructor should always be called with keyword arguments. Arguments are:
group should be None; reserved for future extension when a ThreadGroup class is implemented.
target is the callable object to be invoked by the run() method. Defaults to None, meaning nothing is called.
name is the thread name. By default, a unique name is constructed of the form "Thread-N" where N is a small decimal number.
args is the argument tuple for the target invocation. Defaults to ().
kwargs is a dictionary of keyword arguments for the target invocation. Defaults to {}.
If a subclass overrides the constructor, it must make sure to invoke the base class constructor (Thread.__init__()) before doing anything else to the thread.
68 def get_scan_attr(self) -> dict: 69 """ 70 Get the scan attributes from the scan.xml file. 71 72 Returns 73 ------- 74 dict 75 Dictionary containing the scan number as key and a tuple of retention time and TIC as value. 76 """ 77 from bs4 import BeautifulSoup 78 79 soup = BeautifulSoup(self.scan_attr.open(), "xml") 80 81 list_rt = [float(rt.text) for rt in soup.find_all("minutes")] 82 list_tic = [float(tic.text) for tic in soup.find_all("tic")] 83 list_scan = [int(scan.text) for scan in soup.find_all("count")] 84 85 dict_scan_rt_tic = dict(zip(list_scan, zip(list_rt, list_tic))) 86 87 return dict_scan_rt_tic
Get the scan attributes from the scan.xml file.
Returns
- dict: Dictionary containing the scan number as key and a tuple of retention time and TIC as value.
89 def import_mass_spectra(self) -> None: 90 """ 91 Import the mass spectra from the scan.xml file. 92 """ 93 dict_scan_rt_tic = self.get_scan_attr() 94 95 list_rt, list_tic = ( 96 list(), 97 list(), 98 ) 99 100 list_scans = sorted(list(dict_scan_rt_tic.keys())) 101 102 for scan_number in list_scans: 103 mass_spec = self.get_mass_spectrum(scan_number) 104 105 self.lcms.add_mass_spectrum(mass_spec) 106 107 list_rt.append(dict_scan_rt_tic.get(scan_number)[0]) 108 109 list_tic.append(dict_scan_rt_tic.get(scan_number)[1]) 110 111 self.lcms.retention_time = list_rt 112 self.lcms.tic = list_tic 113 self.lcms.scans_number = list_scans
Import the mass spectra from the scan.xml file.
115 def get_mass_spectrum(self, scan_number: int): 116 """ 117 Get the mass spectrum for a given scan number. 118 119 Parameters 120 ---------- 121 scan_number : int 122 Scan number. 123 124 """ 125 bruker_reader = ReadBrukerSolarix(self.lcms.file_location) 126 127 bruker_transient = bruker_reader.get_transient(scan_number) 128 129 mass_spec = bruker_transient.get_mass_spectrum( 130 plot_result=False, 131 auto_process=self.auto_process, 132 keep_profile=self.keep_profile, 133 ) 134 135 return mass_spec
Get the mass spectrum for a given scan number.
Parameters
- scan_number (int): Scan number.
137 def run(self): 138 """ 139 Run the import_mass_spectra method. 140 """ 141 self.import_mass_spectra()
Run the import_mass_spectra method.
143 def get_lcms_obj(self): 144 """ 145 Get the LCMSBase object. 146 147 Raises 148 ------ 149 Exception 150 If the LCMSBase object is empty. 151 """ 152 if self.lcms: 153 return self.lcms 154 else: 155 raise Exception("Returning an empty LCMSBase class.")
Get the LCMSBase object.
Raises
- Exception: If the LCMSBase object is empty.
Inherited Members
- threading.Thread
- start
- join
- name
- ident
- is_alive
- daemon
- isDaemon
- setDaemon
- getName
- setName
- native_id