source: sasview/calculator/resolution_calculator.py @ 3e41f43

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 3e41f43 was 3be3a80, checked in by Jae Cho <jhjcho@…>, 14 years ago

implemented a resolution estimating tool

  • Property mode set to 100644
File size: 29.0 KB
Line 
1"""
2This object is a small tool to allow user to quickly
3determine the variance in q  from the
4instrumental parameters.
5"""
6from instrument import Sample
7from instrument import Detector
8from instrument import Neutron
9from instrument import Aperture
10# import math stuffs
11from math import pi
12from math import sqrt
13from math import cos
14from math import sin
15from math import atan
16from math import atan2
17from math import pow
18from math import asin
19from math import tan
20
21import numpy
22
23#Plank's constant in cgs unit
24_PLANK_H = 6.62606896E-27
25#Gravitational acc. in cgs unit
26_GRAVITY = 981.0
27
28class ResolutionCalculator(object):
29    """
30    compute resolution in 2D
31    """
32    def __init__(self):
33       
34        # wavelength
35        self.wave = Neutron()
36        # sample
37        self.sample = Sample()
38        # aperture
39        self.aperture = Aperture()
40        # detector
41        self.detector = Detector()
42        # 2d image of the resolution
43        self.image = []
44        # resolutions
45        self.sigma_1 = 0
46        self.sigma_2 = 0
47        self.sigma_1d = 0
48        # q min and max
49        self.qx_min = -0.3
50        self.qx_max = 0.3
51        self.qy_min = -0.3
52        self.qy_max = 0.3
53        # q min and max of the detector
54        self.detector_qx_min = -0.3
55        self.detector_qx_max = 0.3
56        self.detector_qy_min = -0.3
57        self.detector_qy_max = 0.3
58        # plots
59        self.plot = None
60        # instrumental params defaults
61        self.mass = 0
62        self.intensity = 0
63        self.wavelength = 0
64        self.wavelength_spread = 0
65        self.source_aperture_size = []
66        self.source2sample_distance = []
67        self.sample2sample_distance = []
68        self.sample_aperture_size = []
69        self.sample2detector_distance = []
70        self.detector_pix_size = []
71        self.detector_size = []
72        # get all the values of the instrumental parameters
73        self.get_all_instrument_params()
74       
75    def compute_and_plot(self, qx_value, qy_value, qx_min, qx_max, 
76                          qy_min, qy_max, coord = 'polar'):
77        """
78        Compute the resolution
79        : qx_value: x component of q
80        : qy_value: y component of q
81        """
82        # compute 2d resolution
83        _, _, sigma_1, sigma_2 = self.compute(qx_value, qy_value, coord)
84        # make image
85        image = self.get_image(qx_value, qy_value, sigma_1, sigma_2, 
86                               qx_min, qx_max, qy_min, qy_max, coord)
87        # plot image
88        return self.plot_image(image) 
89
90    def compute(self, qx_value, qy_value, coord = 'polar'):
91        """
92        Compute the Q resoltuion in || and + direction of 2D
93        : qx_value: x component of q
94        : qy_value: y component of q
95        """
96        # make sure to update all the variables need.
97        self.get_all_instrument_params()
98        # wavelength
99        lamb = self.wavelength
100   
101        if lamb == 0:
102            msg = "Can't compute the resolution: the wavelength is zero..." 
103            raise RuntimeError, msg
104        # wavelength spread
105        lamb_spread = self.wavelength_spread
106        # Find polar values
107        qr_value, phi = self._get_polar_value(qx_value, qy_value)
108        # vacuum wave transfer
109        knot = 2*pi/lamb
110        # scattering angle theta; always true for plane detector
111        # aligned vertically to the ko direction
112        if qr_value > knot:
113            theta = pi/2
114        else:
115            theta = asin(qr_value/knot)
116        # source aperture size
117        rone = self.source_aperture_size
118        # sample aperture size
119        rtwo = self.sample_aperture_size
120        # detector pixel size
121        rthree = self.detector_pix_size
122        # source to sample(aperture) distance
123        l_ssa = self.source2sample_distance[0]
124        # sample(aperture) to detector distance
125        l_sad = self.sample2detector_distance[0]
126        # sample (aperture) to sample distance
127        l_sas = self.sample2sample_distance[0]
128        # source to sample distance
129        l_one = l_ssa + l_sas
130        # sample to detector distance
131        l_two = l_sad - l_sas
132       
133        # Sample offset correction for l_one and Lp on variance calculation
134        l1_cor = (l_ssa * l_two) / (l_sas + l_two)
135        lp_cor = (l_ssa * l_two) / (l_one + l_two)
136        # the radial distance to the pixel from the center of the detector
137        radius = tan(theta)*l_two
138        #Lp = l_one*l_two/(l_one+l_two)
139        # default polar coordinate
140        comp1 = 'radial'
141        comp2 = 'phi'
142        # in the case of the cartesian coordinate
143        if coord == 'cartesian':
144            comp1 = 'x'
145            comp2 = 'y'
146
147        # sigma in the radial/x direction
148        # for source aperture
149        sigma_1  = self.get_variance(rone, l1_cor, phi, comp1)
150        # for sample apperture
151        sigma_1 += self.get_variance(rtwo, lp_cor, phi, comp1)
152        # for detector pix
153        sigma_1 += self.get_variance(rthree, l_two, phi, comp1)
154        # for gravity term
155        sigma_1 +=  self.get_variance_gravity(l_ssa, l_sad, lamb, lamb_spread, 
156                             phi, comp1, 'on') 
157        # for wavelength spread
158        # reserve for 1d calculation
159        sigma_wave_1 = self.get_variance_wave(radius, l_two, lamb_spread, 
160                                          phi, comp1, 'on')
161        # for 1d
162        variance_1d_1 = sigma_1/2 +sigma_wave_1
163        # normalize
164        variance_1d_1 = knot*knot*variance_1d_1/12
165       
166        # for 2d
167        sigma_1 += sigma_wave_1
168        # normalize
169        sigma_1 = knot*sqrt(sigma_1/12)
170
171        # sigma in the phi/y direction
172        # for source apperture
173        sigma_2  = self.get_variance(rone, l1_cor, phi, comp2)
174        # for sample apperture
175        sigma_2 += self.get_variance(rtwo, lp_cor, phi, comp2)
176        # for detector pix
177        sigma_2 += self.get_variance(rthree, l_two, phi, comp2)
178        # for gravity term
179        sigma_2 +=  self.get_variance_gravity(l_ssa, l_sad, lamb, lamb_spread, 
180                             phi, comp2, 'on')
181        # for wavelength spread
182        # reserve for 1d calculation
183        sigma_wave_2 = self.get_variance_wave(radius, l_two, lamb_spread, 
184                                          phi, comp2, 'on') 
185        # for 1d
186        variance_1d_2 = sigma_2/2 +sigma_wave_2
187        # normalize
188        variance_1d_2 = knot*knot*variance_1d_2/12
189       
190        # for 2d
191        sigma_2 += sigma_wave_2
192        # normalize
193        sigma_2 =  knot*sqrt(sigma_2/12)
194
195        # set sigmas
196        self.sigma_1 = sigma_1
197        self.sigma_2 = sigma_2
198       
199        self.sigma_1d = sqrt(variance_1d_1 + variance_1d_2)
200        return qr_value, phi, sigma_1, sigma_2
201   
202    def get_image(self, qx_value, qy_value, sigma_1, sigma_2,
203                  qx_min, qx_max, qy_min, qy_max, coord = 'polar'): 
204        """
205        Get the resolution in polar coordinate ready to plot
206        : qx_value: qx_value value
207        : qy_value: qy_value value
208        : sigma_1: variance in r direction
209        : sigma_2: variance in phi direction
210        : coord: coordinate system of image, 'polar' or 'cartesian'
211        """
212        # Get  qx_max and qy_max...
213        output = self._get_detector_qxqy_pixels()
214        # Set qx_value/qy_value min/max
215        #qx_min = self.qx_min
216        #qx_max = self.qx_max
217        #qy_min = self.qy_min
218        #qy_max = self.qy_max
219
220        # Find polar values
221        qr_value, phi = self._get_polar_value(qx_value, qy_value)
222       
223        # Check whether the q value is within the detector range
224        msg = "Invalid input: Q value out of the detector range..."
225        if qx_min < self.qx_min:
226            self.qx_min = qx_min
227            #raise ValueError, msg
228        if qx_max > self.qx_max:
229            self.qx_max = qx_max
230            #raise ValueError, msg
231        if qy_min < self.qy_min:
232            self.qy_min = qy_min
233            #raise ValueError, msg
234        if qy_max > self.qy_max:
235            self.qy_max = qy_max
236            #raise ValueError, msg
237
238        # Make an empty graph in the detector scale
239        dx_size = (self.qx_max - self.qx_min) / (1000 - 1)
240        dy_size = (self.qy_max - self.qy_min) / (1000 - 1)
241        x_val = numpy.arange(self.qx_min, self.qx_max, dx_size)
242        y_val = numpy.arange(self.qy_max, self.qy_min, -dy_size)
243        q_1, q_2 = numpy.meshgrid(x_val, y_val)
244
245        # check whether polar or cartesian
246        if coord == 'polar':
247            q_1, q_2 = self._rotate_z(q_1, q_2, phi)
248            qc_1 = qr_value
249            qc_2 = 0.0
250        else:
251            # catesian coordinate
252            # qx_center
253            qc_1 = qx_value
254            # qy_center
255            qc_2 = qy_value
256           
257        # Calculate the 2D Gaussian distribution image
258        image = self._gaussian2d(q_1, q_2, qc_1, qc_2, sigma_1, sigma_2)
259        # Add it if there are more than one inputs.
260        if len(self.image) > 0:
261            self.image += image
262        else:
263            self.image = image
264
265        return self.image
266   
267    def plot_image(self, image):
268        """
269        Plot image using pyplot
270        : image: 2d resolution image
271       
272        : return plt: pylab object
273        """
274        import matplotlib.pyplot as plt
275
276        self.plot = plt
277        plt.xlabel('$\\rm{Q}_{x} [A^{-1}]$')
278        plt.ylabel('$\\rm{Q}_{y} [A^{-1}]$')
279        # Max value of the image
280        max = numpy.max(image)
281        # Image
282        im = plt.imshow(image, 
283                extent = [self.qx_min, self.qx_max, self.qy_min, self.qy_max])
284
285        # bilinear interpolation to make it smoother
286        im.set_interpolation('bilinear')
287
288        return plt
289   
290    def reset_image(self):
291        """
292        Reset image to default (=[])
293        """
294        self.image = []
295       
296    def get_variance(self, size = [], distance = 0, phi = 0, comp = 'radial'):
297        """
298        Get the variance when the slit/pinhole size is given
299        : size: list that can be one(diameter for circular)
300                or two components(lengths for rectangular)
301        : distance: [z, x] where z along the incident beam, x // qx_value
302        : comp: direction of the sigma; can be 'phi', 'y', 'x', and 'radial'
303       
304        : return variance: sigma^2
305        """   
306        # check the length of size (list)
307        len_size = len(size)
308       
309        # define sigma component direction
310        if comp == 'radial':
311            phi_x = cos(phi)
312            phi_y = sin(phi)
313        elif comp == 'phi':
314            phi_x = sin(phi)
315            phi_y = cos(phi)
316        elif comp == 'x':
317            phi_x = 1
318            phi_y = 0
319        elif comp == 'y':
320            phi_x = 0
321            phi_y = 1
322        else:
323            phi_x = 0
324            phi_y = 0
325        # calculate each component
326        # for pinhole w/ radius = size[0]/2
327        if len_size == 1:
328            x_comp = (0.5 * size[0]) * sqrt(3)
329            y_comp = 0
330        # for rectangular slit
331        elif len_size == 2:
332            x_comp = size[0] * phi_x
333            y_comp = size[1] * phi_y
334        # otherwise
335        else:
336            raise ValueError, " Improper input..."
337        # get them squared 
338        sigma  = x_comp * x_comp
339        sigma += y_comp * y_comp
340        # normalize by distance
341        sigma /= (distance * distance)
342
343        return sigma
344
345    def get_variance_wave(self, radius, distance, spread, phi, 
346                          comp = 'radial', switch = 'on'):
347        """
348        Get the variance when the wavelength spread is given
349       
350        : radius: the radial distance from the beam center to the pix of q
351        : distance: sample to detector distance
352        : spread: wavelength spread (ratio)
353        : comp: direction of the sigma; can be 'phi', 'y', 'x', and 'radial'
354       
355        : return variance: sigma^2
356        """
357        if switch.lower() == 'off':
358            return 0
359        # check the singular point
360        if distance == 0 or comp == 'phi':
361            return 0
362        else:
363            # calculate sigma^2
364            sigma = 2 * pow(radius/distance*spread, 2)
365            if comp == 'x':
366                sigma *= (cos(phi)*cos(phi))
367            elif comp == 'y':
368                sigma *= (sin(phi)*sin(phi))
369            else:
370                sigma *= 1         
371               
372            return sigma
373
374    def get_variance_gravity(self, s_distance, d_distance, wavelength, spread, 
375                             phi, comp = 'radial', switch = 'on'):
376        """
377        Get the variance from gravity when the wavelength spread is given
378       
379        : s_distance: source to sample distance
380        : d_distance: sample to detector distance
381        : wavelength: wavelength
382        : spread: wavelength spread (ratio)
383        : comp: direction of the sigma; can be 'phi', 'y', 'x', and 'radial'
384       
385        : return variance: sigma^2
386        """
387        if switch.lower() == 'off':
388            return 0
389        # check the singular point
390        if d_distance == 0 or comp == 'x':
391            return 0
392        else:
393            # neutron mass in cgs unit
394            self.mass = self.get_neutron_mass()
395            # plank constant in cgs unit
396            h_constant = _PLANK_H
397            # gravity in cgs unit
398            gravy = _GRAVITY
399            # m/h
400            m_over_h = self.mass /h_constant
401            # A value
402            a_value = d_distance * (s_distance + d_distance)
403            a_value *= pow(m_over_h / 2, 2)
404            a_value *= gravy
405            # unit correction (1/cm to 1/A) for A and d_distance below
406            a_value *= 1.0E-16
407           
408            # calculate sigma^2
409            sigma = pow(a_value / d_distance, 2)
410            sigma *= pow(wavelength, 4)
411            sigma *= pow(spread, 2)
412            sigma *= 8
413           
414            # only for the polar coordinate
415            if comp == 'radial':
416                sigma *= (sin(phi) * sin(phi))
417            elif comp == 'phi':
418                sigma *= (cos(phi) * cos(phi))
419           
420            return sigma
421       
422    def get_intensity(self):
423        """
424        Get intensity
425        """
426        return self.wave.intensity
427
428    def get_wavelength(self):
429        """
430        Get wavelength
431        """
432        return self.wave.wavelength
433   
434    def get_wavelength_spread(self):
435        """
436        Get wavelength spread
437        """
438        return self.wave.wavelength_spread
439   
440    def get_neutron_mass(self):
441        """
442        Get Neutron mass
443        """
444        return self.wave.mass
445   
446    def get_source_aperture_size(self):
447        """
448        Get source aperture size
449        """
450        return self.aperture.source_size
451   
452    def get_sample_aperture_size(self):
453        """
454        Get sample aperture size
455        """
456        return self.aperture.sample_size
457   
458    def get_detector_pix_size(self):
459        """
460        Get detector pixel size
461        """
462        return self.detector.pix_size
463   
464    def get_detector_size(self):
465        """
466        Get detector size
467        """
468        return self.detector.size   
469     
470    def get_source2sample_distance(self):
471        """
472        Get detector source2sample_distance
473        """
474        return self.aperture.sample_distance
475   
476    def get_sample2sample_distance(self):
477        """
478        Get detector sampleslitsample_distance
479        """
480        return self.sample.distance
481   
482    def get_sample2detector_distance(self):
483        """
484        Get detector sample2detector_distance
485        """
486        return self.detector.distance
487   
488    def set_intensity(self, intensity):
489        """
490        Set intensity
491        """
492        self.wave.set_intensity(intensity)
493   
494    def set_wavelength(self, wavelength):
495        """
496        Set wavelength
497        """
498        self.wave.set_wavelength(wavelength)
499   
500    def set_wavelength_spread(self, wavelength_spread):
501        """
502        Set wavelength spread
503        """
504        self.wave.set_wavelength_spread(wavelength_spread)
505   
506    def set_source_aperture_size(self, size):
507        """
508        Set source aperture size
509       
510        : param size: [dia_value] or [x_value, y_value]
511        """       
512        if len(size) < 1 or len(size) > 2:
513            raise RuntimeError, "The length of the size must be one or two."
514        self.aperture.set_source_size(size)
515       
516    def set_neutron_mass(self, mass):
517        """
518        Set Neutron mass
519        """
520        self.wave.set_mass(mass)
521       
522    def set_sample_aperture_size(self, size):
523        """
524        Set sample aperture size
525       
526        : param size: [dia_value] or [xheight_value, yheight_value]
527        """
528        if len(size) < 1 or len(size) > 2:
529            raise RuntimeError, "The length of the size must be one or two."
530        self.aperture.set_sample_size(size)
531   
532    def set_detector_pix_size(self, size):
533        """
534        Set detector pixel size
535        """
536        self.detector.set_pix_size(size)
537       
538    def set_detector_size(self, size):
539        """
540        Set detector size in number of pixels
541        : param size: [pixel_nums] or [x_pix_num, yx_pix_num]
542        """
543        self.detector.set_size(size)
544       
545    def set_source2sample_distance(self, distance):
546        """
547        Set detector source2sample_distance
548       
549        : param distance: [distance, x_offset]
550        """
551        if len(distance) < 1 or len(distance) > 2:
552            raise RuntimeError, "The length of the size must be one or two."
553        self.aperture.set_sample_distance(distance)
554
555    def set_sample2sample_distance(self, distance):
556        """
557        Set detector sample_slit2sample_distance
558       
559        : param distance: [distance, x_offset]
560        """
561        if len(distance) < 1 or len(distance) > 2:
562            raise RuntimeError, "The length of the size must be one or two."
563        self.sample.set_distance(distance)
564   
565    def set_sample2detector_distance(self, distance):
566        """
567        Set detector sample2detector_distance
568       
569        : param distance: [distance, x_offset]
570        """
571        if len(distance) < 1 or len(distance) > 2:
572            raise RuntimeError, "The length of the size must be one or two."
573        self.detector.set_distance(distance)
574       
575    def get_all_instrument_params(self):
576        """
577        Get all instrumental parameters
578        """
579        self.intensity = self.get_intensity()
580        self.wavelength = self.get_wavelength()
581        self.wavelength_spread = self.get_wavelength_spread()
582        self.mass = self.get_neutron_mass()
583        self.source_aperture_size = self.get_source_aperture_size()
584        self.sample_aperture_size = self.get_sample_aperture_size()
585        self.detector_pix_size = self.get_detector_pix_size()
586        self.detector_size = self.get_detector_size()
587        self.source2sample_distance = self.get_source2sample_distance()
588        self.sample2sample_distance = self.get_sample2sample_distance()
589        self.sample2detector_distance = self.get_sample2detector_distance()
590
591
592         
593    def _rotate_z(self, x_value, y_value, theta= 0.0):
594        """
595        Rotate x-y cordinate around z-axis by theta
596        : x_value: numpy array of x values
597        : y_value: numpy array of y values
598        : theta: angle to rotate by in rad
599       
600        :return: x_prime, y-prime
601        """       
602        # rotate by theta
603        x_prime = x_value * cos(theta) + y_value * sin(theta)
604        y_prime =  -x_value * sin(theta) + y_value * cos(theta)
605   
606        return x_prime, y_prime
607   
608    def _gaussian2d(self, x_val, y_val, x0_val, y0_val, sigma_x, sigma_y):
609        """
610        Calculate 2D Gaussian distribution
611        : x_val: x value
612        : y_val: y value
613        : x0_val: mean value in x-axis
614        : y0_val: mean value in y-axis
615        : sigma_x: variance in x-direction
616        : sigma_y: variance in y-direction
617       
618        : return: gaussian (value)
619        """
620        # call gaussian1d
621        gaussian  = self._gaussian1d(x_val, x0_val, sigma_x)
622        gaussian *= self._gaussian1d(y_val, y0_val, sigma_y)
623 
624        # normalizing factor correction
625        if sigma_x != 0 and sigma_y != 0:
626            gaussian *= sqrt(2 * pi)
627        return gaussian
628
629    def _gaussian1d(self, value, mean, sigma):
630        """
631        Calculate 1D Gaussian distribution
632        : value: value
633        : mean: mean value
634        : sigma: variance
635       
636        : return: gaussian (value)
637        """
638        # default
639        gaussian = 1.0
640        if sigma != 0:
641            # get exponent
642            nu_value = (value - mean) / sigma
643            nu_value *= nu_value
644            nu_value *= -0.5
645            gaussian *= numpy.exp(nu_value)
646            gaussian /= sigma
647            # normalize
648            gaussian /= sqrt(2 * pi)
649           
650        return gaussian
651   
652    def _atan_phi(self, qy_value, qx_value):
653        """
654        Find the angle phi of q on the detector plane for qx_value, qy_value given
655        : qx_value: x component of q
656        : qy_value: y component of q
657       
658        : return phi: the azimuthal angle of q on x-y plane
659        """
660        # default
661        phi = 0
662        # ToDo: This is misterious - sign???
663        #qy_value = -qy_value
664        # Take care of the singular point
665        if qx_value == 0:
666            if qy_value > 0:
667                phi = pi / 2
668            elif qy_value < 0:
669                phi = -pi / 2
670            else:
671                phi = 0
672        else:
673            # the angle
674            phi = atan2(qy_value, qx_value)
675
676        return phi
677
678    def _get_detector_qxqy_pixels(self):
679        """
680        Get the pixel positions of the detector in the qx_value-qy_value space
681        """
682       
683        # update all param values
684        self.get_all_instrument_params()
685       
686        # wavelength
687        wavelength = self.wavelength
688        # Gavity correction
689        delta_y = self._get_beamcenter_drop() # in cm
690       
691        # detector_pix size
692        detector_pix_size = self.detector_pix_size
693        # Square or circular pixel
694        if len(detector_pix_size) == 1:
695            pix_x_size = detector_pix_size[0]
696            pix_y_size = detector_pix_size[0]
697        # rectangular pixel pixel
698        elif len(detector_pix_size) == 2:
699            pix_x_size = detector_pix_size[0]
700            pix_y_size = detector_pix_size[1]
701        else:
702            raise ValueError, " Input value format error..."
703        # Sample to detector distance = sample slit to detector
704        # minus sample offset
705        sample2detector_distance = self.sample2detector_distance[0] - \
706                                    self.sample2sample_distance[0]
707        # detector offset in x-direction
708        detector_offset = 0
709        try:
710            detector_offset = self.sample2detector_distance[1]
711        except:
712            pass
713       
714        # detector size in [no of pix_x,no of pix_y]
715        detector_pix_nums_x = self.detector_size[0]
716       
717        # get pix_y if it exists, otherwse take it from [0]
718        try:
719            detector_pix_nums_y = self.detector_size[1]
720        except:
721            detector_pix_nums_y = self.detector_size[0]
722       
723        # detector offset in pix number
724        offset_x = detector_offset / pix_x_size
725        offset_y = delta_y / pix_y_size
726       
727        # beam center position in pix number (start from 0)
728        center_x, center_y = self._get_beamcenter_position(detector_pix_nums_x,
729                                    detector_pix_nums_y, offset_x, offset_y)
730        # distance [cm] from the beam center on detector plane
731        detector_ind_x = numpy.arange(detector_pix_nums_x)
732        detector_ind_y = numpy.arange(detector_pix_nums_y)
733
734        # shif 0.5 pixel so that pix position is at the center of the pixel
735        detector_ind_x = detector_ind_x + 0.5
736        detector_ind_y = detector_ind_y + 0.5
737
738        # the relative postion from the beam center
739        detector_ind_x = detector_ind_x - center_x
740        detector_ind_y = detector_ind_y - center_y
741       
742        # unit correction in cm
743        detector_ind_x = detector_ind_x * pix_x_size
744        detector_ind_y = detector_ind_y * pix_y_size
745       
746        qx_value = numpy.zeros(len(detector_ind_x))
747        qy_value = numpy.zeros(len(detector_ind_y))
748        i = 0
749
750        for indx in detector_ind_x:
751            qx_value[i] = self._get_qx(indx, sample2detector_distance, wavelength)
752            i += 1
753        i = 0
754        for indy in detector_ind_y:
755            qy_value[i] = self._get_qx(indy, sample2detector_distance, wavelength)
756            i += 1
757           
758        # qx_value and qy_value values in array
759        qx_value = qx_value.repeat(detector_pix_nums_y) 
760        qx_value = qx_value.reshape(detector_pix_nums_x, detector_pix_nums_y)
761        qy_value = qy_value.repeat(detector_pix_nums_x) 
762        qy_value = qy_value.reshape(detector_pix_nums_y, detector_pix_nums_x)
763        qy_value = qy_value.transpose()
764
765        # p min and max values among the center of pixels
766        self.qx_min = numpy.min(qx_value)
767        self.qx_max = numpy.max(qx_value)
768        self.qy_min = numpy.min(qy_value)
769        self.qy_max = numpy.max(qy_value)
770               
771        # Appr. min and max values of the detector display limits
772        # i.e., edges of the last pixels.
773        self.qy_min += self._get_qx(-0.5 * pix_y_size, 
774                                sample2detector_distance, wavelength)
775        self.qy_max += self._get_qx(0.5 * pix_y_size, 
776                                sample2detector_distance, wavelength)
777        #if self.qx_min == self.qx_max:
778        self.qx_min += self._get_qx(-0.5 * pix_x_size, 
779                                sample2detector_distance, wavelength)
780        self.qx_max += self._get_qx(0.5 * pix_x_size, 
781                                    sample2detector_distance, wavelength)
782        # min and max values of detecter
783        self.detector_qx_min = self.qx_min
784        self.detector_qx_max = self.qx_max
785        self.detector_qy_min = self.qy_min
786        self.detector_qy_max = self.qy_max
787
788        # try to set it as a Data2D otherwise pass (not required for now)
789        try:
790            from DataLoader.data_info import Data2D
791            output = Data2D()
792            inten = numpy.zeros_like(qx_value)
793            output.data     = inten
794            output.qx_data  = qx_value
795            output.qy_data  = qy_value           
796        except:
797            pass
798       
799        return output#qx_value,qy_value
800       
801    def _get_qx(self, dx_size, det_dist, wavelength):
802        """
803        :param dx_size: x-distance from beam center [cm]
804        :param det_dist: sample to detector distance [cm]
805       
806        :return: q-value at the given position
807        """
808        # Distance from beam center in the plane of detector
809        plane_dist = dx_size
810        # full scattering angle on the x-axis
811        theta  = atan(plane_dist / det_dist)
812        qx_value     = (2.0 * pi / wavelength) * sin(theta)
813        return qx_value 
814   
815    def _get_polar_value(self, qx_value, qy_value):
816        """
817        Find qr_value and phi from qx_value and qy_value values
818       
819        : return qr_value, phi
820        """
821        # find |q| on detector plane
822        qr_value = sqrt(qx_value*qx_value + qy_value*qy_value)
823        # find angle phi
824        phi = self._atan_phi(qy_value, qx_value)
825       
826        return qr_value, phi
827   
828    def _get_beamcenter_position(self, num_x, num_y, offset_x, offset_y):
829        """
830        :param num_x: number of pixel in x-direction
831        :param num_y: number of pixel in y-direction
832        :param offset: detector offset in x-direction in pix number
833       
834        :return: pix number; pos_x, pos_y in pix index
835        """
836        # beam center position
837        pos_x = num_x / 2
838        pos_y = num_y / 2
839
840        # correction for offset
841        pos_x += offset_x
842        # correction for gravity that is always negative
843        pos_y -= offset_y
844
845        return pos_x, pos_y     
846
847    def _get_beamcenter_drop(self):
848        """
849        Get the beam center drop (delta y) in y diection due to gravity
850       
851        :return delta y: the beam center drop in cm
852        """
853        # Covert unit from A to cm
854        unit_cm = 1e-08
855        # Velocity of neutron in horizontal direction (~ actual velocity)
856        velocity = _PLANK_H / (self.mass * self.wavelength * unit_cm)
857        # Compute delta y
858        delta_y = 0.5
859        delta_y *= _GRAVITY
860        delta_y *= self.sample2detector_distance[0] 
861        delta_y *= (self.source2sample_distance[0] + self.sample2detector_distance[0])
862        delta_y /= (velocity * velocity)
863
864        return delta_y     
865
866   
Note: See TracBrowser for help on using the repository browser.