Peak Search and Fit (readmccd.py)

Module functions

The next documentation comes from the docstring in the header of function or class definition.

readmccd.py

IOimagefile.py

IOimagefile module is made for reading data contained in binary image file fully or partially.

More tools can be found in LaueTools package at sourceforge.net and gitlab.esrf.fr March 2020

LaueTools.IOimagefile.stringint(k, n)[source]

returns string of integer k with n zeros padding (by placing zeros before to have n characters)

Parameters:
  • k – integer to convert
  • n – nb of digits for zero padding
Returns:

string of length n containing integer k

Example: 1 -> ‘0001’ 15 -> ‘0015’

LaueTools.IOimagefile.setfilename(imagefilename, imageindex, nbdigits=4, CCDLabel=None)[source]

reconstruct filename string from imagefilename and update filename index with imageindex

Parameters:
  • imagefilename – filename string (full path or not)
  • imageindex (string) – index in filename
Return filename:
 

input filename with index replaced by input imageindex

Return type:

string

LaueTools.IOimagefile.getIndex_fromfilename(imagefilename, nbdigits=4, CCDLabel=None, stackimageindex=-1)[source]

get integer index from imagefilename string

Parameters:imagefilename – filename string (full path or not)
Returns:file index
LaueTools.IOimagefile.getfilename(dirname, imfileprefix, imfilesuffix=None, numim=None, nbdigits_filename=4)[source]

to get the global file name (name+path) for given components of the name put %4d instead of stringint

LaueTools.IOimagefile.getwildcardstring(CCDlabel)[source]

return smart wildcard to open binary CCD image file with priority of CCD type of CCDlabel

Parameters:CCDlabel (str) – string label defining the CCD type
Returns:string from concatenated strings to be used in wxpython open file dialog box
Return type:str
LaueTools.IOimagefile.getpixelValue(filename, x, y, ccdtypegeometry='edf')[source]

return pixel value at x,y

Warning

Very old function. To be checked. Use better readpixelvalue in plotdip.py

Parameters:
  • filename (str) – path to image file
  • x (int) – x pixel value
  • y (int) – y pixel value
  • ccdtypegeometry (str, optional) – CCD label, defaults to “edf”
Returns:

pixel intensity

Return type:

int

LaueTools.IOimagefile.readheader(filename, offset=4096, CCDLabel='MARCCD165')[source]

return header in a raw format

default offset for marccd image

LaueTools.IOimagefile.read_header_marccd(filename)[source]

return string of parameters found in header in marccd image file .mccd

  • print allsentences displays the header
  • use allsentences.split(‘n’) to get a list
LaueTools.IOimagefile.read_header_marccd2(filename)[source]

return string of parameters comments and exposure time found in header in marccd image file .mccd

  • print allsentences displays the header
  • use allsentences.split(‘n’) to get a list
LaueTools.IOimagefile.read_header_scmos(filename)[source]

return string of parameters comments and exposure time found in header in scmis image file .tif

  • print allsentences displays the header
  • use allsentences.split(‘n’) to get a list
LaueTools.IOimagefile.read_motorsposition_fromheader(filename, CCDLabel='MARCCD165')[source]

return xyzpositions, expo_time from image file header available for “MARCCD165”, “sCMOS”, “sCMOS_fliplr”

LaueTools.IOimagefile.readoneimage_full(filename, frametype='mccd', dirname=None)[source]

too SLOW! reads 1 entire image (marCCD format) :return: PILimage, image object of PIL module (16 bits integer) and arrayofdata: 2D array of intensity #TODO: manage framedim like readoneimage() just below

LaueTools.IOimagefile.readCCDimage(filename, CCDLabel='MARCCD165', dirname=None, stackimageindex=-1, verbose=0)[source]

Read raw data binary image file.

Read raw data binary image file and return pixel intensity 2D array such as to fit the data (2theta, chi) scattering angles representation convention.

