pept.tracking.TimeOfFlight#
- class pept.tracking.TimeOfFlight(temporal_resolution=None, points=False)[source]#
Bases:
LineDataFilter
Compute the positron annihilation locations of each LoR as given by the Time Of Flight (ToF) data of the two LoR timestamps.
Filter signature:
LineData -> TimeOfFlight.fit_sample -> LineData (points = False) LineData -> TimeOfFlight.fit_sample -> PointData (points = True)
Importantly, the input LineData must have at least 8 columns, formatted as [t1, x1, y1, z1, x2, y2, z2, t2] - notice the different timestamps of the two LoR ends.
If points = False (default), the computed ToF points are saved as an extra attribute “tof” in the input LineData; otherwise they are returned directly.
The temporal_resolution should be set to the FWHM of the temporal resolution in the LoR timestamps, in self-consistent units of measure (e.g. m / s or mm / ms, but not mm / s). If it is set, the “temporal_resolution” and “spatial_resolution” extra attributes are set on the ToF points.
New in pept-0.4.2
Examples
Generate 10 random LoRs between (-100, 100) mm and ms with 8 columns for the ToF format.
>>> import numpy as np >>> import pept
>>> rng = np.random.default_rng(0) >>> lors = pept.LineData( >>> rng.uniform(-100, 100, (10, 8)), >>> columns = ["t1", "x1", "y1", "z1", "x2", "y2", "z2", "t2"], >>> ) >>> lors pept.LineData (samples: 1) -------------------------- sample_size = 10 overlap = 0 lines = (rows: 10, columns: 8) [[ 57.4196615 -52.1261114 ... -9.93212667 59.26485406] [-53.8715582 -89.59573979 ... -40.26077344 34.39897559] ... [ 51.59020047 2.55174465 ... -31.13800424 -13.94025361] [ 93.21241616 12.44636845 ... -75.08905883 -42.3338486 ]] columns = ['t1', 'x1', 'y1', 'z1', 'x2', 'y2', 'z2', 't2'] attrs = {}
Compute Time of Flight annihilation locations from the two timestamps in the data above. Assume a temporal resolution of 100 ps - be careful to use self-consistent units; in this case we are using mm and ms:
>>> from pept.tracking import *
>>> temporal_resolution = 100e-12 * 1000 # ms >>> lors_tof = TimeOfFlight(temporal_resolution).fit_sample(lors) >>> lors_tof pept.LineData (samples: 1) -------------------------- sample_size = 10 overlap = 0 lines = (rows: 10, columns: 8) [[ 57.4196615 -52.1261114 ... -9.93212667 59.26485406] [-53.8715582 -89.59573979 ... -40.26077344 34.39897559] ... [ 51.59020047 2.55174465 ... -31.13800424 -13.94025361] [ 93.21241616 12.44636845 ... -75.08905883 -42.3338486 ]] columns = ['t1', 'x1', 'y1', 'z1', 'x2', 'y2', 'z2', 't2'] attrs = { 'tof': pept.PointData (samples: 1) --------------------------- sample_... }
>>> lors_tof.attrs["tof"] pept.PointData (samples: 1) --------------------------- sample_size = 10 overlap = 0 points = (rows: 10, columns: 4) [[ 5.64970655e+01 -3.22092074e+07 2.41767704e+08 -1.30428351e+08] [-9.80068250e+01 -2.48775932e+09 -1.12904720e+10 -6.43480969e+09] ... [ 1.88249731e+01 3.34819602e+09 -8.78848458e+09 2.83529405e+09] [ 2.54392837e+01 1.90343279e+10 -1.92717662e+09 -6.84078611e+09]] columns = ['t', 'x', 'y', 'z'] attrs = { 'temporal_resolution': 1.0000000000000001e-07 'spatial_resolution': 29.9792458 }
Alternatively, you can extract only the ToF points directly:
>>> tof = TimeOfFlight(temporal_resolution, points = True).fit_sample(lors) >>> tof pept.PointData (samples: 1) --------------------------- sample_size = 10 overlap = 0 points = (rows: 10, columns: 4) [[ 5.64970655e+01 -3.22092074e+07 2.41767704e+08 -1.30428351e+08] [-9.80068250e+01 -2.48775932e+09 -1.12904720e+10 -6.43480969e+09] ... [ 1.88249731e+01 3.34819602e+09 -8.78848458e+09 2.83529405e+09] [ 2.54392837e+01 1.90343279e+10 -1.92717662e+09 -6.84078611e+09]] columns = ['t', 'x', 'y', 'z'] attrs = { 'temporal_resolution': 1.0000000000000001e-07 'spatial_resolution': 29.9792458 }
Methods
__init__
([temporal_resolution, points])copy
([deep])Create a deep copy of an instance of this class, including all inner attributes.
fit
(line_data[, executor, max_workers, verbose])Apply self.fit_sample (implemented by subclasses) according to the execution policy.
fit_sample
(sample)load
(filepath)Load a saved / pickled PEPTObject object from filepath.
save
(filepath)Save a PEPTObject instance as a binary pickle object.
- copy(deep=True)#
Create a deep copy of an instance of this class, including all inner attributes.
- fit(line_data, executor='joblib', max_workers=None, verbose=True)#
Apply self.fit_sample (implemented by subclasses) according to the execution policy. Simply return a list of processed samples. If you need a reduction step (e.g. stack all processed samples), apply it in the subclass.
- static load(filepath)#
Load a saved / pickled PEPTObject object from filepath.
Most often the full object state was saved using the .save method.
- Parameters
- filepath
filename
orfile
handle
If filepath is a path (rather than file handle), it is relative to where python is called.
- filepath
- Returns
pept.PEPTObject
subclass
instance
The loaded object.
Examples
Save a LineData instance, then load it back:
>>> lines = pept.LineData([[1, 2, 3, 4, 5, 6, 7]]) >>> lines.save("lines.pickle")
>>> lines_reloaded = pept.LineData.load("lines.pickle")
- save(filepath)#
Save a PEPTObject instance as a binary pickle object.
Saves the full object state, including inner attributes, in a portable binary format. Load back the object using the load method.
- Parameters
- filepath
filename
orfile
handle
If filepath is a path (rather than file handle), it is relative to where python is called.
- filepath
Examples
Save a LineData instance, then load it back:
>>> lines = pept.LineData([[1, 2, 3, 4, 5, 6, 7]]) >>> lines.save("lines.pickle")
>>> lines_reloaded = pept.LineData.load("lines.pickle")