source: sasview/sansdataloader/src/sans/dataloader/manipulations.py @ bd365666

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 bd365666 was f60a8c2, checked in by Mathieu Doucet <doucetm@…>, 12 years ago

Pep-8-ification

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