analysis¶
analysis.beamline¶
- hummingbird.analysis.beamline.averagePhotonEnergy(evt, records, outkey='averagePhotonEnergy')[source]¶
Averages across given photon energies and adds it to evt[“analysis”][outkey].
- Args:
- evt:
The event variable
- records:
A dictionary of photon energy
Record()
objects
- Kwargs:
- outkey(str):
Data key of resulting
Record()
, default is ‘averagePhotonEnergy’
- Authors:
Benedikt J. Daurer
- hummingbird.analysis.beamline.averagePulseEnergy(evt, records, outkey='averagePulseEnergy')[source]¶
Averages across given pulse energies and adds it to evt[“analysis”][outkey].
- Args:
- evt:
The event variable
- records:
A dictionary of pulse energy
Record()
objects
- Kwargs:
- outkey(str):
Data key of resulting
Record()
, default is ‘averagePulseEnergy’
- Authors:
Filipe Maia, Benedikt J. Daurer
analysis.event¶
- hummingbird.analysis.event.printKeys(evt, group=None)[source]¶
prints available keys of Hummingbird events
analysis.hitfinding¶
- hummingbird.analysis.hitfinding.countHits(evt, hit, outkey='nrHits')[source]¶
Counts hits and adds the total nr. of hits to
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- hit:
A boolean (True for hit, False for miss)
- Kwargs:
- outkey(str):
Data key of resulting
Record()
object, default is “nrHits”
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se), Jonas Sellberg, Tomas Ekeberg
- hummingbird.analysis.hitfinding.countHitscore(evt, hitscore, hitscoreThreshold=200, outkey='predef: ')[source]¶
A simple hitfinder that performs a limit test against an already defined hitscore and adds the result to
evt["analysis"][outkey + "isHit"]
, and the hitscore toevt["analysis"][outkey + "hitscore"]
.- Args:
- evt:
The event variable
- hitscore:
A pre-defined hitscore
- Kwargs:
- hitscoreThreshold(int):
Events with hitscore above this threshold are hits, default=200
- outkey(str):
Prefix of data key of resulting
Record()
object, default is “predef: “
- Authors:
Carl Nettelblad (carl.nettelblad@it.uu.se), Benedikt J. Daurer
- hummingbird.analysis.hitfinding.countLitPixels(evt, record, aduThreshold=20, hitscoreThreshold=200, hitscoreDark=0, hitscoreMax=None, mask=None, stack=False, outkey='litpixel: ')[source]¶
A simple hitfinder that counts the number of lit pixels and adds the result to
evt["analysis"][outkey + "isHit"]
,evt["analysis"][outkey + "isMiss"]
, and the hitscore toevt["analysis"][outkey + "hitscore"]
.- Args:
- evt:
The event variable
- record:
A pixel detector
Record()
object
- Kwargs:
- aduThreshold(int):
only pixels above this threshold (in ADUs) are valid, default=20
- hitscoreThreshold(int):
events with hitscore (Nr. of lit pixels) above this threshold are hits, default=200
- hitscoreMax(int):
events with hitscore (Nr. of lit pixels) below this threshold (if not None) are hits, default=None
- hitscoreDark(int):
events with hitscore (Nr. of lit pixels) above this threshold are not darks (so either hit or miss), default=0
- mask(int, bool):
only use masked pixel (mask == True or 1) for counting the nr. of lit pixels
- outkey(str):
Prefix of data key of resulting
Record()
object, default is “litpixel: “
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.hitfinding.countPhotonsAgainstEnergyFunction(evt, photonscore_record, energy_record, energyFunction=<function <lambda>>, outkey='photons: ')[source]¶
A hitfinder that tests given photon score (e.g. photon count) against a predicted photon threshold that is dependent on some given energy and adds a boolean to
evt["analysis"][outkey + "isHit"]
, the hitscore toevt["analysis"][outkey + "hitscore"]
and the limit to ``evt[“analysis”][outkey + “photonLimit”]- Args:
- evt:
The event variable
- photonscore_record:
A
Record()
object containing a photon score, e.g. total photon count
:energy_record:” A
Record()
object containing an energy value, e.g. from gas monitor detector- Kwargs:
- energyFunction(function with double argument):
Function that computes the photon threshold, given the energy
- outkey(str):
Prefix of data key of resulting
Record()
object, default is “photons: “
- Authors:
Carl Nettelblad (carl.nettelblad@it.uu.se)
- hummingbird.analysis.hitfinding.countPhotonsAgainstEnergyPolynomial(evt, photonscore_record, energy_record, energyPolynomial=[200], outkey='photons: ')[source]¶
A hitfinder that tests photon score (e.g. photon count) against a predicted photon threshold that is dependent on some given energy and adds a boolean to
evt["analysis"][outkey + "isHit"]
, the hitscore toevt["analysis"][outkey + "hitscore"]
and the limit to ``evt[“analysis”][outkey + “photonLimit”]- Args:
- evt:
The event variable
- photonscore_record:
A
Record()
object containting a photon score, e.g. total photon count
:energy_record:” A
Record()
object containting an energy value, e.g. from gas monitor detector- Kwargs:
- energyPolynomial:
array_like with polynomial coefficients fed to polyval (polynomial order one less than list length)
- outkey(str):
Prefix of data key of resulting
Record()
object, default is “photons: “
- Authors:
Carl Nettelblad (carl.nettelblad@it.uu.se)
- hummingbird.analysis.hitfinding.countTof(evt, record, signalThreshold=1, minWindow=0, maxWindow=-1, hitscoreThreshold=2, outkey='tof: ')[source]¶
A simple hitfinder that performs a peak counting test on a time-of-flight detector signal, in a specific subwindow and adds the result to
evt["analysis"][outkey + "isHit"]
, and the hitscore toevt["analysis"][outkey + "hitscore"]
.- Args:
- evt:
The event variable
- record:
A ToF detector
Record()
object
- Kwargs:
- signalThreshold(str):
The threshold of the signal, anything above this contributes to the score, default=1
- minWindow(int):
Lower limit of window, default=0
- maxWindow(int):
Upper limit of window, default=1
- hitscoreThreshold(int):
events with hitscore (Nr. of photons) above this threshold are hits, default=2
- outkey(str):
Prefix of data key of resulting
Record()
object, default is “tof: “
- Authors:
Carl Nettelblad (carl.nettelblad@it.uu.se)
- hummingbird.analysis.hitfinding.hitrate(evt, hit, history=100, unit='percent', outkey='hitrate')[source]¶
Counts hits and adds current hit rate to
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- hit:
A boolean (True for hit, False for miss)
- Kwargs:
- history(int):
Buffer length, default = 100
- outkey(str):
Data key of resulting
Record()
object, default is “hitrate”- unit(str):
Unit of hitrate, ‘fraction’ or ‘percent’, default is ‘percent’
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se), Tomas Ekeberg
- hummingbird.analysis.hitfinding.lambda_values(evt, pulse_energy, sum_over_bkg_frames, fit_bkg, sample_params, outkey='')[source]¶
analysis.sizing¶
- hummingbird.analysis.sizing.absolute_error(evt, type_a, key_a, type_b, key_b, out_key=None)[source]¶
Returning the absolute error between two records as a new record.
- hummingbird.analysis.sizing.findCenter(evt, type, key, mask=None, x0=0, y0=0, maxshift=10, threshold=0.5, blur=4)[source]¶
Estimating the center of diffraction based on pair-wise correlation enforcing friedel-symmetry and adding the estimated off center shifts cx and cy to
evt['analysis']['offCenterX']
andevt['analysis']['offCenterX']
.Note
For this function, libspimage needs to be installed.
- Args:
- evt:
The event variable
- type(str):
The event type of detectors, e.g. photonPixelDetectors
- key(str):
The event key of a detector, e.g. CCD
- Kwargs:
- mask(bool or int):
Only valid pixels (mask == True or 1) are used (default: all pixels are valid)
- x0(int):
Initial guess for off center shift in x given in pixels (default = 0)
- y0(int):
Initial guess for off center shift in y given in pixels (default = 0)
- maxshift(int):
Maximum shift (in both directions) in pixels that is used for searching optimum (default = 10)
- threshold(float):
Intensities below this threshold are set to zero (default = 0.5)
- blur(int):
Radius of the blurring kernel used to find the solution quickly (default = 4)
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se), Filipe Maia, Tomas Ekeberg
- hummingbird.analysis.sizing.fitSphere(evt, type, key, mask=None, x0=0, y0=0, d0=100, i0=1.0, wavelength=1.0, pixelsize=110, distance=1000, adu_per_photon=1, quantum_efficiency=1, material='virus', mask_radius=100, downsampling=1, brute_evals=10, photon_counting=True)[source]¶
Estimating the size of particles based on diffraction data using sphere model fitting. Adds results to
evt['analysis'][RESULT]
where RESULT is ‘size’, ‘intensity’, ‘centerx’, ‘centery’, ‘goodness’.Note
For this function, libspimage needs to be installed.
- Args:
- evt:
The event variable
- type(str):
The event type of detectors, e.g. photonPixelDetectors
- key(str):
The event key of a detector, e.g. CCD
- Kwargs:
- x0(int):
Initial guess for off center shift in x (default = 0)
- y0(int):
Initial guess for off center shift in y (default = 0)
- d0(int):
Initial guess for diameter [nm] (default = 100)
- i0(int):
Initial guess for intensity [mJ/um2] (default = 1)
- wavelength(float):
Photon wavelength [nm] (default = 1)
- pixelsize(int):
Side length of a pixel [um] (default=110)
- distance(int):
Distance from interaction to detector [mm] (default = 1000)
- adu_per_photon(int):
ADUs per photon (default = 1)
- quantum_efficiency(float):
Quantum efficiency of the detector (default = 1)
- material(str):
Material of particle, e.g. virus, protein, water, … (default = virus)
- mask_radius(int):
Radius in pixels used for circular mask defining valid pixels for fitting (default=100)
- downsampling(int):
Factor of downsampling, 1 means no downsampling (default = 1)
- brute_evals(int):
Nr. of brute force evaluations for estimating the size (default = 10)
- photon_counting(bool):
If True, Do photon conversion (discretization) before fitting (default = True)
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se), Max Hantke, Filipe Maia
- hummingbird.analysis.sizing.fitSphereRadial(evt, type, radial_distance_key, radial_average_key, mask_r=None, d0=100, i0=1.0, wavelength=1.0, pixelsize=110, distance=1000, adu_per_photon=1, quantum_efficiency=1, material='virus', mask_radius=100, brute_evals=10, photon_counting=True)[source]¶
Estimating the size of particles based on diffraction data using radial sphere model fitting. Adds results to
evt['analysis'][RESULT]
where RESULT is ‘diameter’, ‘intensity’, ‘error’.Note
For this function, libspimage needs to be installed.
- Args:
- evt:
The event variable
- type(str):
The event type of detectors, e.g. photonPixelDetectors
- key(str):
The event key of a detector radial average, e.g. radial average - CCD
- Kwargs:
- d0(int):
Initial guess for diameter [nm] (default = 100)
- i0(int):
Initial guess for intensity [mJ/um2] (default = 1)
- wavelength(float):
Photon wavelength [nm] (default = 1)
- pixelsize(int):
Side length of a pixel [um] (default=110)
- distance(int):
Distance from interaction to detector [mm] (default = 1000)
- adu_per_photon(int):
ADUs per photon (default = 1)
- quantum_efficiency(float):
Quantum efficiency of the detector (default = 1)
- material(str):
Material of particle, e.g. virus, protein, water, … (default = virus)
- mask_radius(int):
Radius in pixels used for circular mask defining valid pixels for fitting (default=100)
- brute_evals(int):
Nr. of brute force evaluations for estimating the size (default = 10)
- photon_counting(bool):
If True, Do photon conversion (discretization) before fitting (default = True)
- Authors:
Max Hantke (hantke@xray.bmc.uu.se), Benedikt J. Daurer (benedikt@xray.bmc.uu.se), Filipe Maia
- hummingbird.analysis.sizing.photon_error(evt, type_data, key_data, type_fit, key_fit, adu_per_photon)[source]¶
- hummingbird.analysis.sizing.sphereModel(evt, type, key_centerx, key_centery, key_diameter, key_intensity, shape, wavelength=1.0, pixelsize=110, distance=1000, adu_per_photon=1, quantum_efficiency=1, material='virus', poisson=False)[source]¶
Return sphere model.
Note
For this function, libspimage needs to be installed.
- Args:
- evt:
The event variable
- type(str):
The event type, e.g. analysis
- key_centerx(str):
The event key of the estimated off center shift in x
- key_centery(str):
The event key of the estimated off center shift in y
- key_diameter(str):
The event key of the estimated diameter
- key_intensity(str):
The event key of the estimated intensity
- shape(tuple):
The shape of the fit
- Kwargs:
- wavelength(float):
Photon wavelength [nm] (default = 1)
- pixelsize(int):
Side length of a pixel [um] (default=110)
- distance(int):
Distance from interaction to detector [mm] (default = 1000)
- adu_per_photon(int):
ADUs per photon (default = 1)
- quantum_efficiency(float):
Quantum efficiency of the detector (default = 1)
- material(str):
Material of particle, e.g. virus, protein, water, … (default = virus)
- poisson(bool):
If True, apply poisson sampling (default = False)
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se), Max Hantke, Filipe Maia
analysis.pixel_detector¶
- hummingbird.analysis.pixel_detector.assemble(evt, type, key, x, y, nx=None, ny=None, subset=None, outkey=None, initkey=None)[source]¶
Asesembles a detector image given some geometry and adds assembled image to
evt["analysis"]["assembled - " + key]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- x(int ndarray):
X coordinates
- y(int ndarray):
Y coordinates
- Kwargs:
- nx(int):
Total width of assembled image (zero padding)
- ny(int):
Total height of assembled image (zero padding)
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.bin(evt, type, key, binning, mask=None)[source]¶
Bin a detector image given a binning factor (and mask). Adds the records
evt["analysis"]["binned image - " + key]
andevt["analysis"]["binned mask - " + key]
.Note
This feature depends on the python package libspimage.
- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- binning(int):
The linear binning factor
- Kwargs:
- mask:
Binary mask, pixels that are masked out are not counted into the binned value.
- Authors:
Max F. Hantke (hantke@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.commonModeCSPAD2x2(evt, type, key, mask=None)[source]¶
Subtraction of common mode using median value of masked pixels (left and right half of detector are treated separately). Adds a record
evt["analysis"]["cm_corrected - " + key]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- Kwargs:
- mask:
Binary mask
- Authors:
Max F. Hantke (hantke@xray.bmc.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.commonModeLines(evt, record, outkey=None, direction='vertical')[source]¶
Common mode correction subtracting the median along lines.
- Args:
- evt:
The event variable
- record:
A pixel detector
Record`
- Kwargs:
- outkey:
The event key for the corrected detecor image, default is “corrected”
- direction:
The direction of the lines across which median is taken, default is vertical
- hummingbird.analysis.pixel_detector.commonModePNCCD(evt, type, key, outkey=None, transpose=False, signal_threshold=None, mask=None, min_nr_pixels_per_median=1)[source]¶
Common mode correction for PNCCDs.
For each row its median value is subtracted (left and right half of detector are treated separately). Adds a record
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- Kwargs:
- outkey(str):
The event key for the corrected image, default is “corrected - “ + key
- transpose(bool):
Apply procedure on transposed image
- signal_threshold(float):
Apply procedure by using only pixels below given value
- mask:
You may provide a boolean mask with values that shall be excluded from common mode correction
- min_nr_pixels_per_median(int):
Common mode correction will be skipped for pixel lines that do not contain as much as this number of values
- Authors:
Max F. Hantke (hantke@xray.bmc.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.commonModePNCCD2(evt, type, key, outkey=None, signal_threshold=None, mask=None, min_nr_pixels_per_median=1)[source]¶
Common mode correction for PNCCDs.
For each row its median value is subtracted (left and right half of detector are treated separately). Adds a record
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- Kwargs:
- outkey(str):
The event key for the corrected image, default is “corrected - “ + key
- transpose(bool):
Apply procedure on transposed image
- signal_threshold(float):
Apply procedure by using only pixels below given value
- mask:
You may provide a boolean mask with values that shall be excluded from common mode correction
- min_nr_pixels_per_median(int):
Common mode correction will be skipped for pixel lines that do not contain as much as this number of values
- Authors:
Max F. Hantke (hantke@xray.bmc.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.cropAndCenter(evt, data_rec, cx=None, cy=None, w=None, h=None, outkey='cropped')[source]¶
- hummingbird.analysis.pixel_detector.getCentral4Asics(evt, type, key)[source]¶
Adds a one-dimensional stack of its 4 centermost asics to
evt["analysis"]["central4Asics"]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- Authors:
Carl Nettelblad (carl.nettelblad@it.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.getSubsetAsics(evt, type, key, subset, output)[source]¶
Adds a one-dimensional stack of an arbitrary subset of asics to
evt["analysis"][output]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- subset(list):
Asic indices
- output:
The output key in analysis
- Authors:
Carl Nettelblad (carl.nettelblad@it.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.maxPhotonValue(evt, record, aduPhoton=1, outkey=None)[source]¶
Estimates the maximum number of photons on one pixel on the detector and adds it to
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- record:
The data record (e.g. evt[‘photonPixelDetectors’][‘CCD’])
- Kwargs:
- aduPhoton(int):
ADU count per photon, default = 1
- outkey(str):
Data key of resulting data record, default is ‘maxPhotons’
- Authors:
Tomas Ekeberg (ekeberg@xray.bmc.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.moveHalf(evt, record, vertical=0, horizontal=0, outkey='half-moved', transpose=False)[source]¶
- hummingbird.analysis.pixel_detector.pnccdGain(evt, record, gainmode)[source]¶
Returns gain (Number of ADUs per photon) based on photon energy record and gain mode.
- Args:
- evt:
The event variable
- record:
A photon energy
Record
given in eV- gainmode:
The gain mode of PNCCD (3,4,5,6) or 0 for no gain
- hummingbird.analysis.pixel_detector.radial(evt, record, mask=None, cx=None, cy=None, key='')[source]¶
Compute the radial average of a detector image given the center position (and a mask). Adds the records
evt["analysis"]["radial average - " + key]
andevt["analysis"]["radial distance - " + key]
.Note
This feature depends on the python package libspimage.
TODO: Updated docstring
- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- Kwargs:
- mask:
Binary mask, pixels that are masked out are not counted into the radial average.
- cx(float):
X-coordinate of the center position. If None the center will be in the middle.
- cy(float):
Y-coordinate of the center position. If None the center will be in the middle.
- Authors:
Max F. Hantke (hantke@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.subtractImage(evt, type, key, image, outkey=None)[source]¶
Subtract an image.
Adds a record
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- type(str):
The event type (e.g. photonPixelDetectors)
- key(str):
The event key (e.g. CCD)
- Kwargs:
- outkey(str):
The event key for the subtracted image, default is “subtraced - “ + key
- Authors:
Max F. Hantke (hantke@xray.bmc.uu.se) Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.threshold(evt, record, threshold, outkey=None)[source]¶
Set all values in an array that are lower than the threshold to zero.
- Args:
- evt:
The event variable
- record:
The data record
- threshold:
Set everything lower than this number to zero
- Kwargs:
- aduPhoton(int):
ADU count per photon, default = 1
- outkey(str):
Data key of resulting data record, default is ‘maxPhotons’
- Authors:
Tomas Ekeberg (ekeberg@xray.bmc.uu.se)
- hummingbird.analysis.pixel_detector.totalNrPhotons(evt, record, aduPhoton=1, aduThreshold=0.5, outkey=None)[source]¶
Estimates the total nr. of photons on the detector and adds it to
evt["analysis"][outkey]
.- Args:
- evt:
The event variable
- record:
The data record (e.g. evt[‘photonPixelDetectors’][‘CCD’])
- Kwargs:
- aduPhoton(int):
ADU count per photon, default = 1
- aduThreshold(int):
only pixels above this threshold given in units of ADUs are valid, default = 0.5
- outkey(str):
Data key of resulting data record, default is ‘nrPhotons’
- Authors:
Benedikt J. Daurer (benedikt@xray.bmc.uu.se)
analysis.template¶
- hummingbird.analysis.template.someAnalysis(evt, type, key, keyword=None)[source]¶
An example for an analysis module. Please document here in the docstring:
what the module is doing
what arguments need to be passed
what the module returns (adds to the event variable)
who the authors are
- Args:
- evt:
The event variable
- type(str):
The event type
- key(str):
The event key
- Kwargs:
- keyword(type):
Kewyword description (default = None)
- Authors:
Name (email), Name (email)