source: sasview/src/sas/dataloader/manipulations.py @ 400155b

ESS_GUIESS_GUI_DocsESS_GUI_batch_fittingESS_GUI_bumps_abstractionESS_GUI_iss1116ESS_GUI_iss879ESS_GUI_iss959ESS_GUI_openclESS_GUI_orderingESS_GUI_sync_sascalccostrafo411magnetic_scattrelease-4.1.1release-4.1.2release-4.2.2release_4.0.1ticket-1009ticket-1094-headlessticket-1242-2d-resolutionticket-1243ticket-1249ticket885unittest-saveload
Last change on this file since 400155b was 400155b, checked in by gonzalezm, 9 years ago

Implementing request from ticket 261 - default number of bins in Annulus [Phi View] is now 36 and the first bin is now centered at 0 degrees

  • Property mode set to 100644
File size: 40.3 KB
Line 
1"""
2Data manipulations for 2D data sets.
3Using the meta data information, various types of averaging
4are performed in Q-space
5"""
6#####################################################################
7#This software was developed by the University of Tennessee as part of the
8#Distributed Data Analysis of Neutron Scattering Experiments (DANSE)
9#project funded by the US National Science Foundation.
10#See the license text in license.txt
11#copyright 2008, University of Tennessee
12######################################################################
13
14#TODO: copy the meta data from the 2D object to the resulting 1D object
15import math
16import numpy
17
18#from data_info import plottable_2D
19from data_info import Data1D
20
21
22def get_q(dx, dy, det_dist, wavelength):
23    """
24    :param dx: x-distance from beam center [mm]
25    :param dy: y-distance from beam center [mm]
26   
27    :return: q-value at the given position
28    """
29    # Distance from beam center in the plane of detector
30    plane_dist = math.sqrt(dx*dx + dy*dy)
31    # Half of the scattering angle
32    theta = 0.5 * math.atan(plane_dist/det_dist)
33    return (4.0 * math.pi/wavelength) * math.sin(theta)
34
35
36def get_q_compo(dx, dy, det_dist, wavelength, compo=None):
37    """
38    This reduces tiny error at very large q.
39    Implementation of this func is not started yet.<--ToDo
40    """
41    if dy == 0:
42        if dx >= 0:
43            angle_xy = 0
44        else:
45            angle_xy = math.pi
46    else:
47        angle_xy = math.atan(dx/dy)
48       
49    if compo == "x":
50        out = get_q(dx, dy, det_dist, wavelength) * math.cos(angle_xy)
51    elif compo == "y":
52        out = get_q(dx, dy, det_dist, wavelength) * math.sin(angle_xy)
53    else:
54        out = get_q(dx, dy, det_dist, wavelength)
55    return out
56
57
58def flip_phi(phi):
59    """
60    Correct phi to within the 0 <= to <= 2pi range
61   
62    :return: phi in >=0 and <=2Pi
63    """
64    Pi = math.pi
65    if phi < 0:
66        phi_out = phi + (2 * Pi)
67    elif phi > (2 * Pi):
68        phi_out = phi - (2 * Pi)
69    else:
70        phi_out = phi
71    return phi_out
72
73
74def reader2D_converter(data2d=None):
75    """
76    convert old 2d format opened by IhorReader or danse_reader
77    to new Data2D format
78   
79    :param data2d: 2d array of Data2D object
80   
81    :return: 1d arrays of Data2D object
82   
83    """
84    if data2d.data == None or data2d.x_bins == None or data2d.y_bins == None:
85        raise ValueError, "Can't convert this data: data=None..."
86   
87    from sas.dataloader.data_info import Data2D
88
89    new_x = numpy.tile(data2d.x_bins, (len(data2d.y_bins), 1))
90    new_y = numpy.tile(data2d.y_bins, (len(data2d.x_bins), 1))
91    new_y = new_y.swapaxes(0, 1)
92
93    new_data = data2d.data.flatten()
94    qx_data = new_x.flatten()
95    qy_data = new_y.flatten()
96    q_data = numpy.sqrt(qx_data*qx_data + qy_data*qy_data)
97    if data2d.err_data == None or numpy.any(data2d.err_data <= 0):
98        new_err_data = numpy.sqrt(numpy.abs(new_data))
99    else:
100        new_err_data = data2d.err_data.flatten()
101    mask = numpy.ones(len(new_data), dtype=bool)
102
103    #TODO: make sense of the following two lines...
104    output = Data2D()
105    output = data2d
106    output.data = new_data
107    output.err_data = new_err_data
108    output.qx_data = qx_data
109    output.qy_data = qy_data
110    output.q_data = q_data
111    output.mask = mask
112
113    return output
114
115
116class _Slab(object):
117    """
118    Compute average I(Q) for a region of interest
119    """
120    def __init__(self, x_min=0.0, x_max=0.0, y_min=0.0,
121                 y_max=0.0, bin_width=0.001):
122        # Minimum Qx value [A-1]
123        self.x_min = x_min
124        # Maximum Qx value [A-1]
125        self.x_max = x_max
126        # Minimum Qy value [A-1]
127        self.y_min = y_min
128        # Maximum Qy value [A-1]
129        self.y_max = y_max
130        # Bin width (step size) [A-1]
131        self.bin_width = bin_width
132        # If True, I(|Q|) will be return, otherwise,
133        # negative q-values are allowed
134        self.fold = False
135       
136    def __call__(self, data2D):
137        return NotImplemented
138       
139    def _avg(self, data2D, maj):
140        """
141        Compute average I(Q_maj) for a region of interest.
142        The major axis is defined as the axis of Q_maj.
143        The minor axis is the axis that we average over.
144         
145        :param data2D: Data2D object
146        :param maj_min: min value on the major axis
147       
148        :return: Data1D object
149        """
150        if len(data2D.detector) != 1:
151            msg = "_Slab._avg: invalid number of "
152            msg += " detectors: %g" % len(data2D.detector)
153            raise RuntimeError, msg
154       
155        # Get data
156        data = data2D.data[numpy.isfinite(data2D.data)]
157        q_data = data2D.q_data[numpy.isfinite(data2D.data)]
158        err_data = data2D.err_data[numpy.isfinite(data2D.data)]
159        qx_data = data2D.qx_data[numpy.isfinite(data2D.data)]
160        qy_data = data2D.qy_data[numpy.isfinite(data2D.data)]
161             
162        # Build array of Q intervals
163        if maj == 'x':
164            if self.fold:
165                x_min = 0
166            else:
167                x_min = self.x_min
168            nbins = int(math.ceil((self.x_max - x_min) / self.bin_width))
169            qbins = self.bin_width * numpy.arange(nbins) + x_min
170        elif maj == 'y':
171            if self.fold:
172                y_min = 0
173            else:
174                y_min = self.y_min
175            nbins = int(math.ceil((self.y_max - y_min)/self.bin_width))
176            qbins = self.bin_width * numpy.arange(nbins) + y_min
177        else:
178            raise RuntimeError, "_Slab._avg: unrecognized axis %s" % str(maj)
179                               
180        x = numpy.zeros(nbins)
181        y = numpy.zeros(nbins)
182        err_y = numpy.zeros(nbins)
183        y_counts = numpy.zeros(nbins)
184
185        # Average pixelsize in q space
186        for npts in range(len(data)):
187            # default frac
188            frac_x = 0
189            frac_y = 0
190            # get ROI
191            if self.x_min <= qx_data[npts] and self.x_max > qx_data[npts]:
192                frac_x = 1
193            if self.y_min <= qy_data[npts] and self.y_max > qy_data[npts]:
194                frac_y = 1
195            frac = frac_x * frac_y
196           
197            if frac == 0:
198                continue
199            # binning: find axis of q
200            if maj == 'x':
201                q_value = qx_data[npts]
202                min = x_min
203            if maj == 'y':
204                q_value = qy_data[npts]
205                min = y_min
206            if self.fold and q_value < 0:
207                q_value = -q_value
208            # bin
209            i_q = int(math.ceil((q_value - min) / self.bin_width)) - 1
210           
211            # skip outside of max bins
212            if i_q < 0 or i_q >= nbins:
213                continue
214           
215            #TODO: find better definition of x[i_q] based on q_data
216            x[i_q] += frac * q_value  # min + (i_q + 1) * self.bin_width / 2.0
217            y[i_q] += frac * data[npts]
218           
219            if err_data == None or err_data[npts] == 0.0:
220                if data[npts] < 0:
221                    data[npts] = -data[npts]
222                err_y[i_q] += frac * frac * data[npts]
223            else:
224                err_y[i_q] += frac * frac * err_data[npts] * err_data[npts]
225            y_counts[i_q] += frac
226               
227        # Average the sums
228        for n in range(nbins):
229            err_y[n] = math.sqrt(err_y[n])
230         
231        err_y = err_y / y_counts
232        y = y / y_counts
233        x = x / y_counts
234        idx = (numpy.isfinite(y) & numpy.isfinite(x))
235       
236        if not idx.any(): 
237            msg = "Average Error: No points inside ROI to average..."
238            raise ValueError, msg
239        #elif len(y[idx])!= nbins:
240        #    msg = "empty bin(s) due to tight binning..."
241        #    print "resulted",nbins- len(y[idx]), msg
242        return Data1D(x=x[idx], y=y[idx], dy=err_y[idx])
243       
244       
245class SlabY(_Slab):
246    """
247    Compute average I(Qy) for a region of interest
248    """
249    def __call__(self, data2D):
250        """
251        Compute average I(Qy) for a region of interest
252         
253        :param data2D: Data2D object
254       
255        :return: Data1D object
256        """
257        return self._avg(data2D, 'y')
258       
259       
260class SlabX(_Slab):
261    """
262    Compute average I(Qx) for a region of interest
263    """
264    def __call__(self, data2D):
265        """
266        Compute average I(Qx) for a region of interest
267         
268        :param data2D: Data2D object
269       
270        :return: Data1D object
271       
272        """
273        return self._avg(data2D, 'x')
274
275
276class Boxsum(object):
277    """
278    Perform the sum of counts in a 2D region of interest.
279    """
280    def __init__(self, x_min=0.0, x_max=0.0, y_min=0.0, y_max=0.0):
281        # Minimum Qx value [A-1]
282        self.x_min = x_min
283        # Maximum Qx value [A-1]
284        self.x_max = x_max
285        # Minimum Qy value [A-1]
286        self.y_min = y_min
287        # Maximum Qy value [A-1]
288        self.y_max = y_max
289
290    def __call__(self, data2D):
291        """
292        Perform the sum in the region of interest
293         
294        :param data2D: Data2D object
295       
296        :return: number of counts, error on number of counts,
297            number of points summed
298       
299        """
300        y, err_y, y_counts = self._sum(data2D)
301       
302        # Average the sums
303        counts = 0 if y_counts == 0 else y
304        error = 0 if y_counts == 0 else math.sqrt(err_y)
305       
306        # Added y_counts to return, SMK & PDB, 04/03/2013
307        return counts, error, y_counts
308       
309    def _sum(self, data2D):
310        """
311        Perform the sum in the region of interest
312       
313        :param data2D: Data2D object
314       
315        :return: number of counts,
316            error on number of counts, number of entries summed
317       
318        """
319        if len(data2D.detector) != 1:
320            msg = "Circular averaging: invalid number "
321            msg += "of detectors: %g" % len(data2D.detector)
322            raise RuntimeError, msg
323        # Get data
324        data = data2D.data[numpy.isfinite(data2D.data)]
325        q_data = data2D.q_data[numpy.isfinite(data2D.data)]
326        err_data = data2D.err_data[numpy.isfinite(data2D.data)]
327        qx_data = data2D.qx_data[numpy.isfinite(data2D.data)]
328        qy_data = data2D.qy_data[numpy.isfinite(data2D.data)]
329   
330        y = 0.0
331        err_y = 0.0
332        y_counts = 0.0
333
334        # Average pixelsize in q space
335        for npts in range(len(data)):
336            # default frac
337            frac_x = 0
338            frac_y = 0
339           
340            # get min and max at each points
341            qx = qx_data[npts]
342            qy = qy_data[npts]
343           
344            # get the ROI
345            if self.x_min <= qx and self.x_max > qx:
346                frac_x = 1
347            if self.y_min <= qy and self.y_max > qy:
348                frac_y = 1
349            #Find the fraction along each directions
350            frac = frac_x * frac_y
351            if frac == 0:
352                continue
353            y += frac * data[npts]
354            if err_data == None or err_data[npts] == 0.0:
355                if data[npts] < 0:
356                    data[npts] = -data[npts]
357                err_y += frac * frac * data[npts]
358            else:
359                err_y += frac * frac * err_data[npts] * err_data[npts]
360            y_counts += frac
361        return y, err_y, y_counts
362
363
364class Boxavg(Boxsum):
365    """
366    Perform the average of counts in a 2D region of interest.
367    """
368    def __init__(self, x_min=0.0, x_max=0.0, y_min=0.0, y_max=0.0):
369        super(Boxavg, self).__init__(x_min=x_min, x_max=x_max,
370                                      y_min=y_min, y_max=y_max)
371
372    def __call__(self, data2D):
373        """
374        Perform the sum in the region of interest
375         
376        :param data2D: Data2D object
377       
378        :return: average counts, error on average counts
379       
380        """
381        y, err_y, y_counts = self._sum(data2D)
382       
383        # Average the sums
384        counts = 0 if y_counts == 0 else y / y_counts
385        error = 0 if y_counts == 0 else math.sqrt(err_y) / y_counts
386       
387        return counts, error
388       
389       
390def get_pixel_fraction_square(x, xmin, xmax):
391    """
392    Return the fraction of the length
393    from xmin to x.::
394   
395     
396           A            B
397       +-----------+---------+
398       xmin        x         xmax
399     
400    :param x: x-value
401    :param xmin: minimum x for the length considered
402    :param xmax: minimum x for the length considered
403   
404    :return: (x-xmin)/(xmax-xmin) when xmin < x < xmax
405   
406    """
407    if x <= xmin:
408        return 0.0
409    if x > xmin and x < xmax:
410        return (x - xmin) / (xmax - xmin)
411    else:
412        return 1.0
413
414
415class CircularAverage(object):
416    """
417    Perform circular averaging on 2D data
418   
419    The data returned is the distribution of counts
420    as a function of Q
421    """
422    def __init__(self, r_min=0.0, r_max=0.0, bin_width=0.0005):
423        # Minimum radius included in the average [A-1]
424        self.r_min = r_min
425        # Maximum radius included in the average [A-1]
426        self.r_max = r_max
427        # Bin width (step size) [A-1]
428        self.bin_width = bin_width
429
430    def __call__(self, data2D, ismask=False):
431        """
432        Perform circular averaging on the data
433       
434        :param data2D: Data2D object
435       
436        :return: Data1D object
437        """
438        # Get data W/ finite values
439        data = data2D.data[numpy.isfinite(data2D.data)]
440        q_data = data2D.q_data[numpy.isfinite(data2D.data)]
441        qx_data = data2D.qx_data[numpy.isfinite(data2D.data)]
442        err_data = data2D.err_data[numpy.isfinite(data2D.data)]
443        mask_data = data2D.mask[numpy.isfinite(data2D.data)]
444       
445        dq_data = None
446       
447        # Get the dq for resolution averaging
448        if data2D.dqx_data != None and data2D.dqy_data != None:
449            # The pinholes and det. pix contribution present
450            # in both direction of the 2D which must be subtracted when
451            # converting to 1D: dq_overlap should calculated ideally at
452            # q = 0. Note This method works on only pinhole geometry.
453            # Extrapolate dqx(r) and dqy(phi) at q = 0, and take an average.
454            z_max = max(data2D.q_data)
455            z_min = min(data2D.q_data)
456            x_max = data2D.dqx_data[data2D.q_data[z_max]]
457            x_min = data2D.dqx_data[data2D.q_data[z_min]]
458            y_max = data2D.dqy_data[data2D.q_data[z_max]]
459            y_min = data2D.dqy_data[data2D.q_data[z_min]]
460            # Find qdx at q = 0
461            dq_overlap_x = (x_min * z_max - x_max * z_min) / (z_max - z_min)
462            # when extrapolation goes wrong
463            if dq_overlap_x > min(data2D.dqx_data):
464                dq_overlap_x = min(data2D.dqx_data)
465            dq_overlap_x *= dq_overlap_x
466            # Find qdx at q = 0
467            dq_overlap_y = (y_min * z_max - y_max * z_min) / (z_max - z_min)
468            # when extrapolation goes wrong
469            if dq_overlap_y > min(data2D.dqy_data):
470                dq_overlap_y = min(data2D.dqy_data)
471            # get dq at q=0.
472            dq_overlap_y *= dq_overlap_y
473
474            dq_overlap = numpy.sqrt((dq_overlap_x + dq_overlap_y) / 2.0)
475            # Final protection of dq
476            if dq_overlap < 0:
477                dq_overlap = y_min
478            dqx_data = data2D.dqx_data[numpy.isfinite(data2D.data)]
479            dqy_data = data2D.dqy_data[numpy.isfinite(data2D.data)] - dq_overlap
480            # def; dqx_data = dq_r dqy_data = dq_phi
481            # Convert dq 2D to 1D here
482            dqx = dqx_data * dqx_data
483            dqy = dqy_data * dqy_data
484            dq_data = numpy.add(dqx, dqy)
485            dq_data = numpy.sqrt(dq_data)
486           
487        #q_data_max = numpy.max(q_data)
488        if len(data2D.q_data) == None:
489            msg = "Circular averaging: invalid q_data: %g" % data2D.q_data
490            raise RuntimeError, msg
491
492        # Build array of Q intervals
493        nbins = int(math.ceil((self.r_max - self.r_min) / self.bin_width))
494        qbins = self.bin_width * numpy.arange(nbins) + self.r_min
495
496        x = numpy.zeros(nbins)
497        y = numpy.zeros(nbins)
498        err_y = numpy.zeros(nbins)
499        err_x = numpy.zeros(nbins)
500        y_counts = numpy.zeros(nbins)
501
502        for npt in range(len(data)):
503           
504            if ismask and not mask_data[npt]:
505                continue
506           
507            frac = 0
508           
509            # q-value at the pixel (j,i)
510            q_value = q_data[npt]
511            data_n = data[npt]
512           
513            ## No need to calculate the frac when all data are within range
514            if self.r_min >= self.r_max:
515                raise ValueError, "Limit Error: min > max"
516           
517            if self.r_min <= q_value and q_value <= self.r_max:
518                frac = 1
519            if frac == 0:
520                continue 
521            i_q = int(math.floor((q_value - self.r_min) / self.bin_width))
522
523            # Take care of the edge case at phi = 2pi.
524            if i_q == nbins:
525                i_q = nbins - 1
526            y[i_q] += frac * data_n
527            # Take dqs from data to get the q_average
528            x[i_q] += frac * q_value
529            if err_data == None or err_data[npt] == 0.0:
530                if data_n < 0:
531                    data_n = -data_n
532                err_y[i_q] += frac * frac * data_n
533            else:
534                err_y[i_q] += frac * frac * err_data[npt] * err_data[npt]
535            if dq_data != None:
536                # To be consistent with dq calculation in 1d reduction,
537                # we need just the averages (not quadratures) because
538                # it should not depend on the number of the q points
539                # in the qr bins.
540                err_x[i_q] += frac * dq_data[npt]
541            else:
542                err_x = None
543            y_counts[i_q] += frac
544       
545        # Average the sums
546        for n in range(nbins):
547            if err_y[n] < 0:
548                err_y[n] = -err_y[n]
549            err_y[n] = math.sqrt(err_y[n])
550            #if err_x != None:
551            #    err_x[n] = math.sqrt(err_x[n])
552           
553        err_y = err_y / y_counts
554        err_y[err_y == 0] = numpy.average(err_y)
555        y = y / y_counts
556        x = x / y_counts
557        idx = (numpy.isfinite(y)) & (numpy.isfinite(x))
558       
559        if err_x != None:
560            d_x = err_x[idx] / y_counts[idx]
561        else:
562            d_x = None
563
564        if not idx.any():
565            msg = "Average Error: No points inside ROI to average..."
566            raise ValueError, msg
567       
568        return Data1D(x=x[idx], y=y[idx], dy=err_y[idx], dx=d_x)
569   
570
571class Ring(object):
572    """
573    Defines a ring on a 2D data set.
574    The ring is defined by r_min, r_max, and
575    the position of the center of the ring.
576   
577    The data returned is the distribution of counts
578    around the ring as a function of phi.
579   
580    Phi_min and phi_max should be defined between 0 and 2*pi
581    in anti-clockwise starting from the x- axis on the left-hand side
582    """
583    #Todo: remove center.
584    def __init__(self, r_min=0, r_max=0, center_x=0, center_y=0, nbins=36):
585        # Minimum radius
586        self.r_min = r_min
587        # Maximum radius
588        self.r_max = r_max
589        # Center of the ring in x
590        self.center_x = center_x
591        # Center of the ring in y
592        self.center_y = center_y
593        # Number of angular bins
594        self.nbins_phi = nbins
595
596       
597    def __call__(self, data2D):
598        """
599        Apply the ring to the data set.
600        Returns the angular distribution for a given q range
601       
602        :param data2D: Data2D object
603       
604        :return: Data1D object
605        """
606        if data2D.__class__.__name__ not in ["Data2D", "plottable_2D"]:
607            raise RuntimeError, "Ring averaging only take plottable_2D objects"
608       
609        Pi = math.pi
610 
611        # Get data
612        data = data2D.data[numpy.isfinite(data2D.data)]
613        q_data = data2D.q_data[numpy.isfinite(data2D.data)]
614        err_data = data2D.err_data[numpy.isfinite(data2D.data)]
615        qx_data = data2D.qx_data[numpy.isfinite(data2D.data)]
616        qy_data = data2D.qy_data[numpy.isfinite(data2D.data)]
617       
618        q_data_max = numpy.max(q_data)
619       
620        # Set space for 1d outputs
621        phi_bins   = numpy.zeros(self.nbins_phi)
622        phi_counts = numpy.zeros(self.nbins_phi)
623        phi_values = numpy.zeros(self.nbins_phi)
624        phi_err    = numpy.zeros(self.nbins_phi)
625       
626        # Shift to apply to calculated phi values in order to center first bin at zero
627        phi_shift = Pi / self.nbins_phi
628
629        for npt in range(len(data)):
630            frac = 0
631            # q-value at the point (npt)
632            q_value = q_data[npt]
633            data_n = data[npt]
634                       
635            # phi-value at the point (npt)
636            phi_value = math.atan2(qy_data[npt], qx_data[npt]) + Pi
637           
638            if self.r_min <= q_value and q_value <= self.r_max:
639                frac = 1
640            if frac == 0:
641                continue
642            # binning           
643            i_phi = int(math.floor((self.nbins_phi) * (phi_value+phi_shift) / (2 * Pi)))
644           
645            # Take care of the edge case at phi = 2pi.
646            if i_phi >= self.nbins_phi:
647                i_phi =  0
648            phi_bins[i_phi] += frac * data[npt]
649           
650            if err_data == None or err_data[npt] == 0.0:
651                if data_n < 0:
652                    data_n = -data_n
653                phi_err[i_phi] += frac * frac * math.fabs(data_n)
654            else:
655                phi_err[i_phi] += frac * frac * err_data[npt] * err_data[npt]
656            phi_counts[i_phi] += frac
657                         
658        for i in range(self.nbins_phi):
659            phi_bins[i] = phi_bins[i] / phi_counts[i]
660            phi_err[i] = math.sqrt(phi_err[i]) / phi_counts[i]
661            phi_values[i] = 2.0 * math.pi / self.nbins_phi * (1.0 * i)
662           
663        idx = (numpy.isfinite(phi_bins))
664
665        if not idx.any():
666            msg = "Average Error: No points inside ROI to average..."
667            raise ValueError, msg
668        #elif len(phi_bins[idx])!= self.nbins_phi:
669        #    print "resulted",self.nbins_phi- len(phi_bins[idx])
670        #,"empty bin(s) due to tight binning..."
671        return Data1D(x=phi_values[idx], y=phi_bins[idx], dy=phi_err[idx])
672   
673   
674def get_pixel_fraction(qmax, q_00, q_01, q_10, q_11):
675    """
676    Returns the fraction of the pixel defined by
677    the four corners (q_00, q_01, q_10, q_11) that
678    has q < qmax.::
679   
680                q_01                q_11
681        y=1         +--------------+
682                    |              |
683                    |              |
684                    |              |
685        y=0         +--------------+
686                q_00                q_10
687       
688                    x=0            x=1
689   
690    """
691    # y side for x = minx
692    x_0 = get_intercept(qmax, q_00, q_01)
693    # y side for x = maxx
694    x_1 = get_intercept(qmax, q_10, q_11)
695   
696    # x side for y = miny
697    y_0 = get_intercept(qmax, q_00, q_10)
698    # x side for y = maxy
699    y_1 = get_intercept(qmax, q_01, q_11)
700   
701    # surface fraction for a 1x1 pixel
702    frac_max = 0
703   
704    if x_0 and x_1:
705        frac_max = (x_0 + x_1) / 2.0
706    elif y_0 and y_1:
707        frac_max = (y_0 + y_1) / 2.0
708    elif x_0 and y_0:
709        if q_00 < q_10:
710            frac_max = x_0 * y_0 / 2.0
711        else:
712            frac_max = 1.0 - x_0 * y_0 / 2.0
713    elif x_0 and y_1:
714        if q_00 < q_10:
715            frac_max = x_0 * y_1 / 2.0
716        else:
717            frac_max = 1.0 - x_0 * y_1 / 2.0
718    elif x_1 and y_0:
719        if q_00 > q_10:
720            frac_max = x_1 * y_0 / 2.0
721        else:
722            frac_max = 1.0 - x_1 * y_0 / 2.0
723    elif x_1 and y_1:
724        if q_00 < q_10:
725            frac_max = 1.0 - (1.0 - x_1) * (1.0 - y_1) / 2.0
726        else:
727            frac_max = (1.0 - x_1) * (1.0 - y_1) / 2.0
728           
729    # If we make it here, there is no intercept between
730    # this pixel and the constant-q ring. We only need
731    # to know if we have to include it or exclude it.
732    elif (q_00 + q_01 + q_10 + q_11) / 4.0 <  qmax:
733        frac_max = 1.0
734
735    return frac_max
736             
737             
738def get_intercept(q, q_0, q_1):
739    """
740    Returns the fraction of the side at which the
741    q-value intercept the pixel, None otherwise.
742    The values returned is the fraction ON THE SIDE
743    OF THE LOWEST Q. ::
744   
745   
746            A           B
747        +-----------+--------+    <--- pixel size
748        0                    1
749        Q_0 -------- Q ----- Q_1   <--- equivalent Q range
750        if Q_1 > Q_0, A is returned
751        if Q_1 < Q_0, B is returned
752        if Q is outside the range of [Q_0, Q_1], None is returned
753         
754    """
755    if q_1 > q_0:
756        if (q > q_0 and q <= q_1):
757            return (q - q_0) / (q_1 - q_0)
758    else:
759        if (q > q_1 and q <= q_0):
760            return (q - q_1) / (q_0 - q_1)
761    return None
762     
763     
764class _Sector:
765    """
766    Defines a sector region on a 2D data set.
767    The sector is defined by r_min, r_max, phi_min, phi_max,
768    and the position of the center of the ring
769    where phi_min and phi_max are defined by the right
770    and left lines wrt central line
771    and phi_max could be less than phi_min.
772   
773    Phi is defined between 0 and 2*pi in anti-clockwise
774    starting from the x- axis on the left-hand side
775    """
776    def __init__(self, r_min, r_max, phi_min=0, phi_max=2*math.pi, nbins=20):
777        self.r_min = r_min
778        self.r_max = r_max
779        self.phi_min = phi_min
780        self.phi_max = phi_max
781        self.nbins = nbins
782       
783    def _agv(self, data2D, run='phi'):
784        """
785        Perform sector averaging.
786       
787        :param data2D: Data2D object
788        :param run:  define the varying parameter ('phi' , 'q' , or 'q2')
789       
790        :return: Data1D object
791        """
792        if data2D.__class__.__name__ not in ["Data2D", "plottable_2D"]:
793            raise RuntimeError, "Ring averaging only take plottable_2D objects"
794        Pi = math.pi
795
796        # Get the all data & info
797        data = data2D.data[numpy.isfinite(data2D.data)]
798        q_data = data2D.q_data[numpy.isfinite(data2D.data)]
799        err_data = data2D.err_data[numpy.isfinite(data2D.data)]
800        qx_data = data2D.qx_data[numpy.isfinite(data2D.data)]
801        qy_data = data2D.qy_data[numpy.isfinite(data2D.data)]
802        dq_data = None
803           
804        # Get the dq for resolution averaging
805        if data2D.dqx_data != None and data2D.dqy_data != None:
806            # The pinholes and det. pix contribution present
807            # in both direction of the 2D which must be subtracted when
808            # converting to 1D: dq_overlap should calculated ideally at
809            # q = 0.
810            # Extrapolate dqy(perp) at q = 0
811            z_max = max(data2D.q_data)
812            z_min = min(data2D.q_data)
813            x_max = data2D.dqx_data[data2D.q_data[z_max]]
814            x_min = data2D.dqx_data[data2D.q_data[z_min]]
815            y_max = data2D.dqy_data[data2D.q_data[z_max]]
816            y_min = data2D.dqy_data[data2D.q_data[z_min]]
817            # Find qdx at q = 0
818            dq_overlap_x = (x_min * z_max - x_max * z_min) / (z_max - z_min)
819            # when extrapolation goes wrong
820            if dq_overlap_x > min(data2D.dqx_data):
821                dq_overlap_x = min(data2D.dqx_data)
822            dq_overlap_x *= dq_overlap_x
823            # Find qdx at q = 0
824            dq_overlap_y = (y_min * z_max - y_max * z_min) / (z_max - z_min)
825            # when extrapolation goes wrong
826            if dq_overlap_y > min(data2D.dqy_data):
827                dq_overlap_y = min(data2D.dqy_data)
828            # get dq at q=0.
829            dq_overlap_y *= dq_overlap_y
830
831            dq_overlap = numpy.sqrt((dq_overlap_x + dq_overlap_y) / 2.0)
832            if dq_overlap < 0:
833                dq_overlap = y_min
834            dqx_data = data2D.dqx_data[numpy.isfinite(data2D.data)]
835            dqy_data = data2D.dqy_data[numpy.isfinite(data2D.data)] - dq_overlap
836            # def; dqx_data = dq_r dqy_data = dq_phi
837            # Convert dq 2D to 1D here
838            dqx = dqx_data * dqx_data
839            dqy = dqy_data * dqy_data
840            dq_data = numpy.add(dqx, dqy)
841            dq_data = numpy.sqrt(dq_data)
842           
843        #set space for 1d outputs
844        x        = numpy.zeros(self.nbins)
845        y        = numpy.zeros(self.nbins)
846        y_err    = numpy.zeros(self.nbins)
847        x_err    = numpy.zeros(self.nbins)
848        y_counts = numpy.zeros(self.nbins)
849                     
850        # Get the min and max into the region: 0 <= phi < 2Pi
851        phi_min = flip_phi(self.phi_min)
852        phi_max = flip_phi(self.phi_max)
853       
854        q_data_max = numpy.max(q_data)
855                     
856        for n in range(len(data)):
857            frac = 0
858           
859            # q-value at the pixel (j,i)
860            q_value = q_data[n]
861            data_n = data[n]
862           
863            # Is pixel within range?
864            is_in = False
865           
866            # phi-value of the pixel (j,i)
867            phi_value = math.atan2(qy_data[n], qx_data[n]) + Pi
868           
869            ## No need to calculate the frac when all data are within range
870            if self.r_min <= q_value and q_value <= self.r_max:
871                frac = 1
872            if frac == 0:
873                continue
874            #In case of two ROIs (symmetric major and minor regions)(for 'q2')
875            if run.lower() == 'q2':
876                ## For minor sector wing
877                # Calculate the minor wing phis
878                phi_min_minor = flip_phi(phi_min - Pi)
879                phi_max_minor = flip_phi(phi_max - Pi)
880                # Check if phis of the minor ring is within 0 to 2pi
881                if phi_min_minor > phi_max_minor:
882                    is_in = (phi_value > phi_min_minor or \
883                              phi_value < phi_max_minor)
884                else:
885                    is_in = (phi_value > phi_min_minor and \
886                              phi_value < phi_max_minor)
887
888            #For all cases(i.e.,for 'q', 'q2', and 'phi')
889            #Find pixels within ROI
890            if phi_min > phi_max:
891                is_in = is_in or (phi_value > phi_min or \
892                                   phi_value < phi_max)
893            else:
894                is_in = is_in or (phi_value >= phi_min  and \
895                                    phi_value < phi_max)
896           
897            if not is_in:
898                frac = 0
899            if frac == 0:
900                continue
901            # Check which type of averaging we need
902            if run.lower() == 'phi':
903                temp_x = (self.nbins) * (phi_value - self.phi_min)
904                temp_y = (self.phi_max - self.phi_min)
905                i_bin = int(math.floor(temp_x / temp_y))
906            else:
907                temp_x = (self.nbins) * (q_value - self.r_min)
908                temp_y = (self.r_max - self.r_min)
909                i_bin = int(math.floor(temp_x / temp_y))
910
911            # Take care of the edge case at phi = 2pi.
912            if i_bin == self.nbins:
913                i_bin = self.nbins - 1
914               
915            ## Get the total y
916            y[i_bin] += frac * data_n
917            x[i_bin] += frac * q_value
918            if err_data[n] == None or err_data[n] == 0.0:
919                if data_n < 0:
920                    data_n = -data_n
921                y_err[i_bin] += frac * frac * data_n
922            else:
923                y_err[i_bin] += frac * frac * err_data[n] * err_data[n]
924               
925            if dq_data != None:
926                # To be consistent with dq calculation in 1d reduction,
927                # we need just the averages (not quadratures) because
928                # it should not depend on the number of the q points
929                # in the qr bins.
930                x_err[i_bin] += frac * dq_data[n]
931            else:
932                x_err = None
933            y_counts[i_bin] += frac
934   
935        # Organize the results
936        for i in range(self.nbins):
937            y[i] = y[i] / y_counts[i]
938            y_err[i] = math.sqrt(y_err[i]) / y_counts[i]
939
940            # The type of averaging: phi,q2, or q
941            # Calculate x[i]should be at the center of the bin
942            if run.lower() == 'phi':
943                x[i] = (self.phi_max - self.phi_min) / self.nbins * \
944                    (1.0 * i + 0.5) + self.phi_min
945            else:
946                # We take the center of ring area, not radius.
947                # This is more accurate than taking the radial center of ring.
948                #delta_r = (self.r_max - self.r_min) / self.nbins
949                #r_inner = self.r_min + delta_r * i
950                #r_outer = r_inner + delta_r
951                #x[i] = math.sqrt((r_inner * r_inner + r_outer * r_outer) / 2)
952                x[i] = x[i] / y_counts[i]
953        y_err[y_err == 0] = numpy.average(y_err)
954        idx = (numpy.isfinite(y) & numpy.isfinite(y_err))
955        if x_err != None:
956            d_x = x_err[idx] / y_counts[idx]
957        else:
958            d_x = None
959        if not idx.any():
960            msg = "Average Error: No points inside sector of ROI to average..."
961            raise ValueError, msg
962        #elif len(y[idx])!= self.nbins:
963        #    print "resulted",self.nbins- len(y[idx]),
964        #"empty bin(s) due to tight binning..."
965        return Data1D(x=x[idx], y=y[idx], dy=y_err[idx], dx=d_x)
966               
967               
968class SectorPhi(_Sector):
969    """
970    Sector average as a function of phi.
971    I(phi) is return and the data is averaged over Q.
972   
973    A sector is defined by r_min, r_max, phi_min, phi_max.
974    The number of bin in phi also has to be defined.
975    """
976    def __call__(self, data2D):
977        """
978        Perform sector average and return I(phi).
979       
980        :param data2D: Data2D object
981        :return: Data1D object
982        """
983        return self._agv(data2D, 'phi')
984   
985   
986class SectorQ(_Sector):
987    """
988    Sector average as a function of Q for both symatric wings.
989    I(Q) is return and the data is averaged over phi.
990   
991    A sector is defined by r_min, r_max, phi_min, phi_max.
992    r_min, r_max, phi_min, phi_max >0.
993    The number of bin in Q also has to be defined.
994    """
995    def __call__(self, data2D):
996        """
997        Perform sector average and return I(Q).
998       
999        :param data2D: Data2D object
1000       
1001        :return: Data1D object
1002        """
1003        return self._agv(data2D, 'q2')
1004
1005
1006class Ringcut(object):
1007    """
1008    Defines a ring on a 2D data set.
1009    The ring is defined by r_min, r_max, and
1010    the position of the center of the ring.
1011   
1012    The data returned is the region inside the ring
1013   
1014    Phi_min and phi_max should be defined between 0 and 2*pi
1015    in anti-clockwise starting from the x- axis on the left-hand side
1016    """
1017    def __init__(self, r_min=0, r_max=0, center_x=0, center_y=0):
1018        # Minimum radius
1019        self.r_min = r_min
1020        # Maximum radius
1021        self.r_max = r_max
1022        # Center of the ring in x
1023        self.center_x = center_x
1024        # Center of the ring in y
1025        self.center_y = center_y
1026
1027    def __call__(self, data2D):
1028        """
1029        Apply the ring to the data set.
1030        Returns the angular distribution for a given q range
1031       
1032        :param data2D: Data2D object
1033       
1034        :return: index array in the range
1035        """
1036        if data2D.__class__.__name__ not in ["Data2D", "plottable_2D"]:
1037            raise RuntimeError, "Ring cut only take plottable_2D objects"
1038
1039        # Get data
1040        qx_data = data2D.qx_data
1041        qy_data = data2D.qy_data
1042        mask = data2D.mask
1043        q_data = numpy.sqrt(qx_data * qx_data + qy_data * qy_data)
1044
1045        # check whether or not the data point is inside ROI
1046        out = (self.r_min <= q_data) & (self.r_max >= q_data)
1047
1048        return (out)
1049       
1050
1051class Boxcut(object):
1052    """
1053    Find a rectangular 2D region of interest.
1054    """
1055    def __init__(self, x_min=0.0, x_max=0.0, y_min=0.0, y_max=0.0):
1056        # Minimum Qx value [A-1]
1057        self.x_min = x_min
1058        # Maximum Qx value [A-1]
1059        self.x_max = x_max
1060        # Minimum Qy value [A-1]
1061        self.y_min = y_min
1062        # Maximum Qy value [A-1]
1063        self.y_max = y_max
1064
1065    def __call__(self, data2D):
1066        """
1067       Find a rectangular 2D region of interest.
1068         
1069       :param data2D: Data2D object
1070       :return: mask, 1d array (len = len(data))
1071           with Trues where the data points are inside ROI, otherwise False
1072        """
1073        mask = self._find(data2D)
1074       
1075        return mask
1076       
1077    def _find(self, data2D):
1078        """
1079        Find a rectangular 2D region of interest.
1080       
1081        :param data2D: Data2D object
1082       
1083        :return: out, 1d array (length = len(data))
1084           with Trues where the data points are inside ROI, otherwise Falses
1085        """
1086        if data2D.__class__.__name__ not in ["Data2D", "plottable_2D"]:
1087            raise RuntimeError, "Boxcut take only plottable_2D objects"
1088        # Get qx_ and qy_data
1089        qx_data = data2D.qx_data
1090        qy_data = data2D.qy_data
1091        mask = data2D.mask
1092       
1093        # check whether or not the data point is inside ROI
1094        outx = (self.x_min <= qx_data) & (self.x_max > qx_data)
1095        outy = (self.y_min <= qy_data) & (self.y_max > qy_data)
1096
1097        return (outx & outy)
1098
1099
1100class Sectorcut(object):
1101    """
1102    Defines a sector (major + minor) region on a 2D data set.
1103    The sector is defined by phi_min, phi_max,
1104    where phi_min and phi_max are defined by the right
1105    and left lines wrt central line.
1106   
1107    Phi_min and phi_max are given in units of radian
1108    and (phi_max-phi_min) should not be larger than pi
1109    """
1110    def __init__(self, phi_min=0, phi_max=math.pi):
1111        self.phi_min = phi_min
1112        self.phi_max = phi_max
1113             
1114    def __call__(self, data2D):
1115        """
1116        Find a rectangular 2D region of interest.
1117       
1118        :param data2D: Data2D object
1119       
1120        :return: mask, 1d array (len = len(data))
1121       
1122        with Trues where the data points are inside ROI, otherwise False
1123        """
1124        mask = self._find(data2D)
1125       
1126        return mask
1127       
1128    def _find(self, data2D):
1129        """
1130        Find a rectangular 2D region of interest.
1131       
1132        :param data2D: Data2D object
1133       
1134        :return: out, 1d array (length = len(data))
1135       
1136        with Trues where the data points are inside ROI, otherwise Falses
1137        """
1138        if data2D.__class__.__name__ not in ["Data2D", "plottable_2D"]:
1139            raise RuntimeError, "Sectorcut take only plottable_2D objects"
1140        Pi = math.pi
1141        # Get data
1142        qx_data = data2D.qx_data
1143        qy_data = data2D.qy_data
1144        phi_data = numpy.zeros(len(qx_data))
1145
1146        # get phi from data
1147        phi_data = numpy.arctan2(qy_data, qx_data)
1148       
1149        # Get the min and max into the region: -pi <= phi < Pi
1150        phi_min_major = flip_phi(self.phi_min + Pi) - Pi
1151        phi_max_major = flip_phi(self.phi_max + Pi) - Pi
1152        # check for major sector
1153        if phi_min_major > phi_max_major:
1154            out_major = (phi_min_major <= phi_data) + (phi_max_major > phi_data)
1155        else:
1156            out_major = (phi_min_major <= phi_data) & (phi_max_major > phi_data)
1157         
1158        # minor sector
1159        # Get the min and max into the region: -pi <= phi < Pi
1160        phi_min_minor = flip_phi(self.phi_min) - Pi
1161        phi_max_minor = flip_phi(self.phi_max) - Pi
1162             
1163        # check for minor sector
1164        if phi_min_minor > phi_max_minor:
1165            out_minor = (phi_min_minor <= phi_data) + \
1166                            (phi_max_minor >= phi_data)
1167        else:
1168            out_minor = (phi_min_minor <= phi_data) & \
1169                            (phi_max_minor >= phi_data)
1170        out = out_major + out_minor
1171       
1172        return out
Note: See TracBrowser for help on using the repository browser.