source: sasview/calculator/resolution_calculator.py @ 4d1dce4

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 4d1dce4 was 9f5b505, checked in by Jae Cho <jhjcho@…>, 13 years ago

small sigma correction in r-direction

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