Parameters:
  • filename (str) – path to image file (fullpath if ` dirname` =None)
  • CCDLabel (str, optional) – label, defaults to “MARCCD165”
  • dirname (str, optional) – folder path, defaults to None
  • stackimageindex (int, optional) – index of images bunch, defaults to -1
  • verbose (int, optional) – 0 or 1, defaults to 0
Raises:

ValueError – if data format and CCD parameters from label are not compatible

Returns:

  • dataimage, 2D array image data pixel intensity properly oriented
  • framedim, iterable of 2 integers shape of dataimage
  • fliprot : string, key for CCD frame transform to orient image

Return type:

tuple of 3 elements

LaueTools.IOimagefile.readoneimage(filename, framedim=(2048, 2048), dirname=None, offset=4096, formatdata='uint16')[source]

returns a 1d array of integers from a binary image file (full data)

Parameters:
  • filename (str) – image file name (full path if dirname=0)
  • framedim (tuple of 2 integers, optional) – detector dimensions, defaults to (2048, 2048)
  • dirname (str, optional) – folder path, defaults to None
  • offset (int, optional) – file header in byte (octet), defaults to 4096
  • formatdata (str, optional) – numpy format of raw binary image pixel value, defaults to “uint16”
Returns:

dataimage : image data pixel intensity

Return type:

1D array

LaueTools.IOimagefile.readoneimage_band(filename, framedim=(2048, 2048), dirname=None, offset=4096, line_startindex=0, line_finalindex=2047, formatdata='uint16')[source]

returns a 1d array of integers from a binary image file. Data located in band according shape of data (framedim)

Parameters:
  • filename – string path to image file (fullpath if `dirname`=None)
  • offset – integer nb of file header bytes
  • framedim – iterable of 2 integers shape of expected 2D data
  • formatdata – string key for numpy dtype to decode binary file
Returns:

dataimage, 1D array, image data pixel intensity

LaueTools.IOimagefile.readoneimage_crop_fast(filename, dirname=None, CCDLabel='MARCCD165', firstElemIndex=0, lastElemIndex=2047)[source]

Returns a 2d array of integers from a binary image file. Data are taken only from a rectangle

with respect to firstElemIndex and lastElemIndex.

Parameters:
  • filename – string, path to image file (fullpath if ` dirname`=None)
  • offset – integer, nb of file header bytes
  • framedim – iterable of 2 integers, shape of expected 2D data
  • formatdata – string, key for numpy dtype to decode binary file
Returns:

dataimage : 1D array image data pixel intensity

LaueTools.IOimagefile.readrectangle_in_image(filename, pixx, pixy, halfboxx, halfboxy, dirname=None, CCDLabel='MARCCD165', verbose=True)[source]

returns a 2d array of integers from a binary image file. Data are taken only from a rectangle centered on pixx, pixy

Returns:dataimage : 2D array, image data pixel intensity
LaueTools.IOimagefile.readoneimage_crop(filename, center, halfboxsize, CCDLabel='PRINCETON', dirname=None)[source]

return a cropped array of data read in an image file

Parameters:
  • filename – string, path to image file (fullpath if ` dirname`=None)
  • center – iterable of 2 integers, (x,y) pixel coordinates
  • halfboxsize – integer or iterable of 2 integers, ROI half size in both directions
Returns:

dataimage : 1D array, image data pixel intensity

Todo

useless?

LaueTools.IOimagefile.readoneimage_manycrops(filename, centers, boxsize, stackimageindex=-1, CCDLabel='MARCCD165', addImax=False, use_data_corrected=None)[source]

reads 1 image and extract many regions centered on center_pixel with xyboxsize dimensions in pixel unit

Parameters:
  • filename – string,fullpath to image file
  • centers – list or array of [int,int] centers (x,y) pixel coordinates
  • use_data_corrected – enter data instead of reading data from file must be a tuple of 3 elements: fulldata, framedim, fliprot where fulldata is a numpy.ndarray as output by readCCDimage()
  • boxsize – iterable 2 elements or integer boxsizes [in x, in y] direction or integer to set a square ROI
