pyimgalgos documentation¶
Class FiberAngles
- contains methods to evaluate angles in fiber diffraction experiments¶
- Usage::
# Import from pyimgalgos.FiberAngles import fraser, calc_phi, calc_beta, funcy
# Fraser transformation: s12rot, s3rot, fraser_img = fraser(arr, beta_deg, dist_pix, center=None, oshape=(1500,1500))
# HBins objects for Fraser transformation: qh_hbins, qv_hbins = fraser_bins(fraser_img, dist_pix)
# Evaluation of fiber tilt angles beta phi (in the image plane) (in transverse to image plane). phi = calc_phi (x1, y1, x2, y2, dist) beta = calc_beta(x1, y1, phi, dist)
#Fit functions yarr = funcy(xarr, phi_deg, bet_deg)
yarr = funcy_l1_v0(xarr, phi_deg, bet_deg, DoR=0.4765, sgnrt=-1.) (depric.) yarr = funcy_l1_v1(xarr, phi_deg, bet_deg, DoR=0.4765, sgnrt=-1.)
yarr2 = funcy2(xarr, a, b, c)
# Conversion methods qx, qy, qz = recipnorm(x, y, z) # returns q/fabs(k) components for 3-d point along k^prime.
xe, ye = rqh_to_xy(re, qh) # Returns reciprocal (xe,ye) coordinates of the qh projection on Ewald sphere of radius re. xe, ye, ze = rqhqv_to_xyz(re, qh, qv) # Returns reciprocal (xe,ye,ze) coordinates of the qh,qv projection on Ewald sphere of radius re.
# Commands to test in the release directory: # python ./pyimgalgos/src/FiberAngles.py <test-id> # where # <test-id> = 1 - test of the Fraser transformation # <test-id> = 2 - test of the phi angle # <test-id> = 3 - test of the beta angle
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
- See:
Created in 2015 by Mikhail Dubrovin
-
FiberAngles.
calc_beta
(xpix, ypix, phi, dist)[source]¶ Evaluates fiber beta angle [rad] for two peaks in equatorial region - xpix - [in] x coordinate of the point - ypix - [in] y coordinate of the point - phi - [in] fiber tilt angle [rad] if the detector plane - dist - [in] distance from sample to detector
-
FiberAngles.
calc_phi
(x1pix, y1pix, x2pix, y2pix, dist)[source]¶ Evaluates fiber phi angle [rad] for two peaks in equatorial region - x1pix - [in] x coordinate of the 1st point - y1pix - [in] y coordinate of the 1st point - x1pix - [in] x coordinate of the 2nd point - y1pix - [in] y coordinate of the 2nd point - dist - [in] distance from sample to detector
-
FiberAngles.
fraser
(arr, beta_deg, z, center=None, oshape=(1500, 1500))[source]¶ Do Fraser correction for an array at angle beta and sample-to-detector distance z, given in number of pixels (109.92um), angle is in degrees. Example: fraser(array,10,909); (10 degrees at 100mm distance)
ASSUMPTION: 1. by default 2-d arr image center corresponds to (x,y) origin - arr - [in] 2-d image array - beta_deg - [in] angle beta in degrees - z - [in] distance from sample to detector given in units of pixel size (110um) - center - [in] center (row,column) location on image, which will be used as (x,y) origin - oshape - [in] ouitput image shape
-
FiberAngles.
fraser_bins
(fraser_img, dist_pix, dqv=0)[source]¶ Returns horizontal and vertical HBins objects for pixels in units of k=1 Fraser imaging array, returned by method fraser(…). Units: sample to detector distance dist_pix given in pixels, dqv - normalized offset for qv (for l=1 etc.)
-
FiberAngles.
fraser_xyz
(x, y, z, beta_deg, k=1.0)[source]¶ Do Fraser transformation for of 3-d points given by x,y,z arrays for angle beta around x and returns horizontal and vertical components of the scattering vector in units of k (1(def)-relative or eV or 1/A). x,y-arrays representing image point coordinates, z-can be scalar - distance from sample to detector. x, y, and z should be in the same units; ex.: number of pixels (109.92um) or [um], angle is in degrees. Example: fraser_xy(x,y,909,10); (10 degrees at 909 pixel (100mm) distance)
ASSUMPTION: - x,y,z - [in] point coordinate arrays with originn in IP - beta_deg - [in] angle beta in degrees - k - [in] scale factor, ex.: wave number in units of (eV or 1/A).
-
FiberAngles.
funcy
(x, phi_deg, bet_deg)¶ Function for parameterization of y(x, phi, beta) of peaks in mediane plane for fiber diffraction ATTENTION!: curve_fit assumes that x and returned y are numpy arrays.
-
FiberAngles.
funcy_l0
(x, phi_deg, bet_deg)[source]¶ Function for parameterization of y(x, phi, beta) of peaks in mediane plane for fiber diffraction ATTENTION!: curve_fit assumes that x and returned y are numpy arrays.
-
FiberAngles.
funcy_l1_v0
(x, phi_deg, bet_deg, DoR=0.474, sgnrt=-1.0)[source]¶ DEPRICATED: D/L - is not a constant as it should be Function for parameterization of y(x, phi, beta). v0: EQUATION FOR D/L=… of peaks in l=1 plane for fiber diffraction DoR = d/R ratio, where d is a distance between l=0 and l=1 on image, DoR = 433/913.27 - default constant R is sample to detector distance ATTENTION!: curve_fit assumes that x and returned y are numpy arrays.
-
FiberAngles.
funcy_l1_v1
(x, phi_deg, bet_deg, DoR=0.4292, sgnrt=1.0)[source]¶ v0: EQUATION FOR D/R. Function for parameterization of y(x, phi, beta) of peaks in l=1 plane for fiber diffraction DoR = d/R ratio, where d is a distance between l=0 and l=1 on image, DoR = 392/913.27 - default R is sample to detector distance ATTENTION!: curve_fit assumes that x and returned y are numpy arrays.
-
FiberAngles.
funcy_v0
(x, phi_deg, bet_deg)¶ Function for parameterization of y(x, phi, beta) of peaks in mediane plane for fiber diffraction ATTENTION!: curve_fit assumes that x and returned y are numpy arrays.
-
FiberAngles.
recipnorm
(x, y, z)[source]¶ Returns normalizd reciprocal space coordinates (qx,qy,qz) of the scattering vector q = (k^prime - k)/abs(k), and assuming that - scattering point is a 3-d space origin, also center of the Ewalds sphere - k points from 3-d space origin to the point with coordinates x, y, z (pixel coordinates relative to IP) - scattering is elastic, no energy loss or gained, abs(k^prime)=abs(k) - reciprocal space origin is in the intersection point of axes z and Ewald’s sphere.
-
FiberAngles.
rotation
(X, Y, angle_deg)[source]¶ For numpy arrays X and Y returns the numpy arrays of Xrot and Yrot rotated by angle_deg
-
FiberAngles.
rotation_cs
(X, Y, C, S)[source]¶ For numpy arrays X and Y returns the numpy arrays of Xrot and Yrot
-
FiberAngles.
rotation_phi_beta
(x, y, z, phi_deg, beta_deg, scale)[source]¶ Returns horizontal and vertical components of the scattering vector in units of scale (k) x, y can be arrays, z-scalar in the same units, ex. scale = k[1/A] or in number of pixels etc.
-
FiberAngles.
rqh_to_xz
(re, qh)[source]¶ Returns reciprocal (qxe,qze) coordinates of the qh projection on Ewald sphere of radius re.
Parameters
- re - (float scalar) Ewald sphere radius (1/A)
- qh - (numpy array) horizontal component of q values (1/A)
- Assumption:
- reciprocal frame origin (0,0) is on Ewald sphere, center of the Ewald sphere is in (-re,0), qh is a length of the Ewald sphere chorde from origin to the point with peak. NOTE: qh, sina, cosa, qxe, qze - can be numpy arrays shaped as qh. Returns: qxe, qze - coordinates of the point on Ewald sphere equivalent to q(re,qh); Ewald sphere frame has an experiment/detector-like definition of axes; - qze - along the beam, - qxe - transverse to the beam in l=0 horizontal plane.
-
FiberAngles.
rqhqv_to_xyz
(re, qh, qv)[source]¶ Returns reciprocal (xe,ye,ze) coordinates of the q(re,qh,qv) projections on Ewald sphere frame. Reciprocal frame has origin on Ewald sphere and experiment/detector-like definition of axes; ze - along the beam, xe, ye - transverse to the beam in horizontal and vertical plane, respectively. Need in this method for l=1 or other 3-d cases.
FiberIndexing
collection of methods for fiber indexing¶
This software was developed in co-operation with Meng for analysis of data cxif5315.
- See:
Created on Oct 7, 2015 by Mikhail Dubrovin
-
class
FiberIndexing.
BinPars
(edges, nbins, vtype=<type 'numpy.float32'>, endpoint=False)[source]¶ Bin parameters storage
-
FiberIndexing.
check_triclinic_primitive_vectors
(a, b, c, alp, bet, gam, vrb=True)[source]¶ - prints input parameters of primitive vectors,
- prints three primitive vectors and their reconstructed parameters,
- reitrive parameters of primitive vectors,
- compare with input and print results of comparison.
-
FiberIndexing.
fill_row
(dr, qv, qh, h, k, l, sigma_ql, sigma_qt, bpq)[source]¶ Returns histogram array (row) for horizontal q component filled by probability to see the peak, modulated by the Gaussian function of dr, where dr is a radial distance between the lattice node and Evald’s sphere.
-
FiberIndexing.
lattice
(b1=(1.0, 0.0, 0.0), b2=(0.0, 1.0, 0.0), b3=(0.0, 0.0, 1.0), hmax=3, kmax=2, lmax=1, cdtype=<type 'numpy.float32'>, hmin=None, kmin=None, lmin=None)[source]¶ returns n-d arrays of 3d coordinates or 2d(if lmax=0) for 3d lattice and Miller hkl indices
-
FiberIndexing.
make_lookup_table
(b1=(1.0, 0.0, 0.0), b2=(0.0, 1.0, 0.0), b3=(0.0, 0.0, 1.0), hmax=3, kmax=2, lmax=1, cdtype=<type 'numpy.float32'>, evald_rad=3, sigma_q=0.001, fout=None, bpq=None, bpomega=None, bpbeta=None, hmin=None, kmin=None, lmin=None)[source]¶ Depricated, see make_lookup_table_v2 with sigma_ql, sigma_qt in stead of sigma_q
-
FiberIndexing.
make_lookup_table_v2
(b1=(1.0, 0.0, 0.0), b2=(0.0, 1.0, 0.0), b3=(0.0, 0.0, 1.0), hmax=3, kmax=2, lmax=1, cdtype=<type 'numpy.float32'>, evald_rad=3, sigma_ql=0.001, sigma_qt=0.001, fout=None, bpq=None, bpomega=None, bpbeta=None, hmin=None, kmin=None, lmin=None)[source]¶ Makes lookup table - crystal lattice nodes information as a function of angle beta and omega, where beta [deg] - fiber axis tilt, omega [deg] - fiber rotation around axis, For each crysal orientation (beta, gamma) lookup table contains info about lattice nodes closest to the Evald’s sphere:
- # beta 0.00 omega 52.50 degree
- # index beta omega h k l dr[1/A] R(h,k,l) qv[1/A] qh[1/A] qt[1/A] ql[1/A] P(omega)
- 106 0.00 52.50 -2 -1 0 -0.002756 0.123157 0.000000 -0.123478 -0.122470 -0.015745 0.165321
- 106 0.00 52.50 1 1 0 -0.000249 0.072533 0.000000 0.072551 0.072347 -0.005436 0.985422
- 106 0.00 52.50 3 5 0 0.000687 0.273564 0.000000 0.273369 0.262250 -0.077171 0.894200
where:
- index - orientation index (just an unique integer number)
- beta, omega [deg] - crystal orientation angles,
- h, k, l - Miller indeces
- dr [1/A] - distance between lattice node and Evald’s sphere
- R(h,k,l) [1/A] - distance between nodes (h,k,l) and (0,0,0)
- qv, qh [1/A] - vertical and horizontal components of scattering vector q
- qt, ql [1/A] - transverse (in horizontal plane) and longitudinal components of vector q
- P(omega) - un-normalized probability (<1) evaluated for dr(omega) using sigma_ql.
File name is generated automatically with current time stamp like lut-cxif5315-r0169-2015-10-23T14:58:36.txt
Input parameters: b1, b2, b3 - reciprocal lattice primitive vectors, hmax, kmax, lmax - lattice node indeces cdtype - data type for lattice node coordinates, evald_rad - Evald’s sphere radius, sigma_ql - expected q resolution along k (due to crystal rotation), sigma_qt - expected qt resolution (in detector plane), fout - open output file object, bpq, bpomega, bpbeta - binning parameters for q, omega, and beta NOTE: Units of b1, b2, b3, evald_rad, and sigma_q should be the same, for example [1/A].
Returns 2-d numpy array for image; summed for all beta probobility(omega vs. q_horizontal).
-
FiberIndexing.
parameters_of_primitive_vectors
(a1, a2, a3)[source]¶ Returns a, b, c, alpha, beta, gamma for three primitive vectors a1, a2, a3.
-
FiberIndexing.
plot_lattice
(b1=(1.0, 0.0, 0.0), b2=(0.0, 1.0, 0.0), b3=(0.0, 0.0, 1.0), hmax=3, kmax=2, lmax=1, cdtype=<type 'numpy.float32'>, evald_rad=0.5, qtol=0.01, prefix='', do_movie=False, delay=400, hmin=None, kmin=None, lmin=None, title_add='')[source]¶ Plots 2-d reciprocal space lattice, evald sphere, generates series of plots for rotated lattice and movie from these plots.
- do_movie = True/False - on/off production of movie
- delay - is a time in msec between movie frames.
-
FiberIndexing.
print_nda
(nda, cmt, fmt=' %8.4f')[source]¶ Prints ndarray and its shape with preceded comment.
-
FiberIndexing.
print_primitive_vectors
(a1, a2, a3, fmt='%10.6f')[source]¶ Prints three primitive vectors and their derived parameters.
-
FiberIndexing.
q_components
(X, Y, Z, evald_rad=0.5)[source]¶ For all defined nodes of the lattice returns dr - distance from evald sphere to the reciprocal lattice node, qv, qh - vertical, horizontal components of the momentum transfer vector. NOTE: X, Y, Z, DX, L, dr, qv, qh, ql, qy, ql are the numpy arrays with shape=(2*lmax+1, 2*kmax+1, 2*hmax+1), evald_rad is a scalar
-
FiberIndexing.
reciprocal_from_bravias
(a1, a2, a3)[source]¶ Returns reciprocal primitive vectors from 3-d Bravias primitive vectors using crystallographer’s definition for conversion as 1/d (2*pi/d - comes from natural physics definition).
-
FiberIndexing.
rotation
(X, Y, angle_deg)[source]¶ For numpy arrays X and Y returns the numpy arrays of Xrot and Yrot rotated by angle_deg
-
FiberIndexing.
rotation_cs
(X, Y, c, s)[source]¶ For numpy arrays X and Y returns the numpy arrays of Xrot and Yrot for specified rotation angle cosine and sine values.
-
FiberIndexing.
round_vzeros
(v, d=10)[source]¶ Returns input vector with rounded to zero components which precision less than requested number of digits.
-
FiberIndexing.
str_omega_drhkl
(ind, beta_deg, omega_deg, dr, r, qv, qh, qt, ql, h, k, l, sigma_ql)[source]¶ Returns the record to save in look-up table or print.
-
FiberIndexing.
triclinic_primitive_vectors
(a=18.36, b=26.65, c=4.81, alpha=90, beta=90, gamma=102.83)[source]¶ Returns 3-d (Bravias) primitive vectors directed along crystal axes (edges) from lattice cell edge lengths [Angstrom or other prefered units] and angles [degree] for triclinic crystal cell parametes:
*----------* / \ / \ / \ / \ / \ gamma \ / *----------* / / / / /alpha / / / *-------/--* c \ / \ beta/ a / \ / \ / \ / *-----b----*
where a, b, c - crystal cell edge lengths, alpha, beta, gamma - interaxial angles around a, b, c edges, respectively’ By design, a1 vector for edge a is directed along x, a2 vector for edge b is in x-y plane, has (x,y,0) components only, a3 vector for edge c has (x,y,z) components.
See geometry details in my logbook p.7-8.
NDArrGenerators
wrapping methods for numpy random array generators¶
Usage:
import pyimgalgos.NDArrGenerators as ag
# Methods
v = ag.prod_of_elements(arr, dtype=np.int)
size = ag.size_from_shape() # use nda.size()
shape = ag.shape_as_2d(sh)
shape = ag.shape_as_3d(sh)
arr2d = ag.reshape_to_2d(nda)
arr3d = ag.reshape_to_3d(nda)
nda = ag.random_standard(shape=(40,60), mu=200, sigma=25, dtype=np.float)
nda = ag.random_exponential(shape=(40,60), a0=100, dtype=np.float)
nda = ag.random_one(shape=(40,60), dtype=np.float)
nda = ag.random_256(shape=(40,60), dtype=np.uint8)
nda = ag.random_xffffffff(shape=(40,60), dtype=np.uint32, add=0xff000000)
nda = ag.aranged_array(shape=(40,60), dtype=np.uint32)
ag.print_ndarr(nda, name='', first=0, last=5)
nda = ag.ring_intensity(r, r0, sigma)
ag.add_ring(arr2d, amp=100, row=4.3, col=5.8, rad=100, sigma=3)
peaks = ag.add_random_peaks(arr2d, npeaks=10, amean=100, arms=50, wmean=2, wrms=0.1)
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created on Nov 23, 2015 by Mikhail Dubrovin
-
NDArrGenerators.
add_random_peaks
(arr2d, npeaks=10, amean=100, arms=50, wmean=2, wrms=0.1)[source]¶ Returns 2-d array with peaks.
-
NDArrGenerators.
add_ring
(arr2d, amp=100, row=4.3, col=5.8, rad=100, sigma=3)[source]¶ Adds peak Gaussian-shaped peak intensity to numpy array arr2d Parameters ———- arr2d : np.array - 2-d numpy array amp : float - ring intensity row : float - ring center row col : float - ring center col rad : float - ring mean radius sigma : float - width of the peak
-
NDArrGenerators.
aranged_array
(shape=(40, 60), dtype=<type 'numpy.uint32'>)[source]¶ Returns numpy array of requested shape and type filling with ascending integer numbers.
-
NDArrGenerators.
print_ndarr
(nda, name='', first=0, last=5)[source]¶ Prints array attributes, title, and a few elements in a single line.
-
NDArrGenerators.
prod_of_elements
(arr, dtype=<type 'int'>)[source]¶ Returns product of sequence elements
-
NDArrGenerators.
random_1
(shape=(40, 60), dtype=<type 'float'>)¶ Returns numpy array of requested shape and type filled with random numbers in the range [0,255].
-
NDArrGenerators.
random_256
(shape=(40, 60), dtype=<type 'numpy.uint8'>)[source]¶ Returns numpy array of requested shape and type filled with random numbers in the range [0,255].
-
NDArrGenerators.
random_exponential
(shape=(40, 60), a0=100, dtype=<type 'float'>)[source]¶ Returns numpy array of requested shape and type filled with exponential distribution for width a0.
-
NDArrGenerators.
random_one
(shape=(40, 60), dtype=<type 'float'>)[source]¶ Returns numpy array of requested shape and type filled with random numbers in the range [0,255].
-
NDArrGenerators.
random_standard
(shape=(40, 60), mu=200, sigma=25, dtype=<type 'float'>)[source]¶ Returns numpy array of requested shape and type filled with normal distribution for mu and sigma.
-
NDArrGenerators.
random_xffffffff
(shape=(40, 60), dtype=<type 'numpy.uint32'>, add=4278190080)[source]¶ Returns numpy array of requested shape and type filled with random numbers in the range [0,0xffffff] with bits 0xff000000 for alpha mask.
-
NDArrGenerators.
ring_intensity
(r, r0, sigma)[source]¶ returns numpy array with ring intensity distribution modulated by Gaussian(r-r0,sigma). Parameters ———- r : np.array - numpy array of radius (i.e. radios for each pixel) r0 : float - radius of the ring sigma : float - width of the ring
-
NDArrGenerators.
shape_as_2d
(sh)[source]¶ Returns 2-d shape for n-d shape if n>2, otherwise returns unchanged shape.
GlobalGraphics
- collection of global graphical methods¶
Usage:
import pyimgalgos.GlobalGraphics as gg
# Methods
fig, axim, axcb = gg.fig_axes(figsize=(13,12), title='Image', dpi=80, win_axim=(0.05, 0.03, 0.87, 0.93), win_axcb=(0.923, 0.03, 0.02, 0.93))
fig, axim, axcb, imsh = gg.fig_axim_axcb_imsh(figsize=(13,12), title='Image', dpi=80, win_axim=(0.05, 0.03, 0.87, 0.93),
win_axcb=(0.923, 0.03, 0.02, 0.93), arr2d=np.zeros((10,10)), origin='upper')
gg.plot_imgcb(fig, axim, axcb, imsh, arr2d, amin=None, amax=None, origin='upper', title=None, cmap='inferno')
gg.plot_img(img, mode=None, amin=None, amax=None, cmap='inferno')
gg.plot_peaks_on_img(peaks, axim, iX, iY, color='w', pbits=0, lw=2)
size = gg.size_of_shape(shape=(2,3,8))
arr = gg.getArrangedImage(shape=(40,60))
arr = gg.getRandomImage(mu=200, sigma=25, shape=(40,60))
hi = gg.getImageAs2DHist(iX,iY,W=None)
img = gg.getImageFromIndexArrays(iX,iY,W=None)
fig, axhi, hi = gg.plotHistogram(arr, amp_range=None, figsize=(6,6), bins=None, title='', window=(0.15, 0.10, 0.78, 0.82))
fig, axhi, hi = gg.hist1d(arr, bins=None, amp_range=None, weights=None, color=None, show_stat=True,
log=False, figsize=(6,5), axwin=(0.15, 0.12, 0.78, 0.80), title=None,
xlabel=None, ylabel=None, titwin=None)
axim = gg.plotImageLarge(arr, img_range=None, amp_range=None, figsize=(12,10), title='Image', origin='upper',
window=(0.05, 0.03, 0.94, 0.94), cmap='inferno')
gg.plotImageAndSpectrum(arr, amp_range=None, cmap='inferno')
fig, ax = gg.plotGraph(x,y, figsize=(5,10), window=(0.15, 0.10, 0.78, 0.86), pfmt='b-', lw=1)
gg.drawCircle(axes, xy0, radius, linewidth=1, color='w', fill=False)
gg.drawCircle(axes, xy0, radius, linewidth=1, color='w', fill=False)
gg.drawCenter(axes, xy0, s=10, linewidth=1, color='w')
gg.drawLine(axes, xarr, yarr, s=10, linewidth=1, color='w')
gg.drawRectangle(axes, xy, width, height, linewidth=1, color='w')
gg.save(fname='img.png', do_save=True, pbits=0377)
gg.save_fig(fig, fname='img.png', do_save=True, pbits=0377)
gg.move(x0=200,y0=100)
gg.move_fig(fig, x0=200, y0=100)
gg.show(mode=None)
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created in 2015 by Mikhail Dubrovin
-
GlobalGraphics.
fig_axes
(figsize=(13, 12), title='Image', dpi=80, win_axim=(0.05, 0.03, 0.87, 0.93), win_axcb=(0.923, 0.03, 0.02, 0.93))[source]¶ Creates and returns figure, and axes for image and color bar
-
GlobalGraphics.
fig_axim_axcb_imsh
(figsize=(13, 12), title='Image', dpi=80, win_axim=(0.05, 0.03, 0.87, 0.93), win_axcb=(0.923, 0.03, 0.02, 0.93), arr2d=array([[ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.], [ 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.]]), origin='upper')[source]¶ Creates and returns figure, axes for image and color bar, imshow object
-
GlobalGraphics.
getImageAs2DHist
(iX, iY, W=None)[source]¶ Makes image from iX, iY coordinate index arrays and associated weights, using np.histogram2d(…).
-
GlobalGraphics.
getImageFromIndexArrays
(iX, iY, W=None)[source]¶ Makes image from iX, iY coordinate index arrays and associated weights, using indexed array.
-
GlobalGraphics.
hist1d
(arr, bins=None, amp_range=None, weights=None, color=None, show_stat=True, log=False, figsize=(6, 5), axwin=(0.15, 0.12, 0.78, 0.8), title=None, xlabel=None, ylabel=None, titwin=None)[source]¶ Makes historgam from input array of values (arr), which are sorted in number of bins (bins) in the range (amp_range=(amin,amax))
-
GlobalGraphics.
plotHistogram
(arr, amp_range=None, figsize=(6, 6), bins=None, title='', window=(0.15, 0.1, 0.78, 0.82))[source]¶ Makes historgam from input array of values (arr), which are sorted in number of bins (bins) in the range (amp_range=(amin,amax))
-
GlobalGraphics.
plot_peaks_on_img
(peaks, axim, iX, iY, color='w', pbits=0, lw=2)[source]¶ Draws peaks on the top of image axes (axim) Plots peaks from array as circles in coordinates of image.
- peaks - 2-d list/tuple of peaks; first 6 values in each peak record should be (s, r, c, amax, atot, npix)
- axim - image axes
- iX - array of x-coordinate indexes for all pixels addressed as [s, r, c] - segment, row, column
- iX - array of y-coordinate indexes for all pixels addressed as [s, r, c] - segment, row, column
- color - peak-ring color
- pbits - verbosity; print 0 - nothing, +1 - peak parameters, +2 - x, y peak coordinate indexes
GlobalUtils
contains collection of global utilities with a single call algorithms¶
Usage:
# Import
# ==============
from pyimgalgos.GlobalUtils import subtract_bkgd, mask_from_windows #, ...
from pyimgalgos.NDArrGenerators import random_standard_array
# Background subtraction
# ======================
# Example for cspad, assuming all nda_*.shape = (32,185,388)
winds_bkgd = [ (s, 10, 100, 270, 370) for s in (4,12,20,28)] # use part of segments 4,12,20, and 28 to subtract bkgd
nda = subtract_bkgd(nda_data, nda_bkgd, mask=nda_mask, winds=winds_bkgd, pbits=0)
# Operations with numpy array shape
# =================================
shape = (32,185,388)
size = size_from_shape(shape) # returns 32*185*388
shape_2d = shape_as_2d(shape) # returns (32*185,388)
arr_2d = reshape_to_2d(nda) # returns re-shaped ndarray
shape = (4,8,185,388)
shape_3d = shape_as_3d(shape) # returns (32,185,388)
arr_3d = reshape_to_3d(nda) # returns re-shaped ndarray
# Make mask n-d numpy array using shape and windows
# =================================================
shape = (2,185,388)
w = 1
winds = [(s, w, 185-w, w, 388-w) for s in (0,1)]
mask = mask_from_windows(shape, winds)
# Make mask as 2,3-d numpy array for a few(width) rows/cols of pixels
# ===================================================================
mask2d = mask_2darr_edges(shape=(185,388), width=2)
mask3d = mask_3darr_edges(shape=(32,185,388), width=5)
# Make mask of local maximal intensity pixels in x-y (ignoring diagonals)
# ===================================================================
# works for 2-d and 3-d arrays only - reshape if needed.
data = random_standard_array(shape=(32,185,388), mu=0, sigma=10)
mask_xy_max = locxymax(data, order=1, mode='clip')
# Get string with time stamp, ex: 2016-01-26T10:40:53
# ===================================================================
ts = str_tstamp(fmt='%Y-%m-%dT%H:%M:%S', time_sec=None)
# Converters for Cheetah
# ======================
runnum, tstamp, tsec, fid = convertCheetahEventName('LCLS_2015_Feb22_r0169_022047_197f7', fmtts='%Y-%m-%dT%H:%M:%S')
table8x8 = table_from_cspad_ndarr(nda_cspad)
nda_cspad = cspad_ndarr_from_table(table8x8)
nda_32x185x388 = cspad_psana_from_cctbx(nda_64x185x194)
nda_64x185x194 = cspad_cctbx_from_psana(nda_32x185x388)
cross_check_cspad_psana_cctbx(nda_32x185x388, nda_64x185x194)
# Safe math
# =========
res = divide_protected(num, den, vsub_zero=0)
# Single line printed for np.array
# ================================
print_ndarr(nda, name='', first=0, last=5)
# Save image in file
# ==================
image = random_standard()
save_image_tiff(image, fname='image.tiff', verb=True) # 16-bit tiff
save_image_file(image, fname='image.png', verb=True) # gif, pdf, eps, png, jpg, jpeg, tiff (8-bit only)
# Create directory
# ==================
create_directory('work-dir')
# Test
# ======================
# is implemented for test numbers from 1 to 9. Command example
# python pyimgalgos/src/GlobalUtils.py 1
See GlobalUtils
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created by Mikhail Dubrovin
-
GlobalUtils.
convertCheetahEventName
(evname, fmtts='%Y-%m-%dT%H:%M:%S')[source]¶ Converts Cheetah event name like ‘LCLS_2015_Feb22_r0169_022047_197f7’ and returns runnum, tstamp, tsec, fid = 169, ‘2015-02-22T02:20:47’, <tsec>, 197f7
-
GlobalUtils.
cross_check_cspad_psana_cctbx
(nda, arr)[source]¶ Apply two-way conversions between psana and cctbx cspad arrays and compare.
-
GlobalUtils.
cspad_cctbx_from_psana
(nda_in)[source]¶ returns cctbx array of ASICs (64, 185, 194) from cspad array (32, 185, 388)
-
GlobalUtils.
cspad_ndarr_from_table
(table8x8)[source]¶ returns cspad array with shape=(32,185,388) generated from table of 2x1s shaped as (8*185, 4*388) in style of Cheetah
-
GlobalUtils.
cspad_ndarr_from_table8x8
(table8x8)¶ returns cspad array with shape=(32,185,388) generated from table of 2x1s shaped as (8*185, 4*388) in style of Cheetah
-
GlobalUtils.
cspad_psana_from_cctbx
(nda_in)[source]¶ returns cspad array (32, 185, 388) from cctbx array of ASICs (64, 185, 194)
-
GlobalUtils.
divide_protected
(num, den, vsub_zero=0)[source]¶ Returns result of devision of numpy arrays num/den with substitution of value vsub_zero for zero den elements.
-
GlobalUtils.
list_of_windarr
(nda, winds=None)[source]¶ Converts 2-d or 3-d numpy array in the list of 2-d numpy arrays for windows @param nda - input 2-d or 3-d numpy array
-
GlobalUtils.
locxymax
(nda, order=1, mode='clip')[source]¶ For 2-d or 3-d numpy array finds mask of local maxima in x and y axes (diagonals are ignored) using scipy.signal.argrelmax and return their product.
@param nda - input ndarray @param order - range to search for local maxima along each dimension @param mode - parameter of scipy.signal.argrelmax of how to treat the boarder
-
GlobalUtils.
mask_2darr_edges
(shape=(185, 388), width=2)[source]¶ Returns mask with masked width rows/colunms of edge pixels for 2-d numpy array.
-
GlobalUtils.
mask_3darr_edges
(shape=(32, 185, 388), width=2)[source]¶ Returns mask with masked width rows/colunms of edge pixels for 3-d numpy array.
-
GlobalUtils.
mask_from_windows
(ashape=(32, 185, 388), winds=None)[source]¶ Makes mask as 2-d or 3-d numpy array defined by the shape with ones in windows. N-d shape for N>3 is converted to 3-d. @param shape - shape of the output numpy array with mask. @param winds - list of windows, each window is a sequence of 5 parameters (segment, rowmin, rowmax, colmin, colmax)
-
GlobalUtils.
mean_of_listwarr
(lst_warr)[source]¶ Evaluates the mean value of the list of 2-d arrays. @lst_warr - list of numpy arrays to evaluate per pixel mean intensity value.
-
GlobalUtils.
print_command_line_parameters
(parser)[source]¶ Prints input arguments and optional parameters
-
GlobalUtils.
save_image_file
(image, fname='image.png', verb=False)[source]¶ Saves files with type by extension gif, pdf, eps, png, jpg, jpeg, tiff (8-bit only), or txt for any other type
-
GlobalUtils.
save_image_tiff
(image, fname='image.tiff', verb=False)[source]¶ Saves image in 16-bit tiff file
-
GlobalUtils.
shape_as_2d
(sh)[source]¶ Returns 2-d shape for n-d shape if n>2, otherwise returns unchanged shape.
-
GlobalUtils.
shape_as_3d
(sh)[source]¶ Returns 3-d shape for n-d shape if n>3, otherwise returns unchanged shape.
-
GlobalUtils.
src_from_rc8x8
(row, col)[source]¶ Converts Cheetah 8x8 ASICs table row and column to seg, row, col coordinates
-
GlobalUtils.
str_tstamp
(fmt='%Y-%m-%dT%H:%M:%S', time_sec=None)[source]¶ Returns string timestamp for specified format and time in sec or current time by default
-
GlobalUtils.
subtract_bkgd
(data, bkgd, mask=None, winds=None, pbits=0)[source]¶ Subtracts numpy array of bkgd from data using normalization in windows for good pixels in mask. Shapes of data, bkgd, and mask numpy arrays should be the same. Each window is specified by 5 parameters: (segment, rowmin, rowmax, colmin, colmax) For 2-d arrays segment index is not used, but still 5 parameters needs to be specified.
Parameters
- data - numpy array for data.
- bkgd - numpy array for background.
- mask - numpy array for mask.
- winds - list of windows, each window is a sequence of 5 parameters.
- pbits - print control bits; =0 - print nothing, !=0 - normalization factor.
-
GlobalUtils.
table8x8_from_cspad_ndarr
(nda_cspad)¶ returns table of 2x1s shaped as (8*185, 4*388) in style of Cheetah generated from cspad array with size=(32*185*388) ordered as in data, shape does not matter.
Graphics
wrapping methods for matplotlib¶
Usage:
import pyimgalgos.Graphics as gr
# Methods
fig = gr.figure(figsize=(13,12), title='Image', dpi=80, facecolor='w', edgecolor='w', frameon=True, move=None)
gr.move_fig(fig, x0=200, y0=100)
gr.move(x0=200, y0=100)
gr.add_axes(fig, axwin=(0.05, 0.03, 0.87, 0.93))
gr.fig_img_cbar_axes(fig=None, win_axim=(0.05, 0.03, 0.87, 0.93), win_axcb=(0.923, 0.03, 0.02, 0.93))
gr.set_win_title(fig, titwin='Image')
gr.add_title_labels_to_axes(axes, title=None, xlabel=None, ylabel=None, fslab=14, fstit=20, color='k')
gr.show(mode=None)
gr.draw()
gr.draw_fig(fig)
gr.save_plt(fname='img.png', verb=True)
gr.save_fig(fig, fname='img.png', verb=True)
hi = gr.hist(axhi, arr, bins=None, amp_range=None, weights=None, color=None, log=False)
imsh = gr.imshow(axim, img, amp_range=None, extent=None, interpolation='nearest', aspect='auto', origin='upper', orientation='horizontal', cmap='inferno')
cbar = gr.colorbar(fig, imsh, axcb, orientation='vertical', amp_range=None)
imsh, cbar = gr.imshow_cbar(fig, axim, axcb, img, amin=None, amax=None, extent=None, interpolation='nearest', aspect='auto', origin='upper', orientation='vertical', cmap='inferno')
- See:
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created in 2015 by Mikhail Dubrovin
-
Graphics.
add_axes
(fig, axwin=(0.05, 0.03, 0.87, 0.93))[source]¶ Add axes to figure from input list of windows.
-
Graphics.
colorbar
(fig, imsh, axcb, orientation='vertical', amp_range=None)[source]¶ orientation = ‘horizontal’ amp_range = (-10,50)
-
Graphics.
fig_img_axes
(fig=None, win_axim=(0.08, 0.05, 0.89, 0.93))[source]¶ Returns figure and image axes
-
Graphics.
fig_img_cbar_axes
(fig=None, win_axim=(0.05, 0.03, 0.87, 0.93), win_axcb=(0.923, 0.03, 0.02, 0.93))[source]¶ Returns figure and axes for image and color bar
-
Graphics.
figure
(figsize=(13, 12), title='Image', dpi=80, facecolor='w', edgecolor='w', frameon=True, move=None)[source]¶ Creates and returns figure
-
Graphics.
hist
(axhi, arr, bins=None, amp_range=None, weights=None, color=None, log=False)[source]¶ Makes historgam from input array of values (arr), which are sorted in number of bins (bins) in the range (amp_range=(amin,amax))
-
Graphics.
imshow
(axim, img, amp_range=None, extent=None, interpolation='nearest', aspect='auto', origin='upper', orientation='horizontal', cmap='inferno')[source]¶ extent - list of four image physical limits for labeling, cmap: ‘jet’, ‘gray_r’, ‘inferno’ #axim.cla()
Class HBins
histogram-style bin parameters holder¶
Usage:
from pyimgalgos.HBins import HBins
# Equal bins constructor
hb = HBins((1,6), nbins=5)
# Variable bins constructor
hb = HBins((1,2,4,6,10))
# Access methods
nbins = hb.nbins() # returns int input parameter - number of bins
edges = hb.edges() # returns np.array input list of bin edges
vmin = hb.vmin() # returns vtype minimal value of bin edges
vmax = hb.vmax() # returns vtype maximal value of bin edges
vtype = hb.vtype() # returns np.dtype - type of bin edge values
equalbins = hb.equalbins() # returns bool True/False for equal/variable size bins
limits = hb.limits() # returns np.array of limits (vmin, vmax)
binedges = hb.binedges() # returns np.array with bin edges of size nbins+1
binedgesleft = hb.binedgesleft() # returns np.array with bin left edges of size nbins
binedgesright = hb.binedgesright() # returns np.array with bin rignt edges of size nbins
bincenters = hb.bincenters() # returns np.array with bin centers of size nbins
binwidth = hb.binwidth() # returns np.array with bin widths of size nbins or scalar bin width for equal bins
halfbinw = hb.halfbinw() # returns np.array with half-bin widths of size nbins or scalar bin half-width for equal bins
strrange = hb.strrange(fmt) # returns str of formatted vmin, vmax, nbins ex: 1-6-5
ind = hb.bin_index(value, edgemode=0) # returns bin index [0,nbins) for value.
indarr = hb.bin_indexes(valarr, edgemode=0) # returns array of bin index [0,nbins) for array of values
hisarr = hb.bin_count(valarr, edgemode=0) # returns array of bin counts [0,nbins) for array of values (histogram value per bin)
# edgemode - defines what to do with underflow overflow indexes;
# = 0 - use indexes 0 and nbins-1 for underflow overflow, respectively
# = 1 - use extended indexes -1 and nbins for underflow overflow, respectively
hb.set_bin_data(data, dtype=np.float) # adds bin data to the HBins object. data size should be equal to hb.nbins()
data = bin_data(dtype=np.float) # returns numpy array of data associated with HBins object.
# Print methods
hb.print_attrs_defined()
hb.print_attrs()
hb.print_attrs_and_methods()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created on 2016-01-15 by Mikhail Dubrovin
-
class
HBins.
HBins
(edges, nbins=None, vtype=<type 'numpy.float32'>)[source]¶ Hystogram-style bin parameters holder
-
__init__
(edges, nbins=None, vtype=<type 'numpy.float32'>)[source]¶ Class constructor for - equal bins, ex: hb = HBins((1,6), nbins=5) - or variable bins, ex: hb = HBins((1,2,4,6,10))
Parameters: - edges - sequence of two or more bin edges - nbins - (int) number of bins for equal size bins - vtype - numpy type of bin values (optional parameter)
-
Class HPolar
- makes 2-d histogram in polar (r-phi) coordinates for imaging detector n-d array data¶
Usage:
# Import
# ------
from pyimgalgos.HPolar import HPolar
# Initialization
# --------------
hp = HPolar(xarr, yarr, mask=None, radedges=None, nradbins=100, phiedges=(0,360), nphibins=32)
# Access methods
# --------------
orb = hp.obj_radbins() # returns HBins object for radial bins
opb = hp.obj_phibins() # returns HBins object for angular bins
rad = hp.pixel_rad()
irad = hp.pixel_irad()
phi0 = hp.pixel_phi0()
phi = hp.pixel_phi()
iphi = hp.pixel_iphi()
iseq = hp.pixel_iseq()
npix = hp.bin_number_of_pixels()
int = hp.bin_intensity(nda)
arr1d = hp.bin_avrg(nda)
arr2d = hp.bin_avrg_rad_phi(nda, do_transp=True)
pixav = hp.pixel_avrg(nda)
pixav = hp.pixel_avrg_interpol(nda, method='linear') # method='nearest' 'cubic'
# Print attributes and n-d arrays
# -------------------------------
hp.print_attrs()
hp.print_ndarrs()
# Global methods
# --------------
from pyimgalgos.HPolar import polarization_factor, divide_protected, cart2polar, polar2cart, bincount
polf = polarization_factor(rad, phi, z)
result = divide_protected(num, den, vsub_zero=0)
r, theta = cart2polar(x, y)
x, y = polar2cart(r, theta)
bin_values = bincount(map_bins, map_weights=None, length=None)
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created in 2015 by Mikhail Dubrovin
-
HPolar.
bincount
(map_bins, map_weights=None, length=None)[source]¶ Wrapper for numpy.bincount with protection of weights and alattening numpy arrays
-
HPolar.
divide_protected
(num, den, vsub_zero=0)[source]¶ Returns result of devision of numpy arrays num/den with substitution of value vsub_zero for zero den elements.
-
HPolar.
polar2cart
(r, theta)[source]¶ For numpy arryys r and theta returns the numpy arrays of x and y
-
HPolar.
polarization_factor
(rad, phi_deg, z)[source]¶ Returns per-pixel polarization factors, assuming that detector is perpendicular to Z.
Class HSpectrum
works with spectral histogram for arbitrary shaped numpy array¶
Usage:
# Import
# ==============
from pyimgalgos.HSpectrum import HSpectrum
# Initialization
# ==============
# 1) for bins of equal size:
range = (vmin, vmax)
nbins = 100
spec = HSpectrum(range, nbins)
# 2) for variable size bins:
bins = (v0, v1, v2, v4, v5, vN) # any number of bin edges
spec = HSpectrum(bins)
# Fill spectrum
# ==============
# nda = ... (get it for each event somehow)
spec.fill(nda)
# Get spectrum
# ==============
histarr, edges, nbins = spec.spectrum()
# Optional
# ==============
spec.print_attrs()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created in 2015 by Mikhail Dubrovin
NDArrSpectrum
- ALIAS FOR HSpectrum - support creation of spectral histogram for arbitrary shaped numpy array.
Usage:
# Import
# ==============
from pyimgalgos.NDArrSpectrum import NDArrSpectrum
# Initialization
# ==============
# 1) for bins of equal size:
range = (vmin, vmax)
nbins = 100
spec = NDArrSpectrum(range, nbins)
# 2) for variable size bins:
bins = (v0, v1, v2, v4, v5, vN) # any number of bin edges
spec = NDArrSpectrum(bins)
# Fill spectrum
# ==============
# nda = ... (get it for each event somehow)
spec.fill(nda)
# Get spectrum
# ==============
histarr, edges, nbins = spec.spectrum()
# Optional
# ==============
spec.print_attrs()
@see pyimgalgos.HSpectrum
,
pyimgalgos.HBins
pyimgalgos.HPolar
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Revision: $Revision$
@version $Id$
@author Mikhail S. Dubrovin
RadialBkgd
- radial background subtraction for imaging detector n-d array data,- extension of the base class
HPolar
for methods subtract_bkgd and subtract_bkgd_interpol.
Usage:
# Import
# ------
from pyimgalgos.RadialBkgd import RadialBkgd
# Initialization
# --------------
rb = RadialBkgd(xarr, yarr, mask=None, radedges=None, nradbins=100, phiedges=(0,360), nphibins=32)
# Access methods
# --------------
orb = rb.obj_radbins() # returns HBins object for radial bins
opb = rb.obj_phibins() # returns HBins object for angular bins
rad = rb.pixel_rad()
irad = rb.pixel_irad()
phi0 = rb.pixel_phi0()
phi = rb.pixel_phi()
iphi = rb.pixel_iphi()
iseq = rb.pixel_iseq()
npix = rb.bin_number_of_pixels()
int = rb.bin_intensity(nda)
arr1d = rb.bin_avrg(nda)
arr2d = rb.bin_avrg_rad_phi(nda, do_transp=True)
pixav = rb.pixel_avrg(nda)
pixav = rb.pixel_avrg_interpol(nda, method='linear') # method='nearest' 'cubic'
cdata = rb.subtract_bkgd(nda)
cdata = rb.subtract_bkgd_interpol(nda, method='linear')
# Print attributes and n-d arrays
# -------------------------------
rb.print_attrs()
rb.print_ndarrs()
# Global methods
# --------------
from pyimgalgos.RadialBkgd import polarization_factor, divide_protected, cart2polar, polar2cart, bincount
polf = polarization_factor(rad, phi, z)
result = divide_protected(num, den, vsub_zero=0)
r, theta = cart2polar(x, y)
x, y = polar2cart(r, theta)
bin_values = bincount(map_bins, map_weights=None, length=None)
@see pyimgalgos.HPolar
pyimgalgos.HBins
See Radial background.
pyimgalgos.HSpectrum
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Revision: $Revision$
@version $Id$
@author Mikhail S. Dubrovin
-
RadialBkgd.
test02
(ntest, prefix='fig-v01')[source]¶ Test for 2-d (default) binning of the rad-phi range of entire image
-
RadialBkgd.
test03
(ntest, prefix='fig-v01')[source]¶ Test for 2-d binning of the restricted rad-phi range of entire image
Class helps to save text file with information about peaks found in peak-finders.
Usage:
# Imports
from pyimgalgos.PeakStore import PeakStore
# Usage with psana types env and evt
pstore = PeakStore(env, 5, prefix='xxx', add_header='TitV1 TitV2 TitV3 ...', pbits=255)
for peak in peaks :
rec = '%s %d %f ...' % (peak[0], peak[5], peak[7],...)
pstore.save_peak(evt, rec)
pstore.close() # is done by default in destructor
# Usage without psana types
pstore = PeakStore('expNNNNN', 5, prefix='xxx', add_header='TitV1 TitV2 TitV3 ...', pbits=255)
for peak in peaks :
rec = '%s %d %f ...' % (peak[0], peak[5], peak[7],...)
pstore.save_peak(peak_rec=rec)
# Print methods
pstore.print_attrs()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
@version $Id$
@author Mikhail S. Dubrovin
Class Quaternion
works with quaternion rotations¶
This software was developed in co-operation with Meng for analysis of data cxif5315.
- See:
Created in April 2016 by Mikhail Dubrovin
-
Quaternion.
normalise_quaternion
(q)[source]¶ q - quaternion
Rescales the quaternion such that its modulus is unity.
Returns: the normalised version of q
-
Quaternion.
quat_rot
(v, q)[source]¶ v - vector (in the form of a “struct rvec”) q - quaternion
Rotates a vector according to a quaternion.
Returns: rotated vector vrot
-
Quaternion.
quaternion_for_angles
(ax, ay, az)[source]¶ Returns: quaternion for three input rotation angles [deg].
-
Quaternion.
quaternion_from_rotmatrix
(m)[source]¶ m - 3-d rotation matrix, class Matrix
Evaluates quaternion from rotation matrix. Implemented as https://en.wikipedia.org/wiki/Rotation_matrix
Returns: normalised quaternion.
-
Quaternion.
quaternion_modulus
(q)[source]¶ q - quaternion
If a quaternion represents a pure rotation, its modulus should be unity. Returns: the modulus of the given quaternion.
-
Quaternion.
quaternion_valid
(q, tol=0.001)[source]¶ q - quaternion
Checks if the given quaternion is normalised.
Returns: 1 if the quaternion is normalised, 0 if not.
-
Quaternion.
record_for_angles
(ax, ay, az)[source]¶ Prints string like: Angles around x,y,z: 72.0 -3.5 176.4 quaternion: w,x,y,z: 0.007459 0.043148 0.586445 0.808804
Entropy.py
- collection of methods to evaluate data array entropy¶
Usage:
# Import
# ==============
from pyimgalgos.Entropy import *
# import for test only
from pyimgalgos.NDArrGenerators import random_standard
from pyimgalgos.GlobalUtils import print_ndarr
arr_float = random_standard(shape=(1000), mu=200, sigma=25, dtype=np.float)
arr_int16 = arr_float.astype(np.int16)
print_ndarr(arr_int16, name='arr_int16', first=0, last=10)
print 'entropy(arr_int16) = %.6f' % entropy(arr_int16)
print 'entropy_v1(arr_int16) = %.6f' % entropy_v1(arr_int16)
print 'entropy_cpo(arr_int16) = %.6f' % entropy_cpo(arr_int16)
See pyimgalgos.GlobalUtils
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
Created by Mikhail Dubrovin
-
Entropy.
entropy
(nda)[source]¶ Evaluates n-d array entropy using formula from https://en.wikipedia.org/wiki/Entropy_%28information_theory%29
-
Entropy.
entropy_cpo
(signal)[source]¶ Entropy evaluation method found by cpo on web
Function returns entropy of a signal, which is 1-D numpy array
TDFileContainer
- text/table data file container - load/hold/provide access to data from text file¶
It is assumed that text data file contains records of the same format and occasionally record-header beginning with character # (hash in [0] position). Example of the file content:
# Exp Run Date Time time(sec) time(nsec) fiduc Evnum Reg Seg Row Col ...
cxif5315 169 2015-02-22 02:20:47 1424600447 486382070 104421 0 EQU 17 153 48 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 494719789 104424 1 EQU 1 161 32 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 494719789 104424 1 EQU 17 170 51 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 503058551 104427 2 EQU 25 170 310 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 503058551 104427 2 EQU 25 180 292 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 511393301 104430 3 EQU 1 162 27 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 536405573 104439 6 ARC 8 11 41 ...
cxif5315 169 2015-02-22 02:20:47 1424600447 536405573 104439 6 ARC 8 10 20 ...
...
Header (without #) should have the same as data number of literal fields separated by spaces. Records in the file should be grupped by unique group-id parameter, for example a group of records may have the same group number or some unique index.
Originaly it is designed to work with text file containing record data generated by peak-finder. It is adopted to work with any other object type beside peak data.
Usage:
# !!! NOTE: None is returned whenever requested information is missing.
# Import
from pyimgalgos.TDFileContainer import TDFileContainer
from pyimgalgos.TDPeakRecord import TDPeakRecord # use it by default in TDFileContainer
from pyimgalgos.TDNodeRecord import TDNodeRecord
from pyimgalgos.TDCheetahPeakRecord import TDCheetahPeakRecord
# Initialization
# for peakfinder records
fname = '/reg/neh/home1/dubrovin/LCLS/rel-mengning/work/pfv2-cxif5315-r0169-2015-09-14T14:28:04.txt'
fc = TDFileContainer(fname, indhdr='Evnum', objtype=TDPeakRecord, pbits=0)
# for index table:
fc = TDFileContainer(fname, indhdr='index', objtype=TDNodeRecord)
# for Cheetah file with peaks:
fc = TDFileContainer(fname, indhdr='frameNumber', objtype=TDCheetahPeakRecord)
gr_nums = fc.group_numbers()
ngrps = fc.number_of_groups()
grnum = fc.current_group_number()
gr_curr = fc.group(grpnum) # returns current or specified group
gr_next = fc.next() # returns next group
gr_prev = fc.previous() # returns previous group
hdr = fc.header()
# Print
fc.print_attrs()
fc.print_content(nlines=None) # prints nline (or all by default) lines from file conteiner
# ____________________________________
# Example of iterations over groups
for grnum in fc.group_num_iterator() :
group = fc.next()
group.print_attrs()
peaks = group.get_objs()
for pk in peaks :
pk.print_short()
# Information available through the TDPeakRecord object pk
# ________________________________________________________
# pk.exp, pk.run, pk.evnum, pk.reg
# pk.date, pk.time, pk.tsec, pk.tnsec, pk.fid
# pk.seg, pk.row, pk.col, pk.amax, pk.atot, pk.npix
# pk.rcent, pk.ccent, pk.rsigma, pk.csigma
# pk.rmin, pk.rmax, pk.cmin, pk.cmax
# pk.bkgd, pk.rms, pk.son
# pk.imrow, pk.imcol
# pk.x, pk.y, pk.r, pk.phi
# pk.sonc
# pk.dphi000
# pk.dphi180
# pk.line
# Example of direct access to group by its number
grpnum = 8 # but grpnum is not necessaraly conecutive number, it should be in fc.group_num_iterator() ...
group = fc.group(grpnum) # returns current or specified group
group.print_attrs()
This software was developed for the LCLS project. If you use all or part of it, please give an appropriate acknowledgment.
- See classes:
pyimgalgos.TDFileContainer
- file records container.pyimgalgos.TDGroup
- holds a list of records associated with a single group.pyimgalgos.TDPeakRecord
- provides access to the peak record.pyimgalgos.TDNodeRecord
- provides access to the look-up table with crystal orientation record.pyimgalgos.TDCheetahPeakRecord
- provides access to the Cheetah peak record.
- See:
Created in 2015 by Mikhail Dubrovin
-
class
TDFileContainer.
TDFileContainer
(fname, indhdr='Evnum', objtype=<class pyimgalgos.TDPeakRecord.TDPeakRecord>, pbits=0)[source]¶ Load and hold record list from file and provide access by group index
-
__init__
(fname, indhdr='Evnum', objtype=<class pyimgalgos.TDPeakRecord.TDPeakRecord>, pbits=0)[source]¶ Constructor Args: fname (str) - text table data file name indhdr (str) - header of the field used for group indexing objtype (TD*Recor) - object type used for data record processing/access pbits (int) - print control bit-word; pbits & 256 - tracking
-
TDGroup - text data event information holder/accessor class.
Works together with TDFileContainer and TDPeakRecord classes.
This software was developed for the LCLS project. If you use all or part of it, please give an appropriate acknowledgment.
@see TDFileContainer - loads/holds text data from class and provides per-event-indexed access. @see TDPeakRecord - holds a list of records associated with a single event.
@version $Id$
@author Mikhail S. Dubrovin
-
class
TDGroup.
TDGroup
(recs, objtype, pbits=0)[source]¶ Gets in constructor a list of text data records and converts them in a list of objects
Class TDMatchRecord
helps to retreive and use peak data in processing¶
Usage:
# import
from pyimgalgos.TDMatchRecord import TDMatchRecord
# make object
rec = TDMatchRecord(line)
# access record attributes
index, beta, omega, h, k, l, dr, R, qv, qh, P =
rec.index, rec.beta, rec.omega, rec.h, rec.k, rec.l, rec.dr, rec.R, rec.qv, rec.qh, rec.P
line = rec.line
# print attributes
rec.print_short()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
- See:
Created in 2015 by Mikhail Dubrovin
Class TDNodeRecord
helps to retreive and use peak data in processing¶
Usage:
# Import
from pyimgalgos.TDNodeRecord import TDNodeRecord
# make object
rec = TDNodeRecord(line)
# access record attributes
index, beta, omega, h, k, l, dr, R, qv, qh, qt, ql, P = rec.index, rec.beta, rec.omega, rec.h, rec.k, rec.l, rec.dr, rec.R, rec.qv, rec.qh, rec.qt, rec.ql, rec.P
line = rec.line
# print attributes
rec.print_short()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
- See:
Created in 2015 by Mikhail Dubrovin
Class TDPeakRecord
helps to retreive and use peak data in processing¶
Usage:
# Imports
from pyimgalgos.TDPeakRecord import TDPeakRecord
# Usage
# make object
pk = TDPeakRecord(line)
# access peak attributes
exp = pk.exp # (str) experiment name
run = pk.run # (int) run number
son = pk.son # (float) S/N for pixel with maximal intensity
sonc = pk.sonc # (float) S/N for all pixels included in the peak
line = pk.line # (str) entire record with peak data
...
# Information available through the TDPeakRecord object pk
# ____________________________________________________
# pk.exp, pk.run, pk.evnum, pk.reg
# pk.date, pk.time, pk.tsec, pk.tnsec, pk.fid
# pk.seg, pk.row, pk.col, pk.amax, pk.atot, pk.npix
# pk.rcent, pk.ccent, pk.rsigma, pk.csigma
# pk.rmin, pk.rmax, pk.cmin, pk.cmax
# pk.bkgd, pk.rms, pk.son
# pk.imrow, pk.imcol
# pk.x, pk.y, pk.r, pk.phi
# pk.sonc
# pk.dphi000
# pk.dphi180
# pk.line
# pk.nsplit
# get evaluated parameters
# pk.peak_signal()
# pk.peak_noise()
# pk.peak_son()
# for peak records with fit information:
# pk.fit_phi, pk.fit_beta
# pk.fit_phi_err, pk.fit_beta_err
# pk.fit_chi2, pk.fit_ndof, pk.fit_prob
# print attributes
pk.print_peak_data()
pk.print_peak_data_short()
pk.print_attrs()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
- See:
Created in 2015 by Mikhail Dubrovin
Class TDNodeRecord
helps to retreive and use peak data in processing¶
Usage:
# Imports
from pyimgalgos.TDCheetahPeakRecord import TDCheetahPeakRecord
# Usage
# make object
rec = TDCheetahPeakRecord(line)
# access peak attributes
rec.line
rec.fields
rec.frameNumber
rec.runnum
rec.tstamp
rec.fid
rec.photonEnergyEv
rec.wavelengthA
rec.GMD
rec.peak_index
rec.peak_x_raw
rec.peak_y_raw
rec.peak_r_assembled
rec.peak_q
rec.peak_resA
rec.nPixels
rec.totalIntensity
rec.maxIntensity
rec.sigmaBG
rec.SNR
rec.tsec
# print attributes
rec.print_peak_data()
rec.print_peak_data_short()
rec.print_attrs()
This software was developed for the SIT project. If you use all or part of it, please give an appropriate acknowledgment.
- See:
Created in 2015 by Mikhail Dubrovin