Returns:

Data, list of 2D array pixel intensity or Data and Imax

LaueTools.IOimagefile.writeimage(outputname, _header, data, dataformat=<class 'numpy.uint16'>)[source]

from data 1d array of integers with header coming from a f.open(‘imagefile’); f.read(headersize);f.close() .. warning:: header contain dimensions for subsequent data. Check before the compatibility of data with header infos(nb of byte per pixel and array dimensions

LaueTools.IOimagefile.write_rawbinary(outputname, data, dataformat=<class 'numpy.uint16'>)[source]

write a binary file without header of a 2D array

used ?

LaueTools.IOimagefile.SumImages(prefixname, suffixname, ind_start, ind_end, dirname=None, plot=0, output_filename=None, CCDLabel=None, nbdigits=0)[source]

sum images and write image with 32 bits per pixel format (4 bytes)

used?

LaueTools.IOimagefile.Add_Images2(prefixname, ind_start, ind_end, plot=0, writefilename=None, CCDLabel='MARCCD165', average=True)[source]

in dev

LaueTools.IOimagefile.Add_Images(prefixname, ind_start, ind_end, plot=0, writefilename=None)[source]

Add continuous sequence of images

Note

Add_Images2 exists

Parameters:
  • prefixname – string, prefix common part of name of files
  • ind_start – int, starting image index
  • ind_end – int, final image index
  • writefilename – string, new image filename where to write datastart (with last image file header read)
Returns:

datastart, array accumulation of 2D data from each image

LaueTools.IOimagefile.get_imagesize(framedim, nbbits_per_pixel, headersize_bytes)[source]

return size of image in byte (= 1 octet = 8 bits)

imageprocessing.py

imageprocessing module is made to modify filter data array

More tools can be found in LaueTools package at sourceforge.net and gitlab.esrf.fr March 2020

LaueTools.imageprocessing.getindices2cropArray(center, halfboxsizeROI, arrayshape, flipxycenter=False)[source]

return array indices limits to crop array data

Parameters:
  • center – iterable of 2 elements (x,y) pixel center of the ROI
  • halfboxsizeROI – integer or iterable of 2 elements half boxsize ROI in two dimensions
  • arrayshape – iterable of 2 integers maximal number of pixels in both directions
  • flipxycenter – boolean True: swap x and y of center with respect to others parameters that remain fixed
Returns:

imin, imax, jmin, jmax : 4 integers 4 indices allowing to slice a 2D np.ndarray

Todo

merge with check_array_indices()

LaueTools.imageprocessing.check_array_indices(imin, imax, jmin, jmax, framedim=None)[source]

Return 4 indices for array slice compatible with framedim

Parameters:
  • imax, jmin, jmax (imin,) – 4 integers mini. and maxi. indices in both directions
  • framedim – iterable of 2 integers shape of the array to be sliced by means of the 4 indices
Returns:

imin, imax, jmin, jmax: 4 integers mini. and maxi. indices in both directions

Todo

merge with getindices2cropArray()

LaueTools.imageprocessing.to8bits(PILimage, normalization_value=None)[source]

convert PIL image (16 bits) in 8 bits PIL image

Returns:
  • [0] 8 bits image
  • [1] corresponding pixels value array

Todo

since not used, may be deleted

LaueTools.imageprocessing.diff_pix(pix, array_pix, radius=1)[source]

returns index in array_pix which is the closest to pix if below the tolerance radius

array_pix: array of 2d pixel points pix: one 2elements pixel point

LaueTools.imageprocessing.minmax(D_array, center, boxsize, framedim=(2048, 2048), withmaxpos=False)[source]

extract min and max from a 2d array in a ROI

Obsolete? Still used in LocalMaxima_ShiftArrays()

Parameters D_array : 2D array

data array
center : iterable of 2 integers
(x,y) pixel center
boxsize : integer or iterable of 2 integers
full boxsize defined in both directions
framedim : iterable of 2 integers
shape of D_array

Return [min, max]: minimium and maximum pixel internsity in ROI [min, max],absolute_max_pos : if withmaxpos is True add in output

the absolute position of the largest pixel

#TODO: replace by scipy.ndimage.extrema # see next functions below # framedim = from dictionary of CCDs D_array shape is flip(framedim)

LaueTools.imageprocessing.getExtrema(data2d, center, boxsize, framedim, ROIcoords=0, flipxycenter=True)[source]

return min max XYposmin, XYposmax values in ROI

Parameters:
  • ROIcoords – 1 in local array indices coordinates 0 in X,Y pixel CCD coordinates
  • flipxycenter – boolean like swap input center coordinates
  • data2d – 2D array data array as read by readCCDimage()
Returns:

min, max, XYposmin, XYposmax: - min : minimum pixel intensity - max : maximum pixel intensity - XYposmin : list of absolute pixel coordinates of lowest pixel - XYposmax : list of absolute pixel coordinates of largest pixel

LaueTools.imageprocessing.getIntegratedIntensity(data2d, center, boxsize, framedim, thresholdlevel=0.2, flipxycenter=True)[source]

return crude estimate of integrated intensity of peak above a given relative threshold

Parameters:
  • ROIcoords – 1 in local array indices coordinates 0 in X,Y pixel CCD coordinates
  • flipxycenter – boolean like swap input center coordinates
  • data2d – 2D array data array as read by readCCDimage()
  • Thresholdlevel – relative level above which pixel intensity must be taken into account I(p)- minimum> Thresholdlevel* (maximum-minimum)
Returns:

integrated intensity, minimum absolute intensity, nbpixels used for the summation

LaueTools.imageprocessing.getMinMax(data2d, center, boxsize, framedim)[source]

return min and max values in ROI

Parameters:

data2d : 2D array
array as read by readCCDimage
LaueTools.imageprocessing.minmax_fast(D_array, centers, boxsize=(25, 25))[source]

extract min (considered as background in boxsize) and intensity at center from a 2d array at different places (centers)

centers is tuple a two array ( array([slow indices]), array([fast indices]))

return:

[0] background values [1] intensity value

used?

LaueTools.imageprocessing.normalize_shape(shape)[source]

return shape in case a scalar was given: return (shape,)

LaueTools.imageprocessing.LoG(r, sigma=None, dim=1, r0=None, peakVal=None)[source]

note: returns negative Laplacian-of-Gaussian (aka. mexican hat) zero-point will be at sqrt(dim)*sigma integral is _always_ 0 if peakVal is None: uses “mathematical” “gaussian derived” norm if r0 is not None: specify radius of zero-point (IGNORE sigma !!)

LaueTools.imageprocessing.LoGArr(shape=(256, 256), r0=None, sigma=None, peakVal=None, orig=None, wrap=0, dtype=<class 'numpy.float32'>)[source]

returns n-dim Laplacian-of-Gaussian (aka. mexican hat) if peakVal is not None

result max is peakVal

if r0 is not None: specify radius of zero-point (IGNORE sigma !!)

credits: “Sebastian Haase <haase@msg.ucsf.edu>”

LaueTools.imageprocessing.radialArr(shape, func, orig=None, wrap=False, dtype=<class 'numpy.float32'>)[source]

generates and returns radially symmetric function sampled in volume(image) of shape shape if orig is None the origin defaults to the center func is a 1D function with 1 paramater: r

if shape is a scalar uses implicitely (shape,) wrap tells if functions is continued wrapping around image boundaries wrap can be True or False or a tuple same length as shape:

then wrap is given for each axis sperately
LaueTools.imageprocessing.LocalMaxima_ndimage(Data, peakVal=4, boxsize=5, central_radius=2, threshold=1000, connectivity=1, returnfloatmeanpos=0, autothresholdpercentage=None)[source]

returns (float) i,j positions in array of each blob (peak, spot, assembly of hot pixels or whatever)

Note

used only in LocalMaxima_KernelConvolution

inputs

peakVal, boxsize, central_radius:
parameters for numerical convolution with a mexican-hat-like kernel
threshold:
intensity threshold of filtered Data (by convolution with the kernel) above which blob signal will be considered if = 0 : take all blobs at the expense of processing time
connectivity :
1 for filled square 3*3 connectivity 0 for 3*3 star like connectivity
autothresholdpercentage :
threshold in filtered image with respect to the maximum intensity in filtered image

output: array (n,2): array of 2 indices

LaueTools.imageprocessing.ConvolvebyKernel(Data, peakVal=4, boxsize=5, central_radius=2)[source]

Convolve Data array witn mexican-hat kernel

inputs: Data : 2D array containing pixel intensities peakVal > central_radius : defines pixel distance from box center where weights are positive

(in the middle) and negative farther to converge back to zero

boxsize : size of the box

ouput: array (same shape as Data)

LaueTools.imageprocessing.LocalMaxima_KernelConvolution(Data, framedim=(2048, 2048), peakValConvolve=4, boxsizeConvolve=5, central_radiusConvolve=2, thresholdConvolve=1000, connectivity=1, IntensityThreshold=500, boxsize_for_probing_minimal_value_background=30, return_nb_raw_blobs=0, peakposition_definition='max')[source]

return local maxima (blobs) position and amplitude in Data by using convolution with a mexican hat like kernel.

Two Thresholds are used sequently:
  • thresholdConvolve : level under which intensity of kernel-convolved array is discarded
  • IntensityThreshold : level under which blob whose local intensity amplitude in raw array is discarded
Parameters:
  • Data – 2D array containing pixel intensities
  • boxsizeConvolve, central_radiusConvolve (peakValConvolve,) – convolution kernel parameters
  • thresholdConvolve – minimum threshold (expressed in unit of convolved array intensity) under which convoluted blob is rejected.It can be zero (all blobs are accepted but time consuming)
  • connectivity

    shape of connectivity pattern to consider pixels belonging to the same blob.

    • 1: filled square (1 pixel connected to 8 neighbours)
    • 0: star (4 neighbours in vertical and horizontal direction)
  • IntensityThreshold – minimum local blob amplitude to accept
  • boxsize_for_probing_minimal_value_background – boxsize to evaluate the background and the blob amplitude
  • peakposition_definition – string (‘max’ or ‘center’) key to assign to the blob position its hottest pixel position or its center (no weight)
Returns:

peakslist : array like (n,2)

list of peaks position (pixel)

Ipixmax : array like (n,1) of integer

list of highest pixel intensity in the vicinity of each peak

npeaks : integer

nb of peaks (if return_nb_raw_blobs =1)

LaueTools.imageprocessing.LocalMaxima_ShiftArrays(Data, framedim=(2048, 2048), IntensityThreshold=500, Saturation_value=65535, boxsize_for_probing_minimal_value_background=30, nb_of_shift=25, pixeldistance_remove_duplicates=25, verbose=0)[source]

blob search or local maxima search by shift array method (kind of derivative)

Warning

Flat peak (= two neighbouring pixel with rigourouslty the same intensity) is not detected

LaueTools.imageprocessing.shiftarrays_accum(Data_array, n, dimensions=1, diags=0)[source]

idem than shiftarrays() but with all intermediate shifted arrays 1D returns 3 arrays corresponding to shifted arrays by n in two directions and original one 2D returns 5 arrays corresponding to shifted arrays by n in two directions and original one

these arrays are ready for comparison with eg np.greater

Data_array must have shape (slowdim,fastdim) so that slowdim-2*n>=1 and fastdim-2*n>=1 (ie central array with zero shift has some elements)

TODO: replace append by a pre allocated array

Note

readmccd.localmaxima is better

LaueTools.imageprocessing.LocalMaxima_from_thresholdarray(Data, IntensityThreshold=400, rois=None, framedim=None, verbose=False)[source]

return center of mass of each blobs composes by pixels above IntensityThreshold

if Centers = list of (x,y, halfboxsizex, halfboxsizey) perform only blob search in theses ROIs

Warning

center of mass of blob where all intensities are set to 1

LaueTools.imageprocessing.localmaxima(DataArray, n, diags=1, verbose=0)[source]

from DataArray 2D returns (array([i1,i2,…,ip]),array([j1,j2,…,jp])) of indices where pixels value is higher in two direction up to n pixels

this tuple can be easily used after in the following manner: DataArray[tupleresult] is an array of the intensity of the hottest pixels in array

in similar way with only four cardinal directions neighbouring (found in the web): import numpy as N def local_minima(array2d):

return ((array2d <= np.roll(array2d, 1, 0)) &
(array2d <= np.roll(array2d, -1, 0)) & (array2d <= np.roll(array2d, 1, 1)) & (array2d <= np.roll(array2d, -1, 1)))

WARNING: flat top peak are not detected !!

LaueTools.imageprocessing.gauss_kern(size, sizey=None)[source]

Returns a normalized 2D gauss kernel array for convolutions

LaueTools.imageprocessing.blur_image(im, n, ny=None)[source]
blurs the image by convolving with a gaussian kernel of typical
size n. The optional keyword argument ny allows for a different size in the y direction.
LaueTools.imageprocessing.blurCCD(im, n)[source]

apply a blur filter to image ndarray

LaueTools.imageprocessing.circularMask(center, radius, arrayshape)[source]

return a boolean ndarray of elem in array inside a mask

LaueTools.imageprocessing.compute_autobackground_image(dataimage, boxsizefilter=10)[source]

return 2D array of filtered data array :param dataimage: array of image data :type dataimage: 2D array

LaueTools.imageprocessing.computefilteredimage(dataimage, bkg_image, CCDlabel, kernelsize=5, formulaexpression='A-B', usemask=True)[source]

return 2D array of initial image data without background given by bkg_image data

usemask : True then substract bkg image on masked raw data
False apply formula on all pixels (no mask)
Parameters:
  • dataimage (2D array) – array of image data
  • bkg_image (2D array) – array of filtered image data (background)
  • CCDlabel (string) – key for CCD dictionary
LaueTools.imageprocessing.filterimage(image_array, framedim, blurredimage=None, kernelsize=5, mask_parameters=None, clipvalues=None, imageformat=<class 'numpy.uint16'>)[source]

compute a difference of images inside a region defined by a mask

Parameters:
  • blurredimage – ndarray image to substract to image_array
  • kernelsize – pixel size of gaussian kernel if blurredimage is None
  • mask_parameters – circular mask parameter: center=(x,y), radius, value outside mask
LaueTools.imageprocessing.rebin2Darray(inputarray, bin_dims=(2, 2), operator='mean')[source]

rebin 2D array by applying an operator to define the value of one element from the other

Parameters:
  • operator – mean, min, max, sum
  • bin_dims – side sizes of binning. (2,3) means 2X3
LaueTools.imageprocessing.blurCCD_with_binning(im, n, binsize=(2, 2))[source]

blur the array by rebinning before and after aplying the filter

LaueTools.imageprocessing.filter_minimum(im, boxsize=10)[source]

return filtered image using minimum filter

LaueTools.imageprocessing.remove_minimum_background(im, boxsize=10)[source]

remove to image array the array resulting from minimum_filter

LaueTools.imageprocessing.plot_image_markers(image, markerpos, position_definition=1)[source]

plot 2D array (image) with markers at first two columns of (markerpos)

Note

used in LaueHDF5. Could be better implementation in some notebooks

LaueTools.imageprocessing.applyformula_on_images(A, B, formulaexpression='A-B', SaturationLevel=None, clipintensities=True)[source]

calculate image data array from math expression

Parameters:
  • B (A,) – ndarray of the same shape
  • SaturationLevel – saturation level of intensity
  • clipintensities – clip resulting intensities to zero and saturation value