source: sasview/sanscalculator/src/sans/calculator/resolution_calculator.py @ 9a945a2

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 9a945a2 was 5f3164c, checked in by Jae Cho <jhjcho@…>, 13 years ago

updated the gravitational effect based on new calculation and fixed lines of detector limits

  • Property mode set to 100644
File size: 41.7 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 TOF as 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        self.image_lam = []
39        # resolutions
40        # lamda in r-direction
41        self.sigma_lamd = 0
42        # x-dir (no lamda)
43        self.sigma_1 = 0
44        #y-dir (no lamda)
45        self.sigma_2 = 0
46        # 1D total
47        self.sigma_1d = 0
48        self.gravity_phi = None
49        # q min and max
50        self.qx_min = -0.3
51        self.qx_max = 0.3
52        self.qy_min = -0.3
53        self.qy_max = 0.3
54        # q min and max of the detector
55        self.detector_qx_min = -0.3
56        self.detector_qx_max = 0.3
57        self.detector_qy_min = -0.3
58        self.detector_qy_max = 0.3
59        # possible max qrange
60        self.qxmin_limit = 0
61        self.qxmax_limit = 0
62        self.qymin_limit = 0
63        self.qymax_limit = 0
64
65        # plots
66        self.plot = None
67        # instrumental params defaults
68        self.mass = 0
69        self.intensity = 0
70        self.wavelength = 0
71        self.wavelength_spread = 0
72        self.source_aperture_size = []
73        self.source2sample_distance = []
74        self.sample2sample_distance = []
75        self.sample_aperture_size = []
76        self.sample2detector_distance = []
77        self.detector_pix_size = []
78        self.detector_size = []
79        # get all the values of the instrumental parameters
80        #self.intensity = self.get_intensity()
81        #self.wavelength = self.get_wavelength()
82        #self.wavelength_spread = self.get_wavelength_spread()
83        self.get_all_instrument_params()
84        # max q range for all lambdas
85        self.qxrange = []
86        self.qyrange = []
87       
88    def compute_and_plot(self, qx_value, qy_value, qx_min, qx_max, 
89                          qy_min, qy_max, coord = 'cartesian'):
90        """
91        Compute the resolution
92        : qx_value: x component of q
93        : qy_value: y component of q
94        """
95        # make sure to update all the variables need.
96        # except lambda, dlambda, and intensity
97        self.get_all_instrument_params()
98        # wavelength etc.
99        lamda_list, dlamb_list = self.get_wave_list()
100        intens_list = []#self.get_intensity_list()
101
102        sig1_list = []
103        sig2_list = []
104        sigr_list = []
105        sigma1d_list = []
106        num_lamda = len(lamda_list)
107        for num in range(num_lamda):
108            lam = lamda_list[num]
109            # wavelength spread
110            dlam = dlamb_list[num]
111            intens = self.setup_tof(lam, dlam)
112            intens_list.append(intens)
113            # cehck if tof
114            if num_lamda > 1:
115                tof = True
116            else:
117                tof = False
118            # compute 2d resolution
119            _, _, sigma_1, sigma_2, sigma_r, sigma1d = \
120                        self.compute(lam, dlam, qx_value, qy_value, coord, tof)
121            # make image
122            image = self.get_image(qx_value, qy_value, sigma_1, sigma_2, 
123                            sigma_r, qx_min, qx_max, qy_min, qy_max,
124                            coord, False)
125           
126            # Non tof mode to be speed up
127            #if num_lamda < 2:
128            #    return self.plot_image(image)
129           
130            if qx_min > self.qx_min:
131                qx_min = self.qx_min
132            if qx_max < self.qx_max:
133                qx_max = self.qx_max
134            if qy_min > self.qy_min:
135                qy_min = self.qy_min
136            if qy_max < self.qy_max:
137                qy_max = self.qy_max
138               
139            # set max qranges
140            self.qxrange = [qx_min, qx_max]
141            self.qyrange = [qy_min, qy_max]
142            #print qy_max+qy_min,qy_max,qy_min
143            sig1_list.append(sigma_1)
144            sig2_list.append(sigma_2)
145            sigr_list.append(sigma_r)
146            sigma1d_list.append(sigma1d)
147        # redraw image in global 2d q-space.   
148        self.image_lam = []
149        total_intensity = 0
150        sigma_1 = 0
151        sigma_r = 0
152        sigma_2 = 0
153        sigma1d = 0
154        for ind in range(num_lamda):
155            lam = lamda_list[ind]
156            dlam = dlamb_list[ind]
157            intens = self.setup_tof(lam, dlam)
158            out = self.get_image(qx_value, qy_value, sig1_list[ind], 
159                                   sig2_list[ind], sigr_list[ind], 
160                                   qx_min, qx_max, qy_min, qy_max, coord)
161            # this is the case of q being outside the detector
162            #if numpy.all(out==0.0):
163            #    continue
164            image = out
165            # set variance as sigmas
166            sigma_1 += sig1_list[ind] *  sig1_list[ind] * self.intensity
167            sigma_r += sigr_list[ind] * sigr_list[ind] * self.intensity
168            sigma_2 += sig2_list[ind] * sig2_list[ind] * self.intensity
169            sigma1d += sigma1d_list[ind] * sigma1d_list[ind] * self.intensity
170            total_intensity += self.intensity
171       
172        if total_intensity != 0:
173            # average variance
174            image_out = image / total_intensity
175            sigma_1 = sigma_1 / total_intensity
176            sigma_r = sigma_r / total_intensity
177            sigma_2 = sigma_2 / total_intensity
178            sigma1d = sigma1d / total_intensity
179            # set sigmas
180            self.sigma_1 = sqrt(sigma_1)
181            self.sigma_lamd = sqrt(sigma_r)
182            self.sigma_2 = sqrt(sigma_2)
183            self.sigma_1d = sqrt(sigma1d)
184            # rescale
185            max_im_val = 1 #image_out.max()
186            if max_im_val > 0:
187                image_out /= max_im_val
188        else:
189            image_out = image * 0.0
190            # Don't calculate sigmas nor set self.sigmas!
191            sigma_1 = 0
192            sigma_r = 0
193            sigma_2 = 0
194            sigma1d = 0
195        if len(self.image) > 0:
196            self.image += image_out
197        else:
198            self.image = image_out
199       
200        # plot image
201        return self.plot_image(self.image) 
202   
203    def setup_tof(self, wavelength, wavelength_spread):
204        """
205        Setup all parameters in instrument
206       
207        : param ind: index of lambda, etc
208        """
209
210        # set wave.wavelength
211        self.set_wavelength(wavelength)
212        self.set_wavelength_spread(wavelength_spread)
213        self.intensity = self.wave.get_intensity()
214       
215        if wavelength == 0:
216            msg = "Can't compute the resolution: the wavelength is zero..." 
217            raise RuntimeError, msg
218        return self.intensity
219       
220    def compute(self, wavelength, wavelength_spread, qx_value, qy_value, 
221                coord = 'cartesian', tof=False):
222        """
223        Compute the Q resoltuion in || and + direction of 2D
224        : qx_value: x component of q
225        : qy_value: y component of q
226        """
227        coord = 'cartesian'
228        lamb = wavelength
229        lamb_spread = wavelength_spread
230        # the shape of wavelength distribution
231
232        if tof:
233            # rectangular
234            tof_factor =  2
235        else:
236            # triangular
237            tof_factor =  1
238        # Find polar values
239        qr_value, phi = self._get_polar_value(qx_value, qy_value)
240        # vacuum wave transfer
241        knot = 2*pi/lamb
242        # scattering angle theta; always true for plane detector
243        # aligned vertically to the ko direction
244        if qr_value > knot:
245            theta = pi/2
246        else:
247            theta = math.asin(qr_value/knot)
248        # source aperture size
249        rone = self.source_aperture_size
250        # sample aperture size
251        rtwo = self.sample_aperture_size
252        # detector pixel size
253        rthree = self.detector_pix_size
254        # source to sample(aperture) distance
255        l_ssa = self.source2sample_distance[0]
256        # sample(aperture) to detector distance
257        l_sad = self.sample2detector_distance[0]
258        # sample (aperture) to sample distance
259        l_sas = self.sample2sample_distance[0]
260        # source to sample distance
261        l_one = l_ssa + l_sas
262        # sample to detector distance
263        l_two = l_sad - l_sas
264       
265        # Sample offset correction for l_one and Lp on variance calculation
266        l1_cor = (l_ssa * l_two) / (l_sas + l_two)
267        lp_cor = (l_ssa * l_two) / (l_one + l_two)
268        # the radial distance to the pixel from the center of the detector
269        radius = math.tan(theta) * l_two
270        #Lp = l_one*l_two/(l_one+l_two)
271        # default polar coordinate
272        comp1 = 'radial'
273        comp2 = 'phi'
274        # in the case of the cartesian coordinate
275        if coord == 'cartesian':
276            comp1 = 'x'
277            comp2 = 'y'
278
279        # sigma in the radial/x direction
280        # for source aperture
281        sigma_1  = self.get_variance(rone, l1_cor, phi, comp1)
282        # for sample apperture
283        sigma_1 += self.get_variance(rtwo, lp_cor, phi, comp1)
284        # for detector pix
285        sigma_1 += self.get_variance(rthree, l_two, phi, comp1)
286        # for gravity term for 1d
287        sigma_1grav1d =  self.get_variance_gravity(l_ssa, l_sad, lamb, lamb_spread, 
288                             phi, comp1, 'on') / tof_factor
289        # for wavelength spread
290        # reserve for 1d calculation
291        A_value = self._cal_A_value(lamb, l_ssa, l_sad)
292        sigma_wave_1, sigma_wave_1_1d = self.get_variance_wave(A_value, 
293                                          radius, l_two, lamb_spread, 
294                                          phi, 'radial', 'on')
295        sigma_wave_1 /= tof_factor
296        sigma_wave_1_1d /=  tof_factor
297        # for 1d
298        variance_1d_1 = (sigma_1 + sigma_1grav1d) /2 + sigma_wave_1_1d
299        # normalize
300        variance_1d_1 = knot * knot * variance_1d_1 / 12
301       
302        # for 2d
303        #sigma_1 += sigma_wave_1
304        # normalize
305        sigma_1 = knot*sqrt(sigma_1 / 12)
306        sigma_r = knot*sqrt(sigma_wave_1 / (tof_factor *12))
307        # sigma in the phi/y direction
308        # for source apperture
309        sigma_2  = self.get_variance(rone, l1_cor, phi, comp2)
310
311        # for sample apperture
312        sigma_2 += self.get_variance(rtwo, lp_cor, phi, comp2)
313
314        # for detector pix
315        sigma_2 += self.get_variance(rthree, l_two, phi, comp2)
316
317        # for gravity term for 1d
318        sigma_2grav1d =  self.get_variance_gravity(l_ssa, l_sad, lamb, lamb_spread, 
319                             phi, comp2, 'on') / tof_factor
320
321       
322        # for wavelength spread
323        # reserve for 1d calculation
324        sigma_wave_2, sigma_wave_2_1d = self.get_variance_wave(A_value, 
325                                          radius, l_two, lamb_spread, 
326                                          phi, 'phi', 'on') 
327        sigma_wave_2 /=  tof_factor
328        sigma_wave_2_1d /=  tof_factor
329        # for 1d
330        variance_1d_2 = (sigma_2 + sigma_2grav1d) / 2 + sigma_wave_2_1d
331        # normalize
332        variance_1d_2 = knot*knot*variance_1d_2 / 12
333       
334        # for 2d
335        #sigma_2 =  knot*sqrt(sigma_2/12)
336        #sigma_2 += sigma_wave_2
337        # normalize
338        sigma_2 =  knot * sqrt(sigma_2 / 12)
339        sigma1d = sqrt(variance_1d_1 + variance_1d_2)
340        # set sigmas
341        self.sigma_1 = sigma_1
342        self.sigma_lamd = sigma_r
343        self.sigma_2 = sigma_2
344        self.sigma_1d = sigma1d
345        return qr_value, phi, sigma_1, sigma_2, sigma_r, sigma1d
346   
347    def _within_detector_range(self,qx_value, qy_value):
348        """
349        check if qvalues are within detector range
350        """
351        # detector range
352        detector_qx_min = self.detector_qx_min
353        detector_qx_max = self.detector_qx_max
354        detector_qy_min = self.detector_qy_min
355        detector_qy_max = self.detector_qy_max
356        if self.qxmin_limit > detector_qx_min:
357            self.qxmin_limit = detector_qx_min
358        if self.qxmax_limit < detector_qx_max:
359            self.qxmax_limit = detector_qx_max
360        if self.qymin_limit > detector_qy_min:
361            self.qymin_limit = detector_qy_min
362        if self.qymax_limit < detector_qy_max:
363            self.qymax_limit = detector_qy_max
364        if qx_value < detector_qx_min or qx_value > detector_qx_max:
365            return False
366        if qy_value < detector_qy_min or qy_value > detector_qy_max:
367            return False
368        return True
369   
370    def get_image(self, qx_value, qy_value, sigma_1, sigma_2, sigma_r,
371                  qx_min, qx_max, qy_min, qy_max, 
372                  coord = 'cartesian', full_cal=True): 
373        """
374        Get the resolution in polar coordinate ready to plot
375        : qx_value: qx_value value
376        : qy_value: qy_value value
377        : sigma_1: variance in r direction
378        : sigma_2: variance in phi direction
379        : coord: coordinate system of image, 'polar' or 'cartesian'
380        """
381        # Get  qx_max and qy_max...
382        output = self._get_detector_qxqy_pixels()
383       
384        qr_value, phi = self._get_polar_value(qx_value, qy_value)
385
386        # Check whether the q value is within the detector range
387        msg = "Invalid input: Q value out of the detector range..."
388        if qx_min < self.qx_min:
389            self.qx_min = qx_min
390            #raise ValueError, msg
391        if qx_max > self.qx_max:
392            self.qx_max = qx_max
393            #raise ValueError, msg
394        if qy_min < self.qy_min:
395            self.qy_min = qy_min
396            #raise ValueError, msg
397        if qy_max > self.qy_max:
398            self.qy_max = qy_max
399            #raise ValueError, msg
400        if not full_cal:
401            return None
402 
403        # Make an empty graph in the detector scale
404        dx_size = (self.qx_max - self.qx_min) / (1000 - 1)
405        dy_size = (self.qy_max - self.qy_min) / (1000 - 1)
406        x_val = numpy.arange(self.qx_min, self.qx_max, dx_size)
407        y_val = numpy.arange(self.qy_max, self.qy_min, -dy_size)
408        q_1, q_2 = numpy.meshgrid(x_val, y_val)
409        #q_phi = numpy.arctan(q_1,q_2)
410        # check whether polar or cartesian
411        if coord == 'polar':
412            # Find polar values
413            qr_value, phi = self._get_polar_value(qx_value, qy_value)
414            q_1, q_2 = self._rotate_z(q_1, q_2, phi)
415            qc_1 = qr_value
416            qc_2 = 0.0
417            # Calculate the 2D Gaussian distribution image
418            image = self._gaussian2d_polar(q_1, q_2, qc_1, qc_2, 
419                                 sigma_1, sigma_2, sigma_r)
420        else:
421            # catesian coordinate
422            # qx_center
423            qc_1 = qx_value
424            # qy_center
425            qc_2 = qy_value
426           
427            # Calculate the 2D Gaussian distribution image
428            image = self._gaussian2d(q_1, q_2, qc_1, qc_2, 
429                                     sigma_1, sigma_2, sigma_r)
430        # out side of detector
431        if not self._within_detector_range(qx_value, qy_value):
432            image *= 0.0
433            self.intensity = 0.0
434            #return self.image
435
436        # Add it if there are more than one inputs.
437        if len(self.image_lam) > 0:
438            self.image_lam += image * self.intensity
439        else:
440            self.image_lam = image * self.intensity
441       
442        return self.image_lam
443   
444    def plot_image(self, image):
445        """
446        Plot image using pyplot
447        : image: 2d resolution image
448       
449        : return plt: pylab object
450        """
451        import matplotlib.pyplot as plt
452
453        self.plot = plt
454        plt.xlabel('$\\rm{Q}_{x} [A^{-1}]$')
455        plt.ylabel('$\\rm{Q}_{y} [A^{-1}]$')
456        # Max value of the image
457        max = numpy.max(image)
458        qx_min, qx_max, qy_min, qy_max = self.get_detector_qrange()
459
460        # Image
461        im = plt.imshow(image, 
462                extent = [qx_min, qx_max, qy_min, qy_max])
463
464        # bilinear interpolation to make it smoother
465        im.set_interpolation('bilinear')
466
467        return plt
468   
469    def reset_image(self):
470        """
471        Reset image to default (=[])
472        """
473        self.image = []
474       
475    def get_variance(self, size = [], distance = 0, phi = 0, comp = 'radial'):
476        """
477        Get the variance when the slit/pinhole size is given
478        : size: list that can be one(diameter for circular)
479                or two components(lengths for rectangular)
480        : distance: [z, x] where z along the incident beam, x // qx_value
481        : comp: direction of the sigma; can be 'phi', 'y', 'x', and 'radial'
482       
483        : return variance: sigma^2
484        """   
485        # check the length of size (list)
486        len_size = len(size)
487       
488        # define sigma component direction
489        if comp == 'radial':
490            phi_x = math.cos(phi)
491            phi_y = math.sin(phi)
492        elif comp == 'phi':
493            phi_x = math.sin(phi)
494            phi_y = math.cos(phi)
495        elif comp == 'x':
496            phi_x = 1
497            phi_y = 0
498        elif comp == 'y':
499            phi_x = 0
500            phi_y = 1
501        else:
502            phi_x = 0
503            phi_y = 0
504        # calculate each component
505        # for pinhole w/ radius = size[0]/2
506        if len_size == 1:
507            x_comp = (0.5 * size[0]) * sqrt(3)
508            y_comp = 0
509        # for rectangular slit
510        elif len_size == 2:
511            x_comp = size[0] * phi_x
512            y_comp = size[1] * phi_y
513        # otherwise
514        else:
515            raise ValueError, " Improper input..."
516        # get them squared 
517        sigma  = x_comp * x_comp
518        sigma += y_comp * y_comp
519        # normalize by distance
520        sigma /= (distance * distance)
521
522        return sigma
523
524    def get_variance_wave(self, A_value, radius, distance, spread, phi, 
525                          comp = 'radial', switch = 'on'):
526        """
527        Get the variance when the wavelength spread is given
528       
529        : radius: the radial distance from the beam center to the pix of q
530        : distance: sample to detector distance
531        : spread: wavelength spread (ratio)
532        : comp: direction of the sigma; can be 'phi', 'y', 'x', and 'radial'
533       
534        : return variance: sigma^2 for 2d, sigma^2 for 1d [tuple]
535        """
536        if switch.lower() == 'off':
537            return 0, 0
538        # check the singular point
539        if distance == 0 or comp == 'phi':
540            return 0, 0
541        else:
542            # calculate sigma^2 for 1d
543            sigma1d = 2 * math.pow(radius/distance*spread, 2)
544            if comp == 'x':
545                sigma1d *= (math.cos(phi)*math.cos(phi))
546            elif comp == 'y':
547                sigma1d *= (math.sin(phi)*math.sin(phi))
548            else:
549                sigma1d *= 1 
550            # sigma^2 for 2d 
551            # shift the coordinate due to the gravitational shift
552            rad_x = radius * math.cos(phi)
553            rad_y = A_value - radius * math.sin(phi)
554            radius = math.sqrt(rad_x * rad_x + rad_y * rad_y)
555            # new phi
556            phi = math.atan2(-rad_y, rad_x)
557            self.gravity_phi = phi
558            # calculate sigma^2
559            sigma = 2 * math.pow(radius/distance*spread, 2)
560            if comp == 'x':
561                sigma *= (math.cos(phi)*math.cos(phi))
562            elif comp == 'y':
563                sigma *= (math.sin(phi)*math.sin(phi))
564            else:
565                sigma *= 1         
566               
567            return sigma, sigma1d
568
569    def get_variance_gravity(self, s_distance, d_distance, wavelength, spread, 
570                             phi, comp = 'radial', switch = 'on'):
571        """
572        Get the variance from gravity when the wavelength spread is given
573       
574        : s_distance: source to sample distance
575        : d_distance: sample to detector distance
576        : wavelength: wavelength
577        : spread: wavelength spread (ratio)
578        : comp: direction of the sigma; can be 'phi', 'y', 'x', and 'radial'
579       
580        : return variance: sigma^2
581        """
582        if switch.lower() == 'off':
583            return 0
584        if self.mass == 0.0:
585            return 0
586        # check the singular point
587        if d_distance == 0 or comp == 'x':
588            return 0
589        else:
590            a_value = self._cal_A_value(None, s_distance, d_distance)
591            # calculate sigma^2
592            sigma = math.pow(a_value / d_distance, 2)
593            sigma *= math.pow(wavelength, 4)
594            sigma *= math.pow(spread, 2)
595            sigma *= 8
596           
597            # only for the polar coordinate
598            #if comp == 'radial':
599            #    sigma *= (math.sin(phi) * math.sin(phi))
600            #elif comp == 'phi':
601            #    sigma *= (math.cos(phi) * math.cos(phi))
602           
603            return sigma
604   
605    def _cal_A_value(self, lamda, s_distance, d_distance):   
606        """
607        Calculate A value for gravity
608       
609        : s_distance: source to sample distance
610        : d_distance: sample to detector distance
611        """
612        # neutron mass in cgs unit
613        self.mass = self.get_neutron_mass()
614        # plank constant in cgs unit
615        h_constant = _PLANK_H
616        # gravity in cgs unit
617        gravy = _GRAVITY
618        # m/h
619        m_over_h = self.mass /h_constant
620        # A value
621        a_value = d_distance * (s_distance + d_distance)
622        a_value *= math.pow(m_over_h / 2, 2)
623        a_value *= gravy
624        # unit correction (1/cm to 1/A) for A and d_distance below
625        a_value *= 1.0E-16
626        # if lamda is give (broad meanning of A)  return 2* lamda^2 * A
627        if lamda != None:
628            a_value *= (4 * lamda * lamda)
629        return a_value
630   
631    def get_intensity(self):
632        """
633        Get intensity
634        """
635        return self.wave.intensity
636
637    def get_wavelength(self):
638        """
639        Get wavelength
640        """
641        return self.wave.wavelength
642   
643    def get_spectrum(self):
644        """
645        Get spectrum
646        """
647        return self.wave.spectrum   
648   
649    def get_default_spectrum(self):
650        """
651        Get default_spectrum
652        """
653        return self.wave.get_default_spectrum()   
654   
655    def get_spectrum(self):
656        """
657        Get _spectrum
658        """
659        return self.wave.get_spectrum() 
660         
661    def get_wavelength_spread(self):
662        """
663        Get wavelength spread
664        """
665        return self.wave.wavelength_spread
666   
667    def get_neutron_mass(self):
668        """
669        Get Neutron mass
670        """
671        return self.wave.mass
672   
673    def get_source_aperture_size(self):
674        """
675        Get source aperture size
676        """
677        return self.aperture.source_size
678   
679    def get_sample_aperture_size(self):
680        """
681        Get sample aperture size
682        """
683        return self.aperture.sample_size
684   
685    def get_detector_pix_size(self):
686        """
687        Get detector pixel size
688        """
689        return self.detector.pix_size
690       
691    def get_detector_size(self):
692        """
693        Get detector size
694        """
695        return self.detector.size   
696     
697    def get_source2sample_distance(self):
698        """
699        Get detector source2sample_distance
700        """
701        return self.aperture.sample_distance
702   
703    def get_sample2sample_distance(self):
704        """
705        Get detector sampleslitsample_distance
706        """
707        return self.sample.distance
708   
709    def get_sample2detector_distance(self):
710        """
711        Get detector sample2detector_distance
712        """
713        return self.detector.distance
714   
715    def set_intensity(self, intensity):
716        """
717        Set intensity
718        """
719        self.wave.set_intensity(intensity)
720       
721    def set_wave(self, wavelength):
722        """
723        Set wavelength list or wavelength
724        """
725        if wavelength.__class__.__name__ == 'list':
726            self.wave.set_wave_list(wavelength)
727        elif wavelength.__class__.__name__ == 'float':
728            self.wave.set_wave_list([wavelength])
729            #self.set_wavelength(wavelength)
730        else:
731            raise
732   
733    def set_wave_spread(self, wavelength_spread):
734        """
735        Set wavelength spread  or wavelength spread
736        """
737        if wavelength_spread.__class__.__name__ == 'list':
738            self.wave.set_wave_spread_list(wavelength_spread)
739        elif wavelength_spread.__class__.__name__ == 'float':
740            self.wave.set_wave_spread_list([wavelength_spread])
741            #self.set_wavelength_spread(wavelength_spread)
742        else:
743            raise
744       
745    def set_wavelength(self, wavelength):
746        """
747        Set wavelength
748        """
749        self.wavelength = wavelength
750        self.wave.set_wavelength(wavelength)
751       
752    def set_spectrum(self, spectrum):
753        """
754        Set spectrum
755        """
756        self.spectrum = spectrum
757        self.wave.set_spectrum(spectrum) 
758         
759    def set_wavelength_spread(self, wavelength_spread):
760        """
761        Set wavelength spread
762        """
763        self.wavelength_spread = wavelength_spread
764        self.wave.set_wavelength_spread(wavelength_spread)
765       
766    def set_wave_list(self, wavelength_list, wavelengthspread_list):
767        """
768        Set wavelength and its spread list
769        """
770        self.wave.set_wave_list(wavelength_list, wavelengthspread_list)
771   
772    def get_wave_list(self):
773        """
774        Set wavelength spread
775        """
776        return self.wave.get_wave_list()
777   
778    def get_intensity_list(self):
779        """
780        Set wavelength spread
781        """
782        return self.wave.get_intensity_list()
783           
784    def set_source_aperture_size(self, size):
785        """
786        Set source aperture size
787       
788        : param size: [dia_value] or [x_value, y_value]
789        """       
790        if len(size) < 1 or len(size) > 2:
791            raise RuntimeError, "The length of the size must be one or two."
792        self.aperture.set_source_size(size)
793       
794    def set_neutron_mass(self, mass):
795        """
796        Set Neutron mass
797        """
798        self.wave.set_mass(mass)
799        self.mass = mass
800       
801    def set_sample_aperture_size(self, size):
802        """
803        Set sample aperture size
804       
805        : param size: [dia_value] or [xheight_value, yheight_value]
806        """
807        if len(size) < 1 or len(size) > 2:
808            raise RuntimeError, "The length of the size must be one or two."
809        self.aperture.set_sample_size(size)
810   
811    def set_detector_pix_size(self, size):
812        """
813        Set detector pixel size
814        """
815        self.detector.set_pix_size(size)
816       
817    def set_detector_size(self, size):
818        """
819        Set detector size in number of pixels
820        : param size: [pixel_nums] or [x_pix_num, yx_pix_num]
821        """
822        self.detector.set_size(size)
823       
824    def set_source2sample_distance(self, distance):
825        """
826        Set detector source2sample_distance
827       
828        : param distance: [distance, x_offset]
829        """
830        if len(distance) < 1 or len(distance) > 2:
831            raise RuntimeError, "The length of the size must be one or two."
832        self.aperture.set_sample_distance(distance)
833
834    def set_sample2sample_distance(self, distance):
835        """
836        Set detector sample_slit2sample_distance
837       
838        : param distance: [distance, x_offset]
839        """
840        if len(distance) < 1 or len(distance) > 2:
841            raise RuntimeError, "The length of the size must be one or two."
842        self.sample.set_distance(distance)
843   
844    def set_sample2detector_distance(self, distance):
845        """
846        Set detector sample2detector_distance
847       
848        : param distance: [distance, x_offset]
849        """
850        if len(distance) < 1 or len(distance) > 2:
851            raise RuntimeError, "The length of the size must be one or two."
852        self.detector.set_distance(distance)
853       
854    def get_all_instrument_params(self):
855        """
856        Get all instrumental parameters
857        """
858        #self.intensity = self.get_intensity()
859        #self.wavelength = self.get_wavelength()
860        #self.wavelength_spread = self.get_wavelength_spread()
861        self.mass = self.get_neutron_mass()
862        self.spectrum = self.get_spectrum()
863        self.source_aperture_size = self.get_source_aperture_size()
864        self.sample_aperture_size = self.get_sample_aperture_size()
865        self.detector_pix_size = self.get_detector_pix_size()
866        self.detector_size = self.get_detector_size()
867        self.source2sample_distance = self.get_source2sample_distance()
868        self.sample2sample_distance = self.get_sample2sample_distance()
869        self.sample2detector_distance = self.get_sample2detector_distance()
870   
871    def get_detector_qrange(self):
872        """
873        get max detector q ranges
874       
875        : return: qx_min, qx_max, qy_min, qy_max tuple
876        """
877        if len(self.qxrange) != 2 or len(self.qyrange) != 2:
878            return None
879        qx_min = self.qxrange[0]
880        qx_max = self.qxrange[1]
881        qy_min = self.qyrange[0]
882        qy_max = self.qyrange[1]
883       
884        return qx_min, qx_max, qy_min, qy_max
885   
886    def _rotate_z(self, x_value, y_value, theta= 0.0):
887        """
888        Rotate x-y cordinate around z-axis by theta
889        : x_value: numpy array of x values
890        : y_value: numpy array of y values
891        : theta: angle to rotate by in rad
892       
893        :return: x_prime, y-prime
894        """       
895        # rotate by theta
896        x_prime = x_value * math.cos(theta) + y_value * math.sin(theta)
897        y_prime =  -x_value * math.sin(theta) + y_value * math.cos(theta)
898   
899        return x_prime, y_prime
900   
901    def _gaussian2d(self, x_val, y_val, x0_val, y0_val, 
902                    sigma_x, sigma_y, sigma_r):
903        """
904        Calculate 2D Gaussian distribution
905        : x_val: x value
906        : y_val: y value
907        : x0_val: mean value in x-axis
908        : y0_val: mean value in y-axis
909        : sigma_x: variance in x-direction
910        : sigma_y: variance in y-direction
911       
912        : return: gaussian (value)
913        """
914        # phi values at each points (not at the center)
915        x_value = x_val - x0_val
916        y_value = y_val - y0_val
917        phi_i = numpy.arctan2(y_val, x_val)
918       
919        # phi correction due to the gravity shift (in phi)
920        phi_0 = math.atan2(y0_val, x0_val)
921        phi_i = phi_i - phi_0 + self.gravity_phi
922
923        sin_phi = numpy.sin(self.gravity_phi)
924        cos_phi = numpy.cos(self.gravity_phi)
925       
926        x_p = x_value * cos_phi + y_value * sin_phi
927        y_p = -x_value * sin_phi + y_value * cos_phi
928       
929        new_sig_x = sqrt(sigma_r * sigma_r / (sigma_x * sigma_x ) + 1)
930        new_sig_y = sqrt(sigma_r * sigma_r / (sigma_y * sigma_y ) + 1)
931        new_x = x_p * cos_phi / new_sig_x - y_p * sin_phi
932        new_x /= sigma_x
933        new_y = x_p * sin_phi / new_sig_y + y_p * cos_phi
934        new_y /= sigma_y
935
936        nu_value = -0.5 *(new_x * new_x + new_y * new_y)
937
938        gaussian = numpy.exp(nu_value)
939        # normalizing factor correction
940        gaussian /= gaussian.sum()
941
942        return gaussian
943
944    def _gaussian2d_polar(self, x_val, y_val, x0_val, y0_val, 
945                        sigma_x, sigma_y, sigma_r):
946        """
947        Calculate 2D Gaussian distribution for polar coodinate
948        : x_val: x value
949        : y_val: y value
950        : x0_val: mean value in x-axis
951        : y0_val: mean value in y-axis
952        : sigma_x: variance in r-direction
953        : sigma_y: variance in phi-direction
954        : sigma_r: wavelength variance in r-direction
955       
956        : return: gaussian (value)
957        """
958        sigma_x = sqrt(sigma_x * sigma_x + sigma_r * sigma_r)
959        # call gaussian1d
960        gaussian  = self._gaussian1d(x_val, x0_val, sigma_x)
961        gaussian *= self._gaussian1d(y_val, y0_val, sigma_y)
962 
963        # normalizing factor correction
964        if sigma_x != 0 and sigma_y != 0:
965            gaussian *= sqrt(2 * pi)
966        return gaussian
967   
968    def _gaussian1d(self, value, mean, sigma):
969        """
970        Calculate 1D Gaussian distribution
971        : value: value
972        : mean: mean value
973        : sigma: variance
974       
975        : return: gaussian (value)
976        """
977        # default
978        gaussian = 1.0
979        if sigma != 0:
980            # get exponent
981            nu_value = (value - mean) / sigma
982            nu_value *= nu_value
983            nu_value *= -0.5
984            gaussian *= numpy.exp(nu_value)
985            gaussian /= sigma
986            # normalize
987            gaussian /= sqrt(2 * pi)
988           
989        return gaussian
990   
991    def _atan_phi(self, qy_value, qx_value):
992        """
993        Find the angle phi of q on the detector plane for qx_value, qy_value given
994        : qx_value: x component of q
995        : qy_value: y component of q
996       
997        : return phi: the azimuthal angle of q on x-y plane
998        """
999        phi = math.atan2(qy_value, qx_value)
1000        return phi
1001        # default
1002        phi = 0
1003        # ToDo: This is misterious - sign???
1004        #qy_value = -qy_value
1005        # Take care of the singular point
1006        if qx_value == 0:
1007            if qy_value > 0:
1008                phi = pi / 2
1009            elif qy_value < 0:
1010                phi = -pi / 2
1011            else:
1012                phi = 0
1013        else:
1014            # the angle
1015            phi = math.atan2(qy_value, qx_value)
1016
1017        return phi
1018
1019    def _get_detector_qxqy_pixels(self):
1020        """
1021        Get the pixel positions of the detector in the qx_value-qy_value space
1022        """
1023       
1024        # update all param values
1025        self.get_all_instrument_params()
1026       
1027        # wavelength
1028        wavelength = self.wave.wavelength
1029        # Gavity correction
1030        delta_y = self._get_beamcenter_drop() # in cm
1031       
1032        # detector_pix size
1033        detector_pix_size = self.detector_pix_size
1034        # Square or circular pixel
1035        if len(detector_pix_size) == 1:
1036            pix_x_size = detector_pix_size[0]
1037            pix_y_size = detector_pix_size[0]
1038        # rectangular pixel pixel
1039        elif len(detector_pix_size) == 2:
1040            pix_x_size = detector_pix_size[0]
1041            pix_y_size = detector_pix_size[1]
1042        else:
1043            raise ValueError, " Input value format error..."
1044        # Sample to detector distance = sample slit to detector
1045        # minus sample offset
1046        sample2detector_distance = self.sample2detector_distance[0] - \
1047                                    self.sample2sample_distance[0]
1048        # detector offset in x-direction
1049        detector_offset = 0
1050        try:
1051            detector_offset = self.sample2detector_distance[1]
1052        except:
1053            pass
1054       
1055        # detector size in [no of pix_x,no of pix_y]
1056        detector_pix_nums_x = self.detector_size[0]
1057       
1058        # get pix_y if it exists, otherwse take it from [0]
1059        try:
1060            detector_pix_nums_y = self.detector_size[1]
1061        except:
1062            detector_pix_nums_y = self.detector_size[0]
1063       
1064        # detector offset in pix number
1065        offset_x = detector_offset / pix_x_size
1066        offset_y = delta_y / pix_y_size
1067       
1068        # beam center position in pix number (start from 0)
1069        center_x, center_y = self._get_beamcenter_position(detector_pix_nums_x,
1070                                    detector_pix_nums_y, offset_x, offset_y)
1071        # distance [cm] from the beam center on detector plane
1072        detector_ind_x = numpy.arange(detector_pix_nums_x)
1073        detector_ind_y = numpy.arange(detector_pix_nums_y)
1074
1075        # shif 0.5 pixel so that pix position is at the center of the pixel
1076        detector_ind_x = detector_ind_x + 0.5
1077        detector_ind_y = detector_ind_y + 0.5
1078
1079        # the relative postion from the beam center
1080        detector_ind_x = detector_ind_x - center_x
1081        detector_ind_y = detector_ind_y - center_y
1082       
1083        # unit correction in cm
1084        detector_ind_x = detector_ind_x * pix_x_size
1085        detector_ind_y = detector_ind_y * pix_y_size
1086       
1087        qx_value = numpy.zeros(len(detector_ind_x))
1088        qy_value = numpy.zeros(len(detector_ind_y))
1089        i = 0
1090
1091        for indx in detector_ind_x:
1092            qx_value[i] = self._get_qx(indx, sample2detector_distance, wavelength)
1093            i += 1
1094        i = 0
1095        for indy in detector_ind_y:
1096            qy_value[i] = self._get_qx(indy, sample2detector_distance, wavelength)
1097            i += 1
1098           
1099        # qx_value and qy_value values in array
1100        qx_value = qx_value.repeat(detector_pix_nums_y) 
1101        qx_value = qx_value.reshape(detector_pix_nums_x, detector_pix_nums_y)
1102        qy_value = qy_value.repeat(detector_pix_nums_x) 
1103        qy_value = qy_value.reshape(detector_pix_nums_y, detector_pix_nums_x)
1104        qy_value = qy_value.transpose()
1105
1106        # p min and max values among the center of pixels
1107        self.qx_min = numpy.min(qx_value)
1108        self.qx_max = numpy.max(qx_value)
1109        self.qy_min = numpy.min(qy_value)
1110        self.qy_max = numpy.max(qy_value)
1111               
1112        # Appr. min and max values of the detector display limits
1113        # i.e., edges of the last pixels.
1114        self.qy_min += self._get_qx(-0.5 * pix_y_size, 
1115                                sample2detector_distance, wavelength)
1116        self.qy_max += self._get_qx(0.5 * pix_y_size, 
1117                                sample2detector_distance, wavelength)
1118        #if self.qx_min == self.qx_max:
1119        self.qx_min += self._get_qx(-0.5 * pix_x_size, 
1120                                sample2detector_distance, wavelength)
1121        self.qx_max += self._get_qx(0.5 * pix_x_size, 
1122                                    sample2detector_distance, wavelength)
1123       
1124        # min and max values of detecter
1125        self.detector_qx_min = self.qx_min
1126        self.detector_qx_max = self.qx_max
1127        self.detector_qy_min = self.qy_min
1128        self.detector_qy_max = self.qy_max
1129       
1130       
1131
1132        # try to set it as a Data2D otherwise pass (not required for now)
1133        try:
1134            from sans.dataloader.data_info import Data2D
1135            output = Data2D()
1136            inten = numpy.zeros_like(qx_value)
1137            output.data     = inten
1138            output.qx_data  = qx_value
1139            output.qy_data  = qy_value           
1140        except:
1141            pass
1142       
1143        return output#qx_value,qy_value
1144       
1145    def _get_qx(self, dx_size, det_dist, wavelength):
1146        """
1147        :param dx_size: x-distance from beam center [cm]
1148        :param det_dist: sample to detector distance [cm]
1149       
1150        :return: q-value at the given position
1151        """
1152        # Distance from beam center in the plane of detector
1153        plane_dist = dx_size
1154        # full scattering angle on the x-axis
1155        theta  = numpy.arctan(plane_dist / det_dist)
1156        qx_value     = (2.0 * pi / wavelength) * numpy.sin(theta)
1157        return qx_value 
1158   
1159    def _get_polar_value(self, qx_value, qy_value):
1160        """
1161        Find qr_value and phi from qx_value and qy_value values
1162       
1163        : return qr_value, phi
1164        """
1165        # find |q| on detector plane
1166        qr_value = sqrt(qx_value*qx_value + qy_value*qy_value)
1167        # find angle phi
1168        phi = self._atan_phi(qy_value, qx_value)
1169       
1170        return qr_value, phi
1171   
1172    def _get_beamcenter_position(self, num_x, num_y, offset_x, offset_y):
1173        """
1174        :param num_x: number of pixel in x-direction
1175        :param num_y: number of pixel in y-direction
1176        :param offset: detector offset in x-direction in pix number
1177       
1178        :return: pix number; pos_x, pos_y in pix index
1179        """
1180        # beam center position
1181        pos_x = num_x / 2
1182        pos_y = num_y / 2
1183
1184        # correction for offset
1185        pos_x += offset_x
1186        # correction for gravity that is always negative
1187        pos_y -= offset_y
1188
1189        return pos_x, pos_y     
1190
1191    def _get_beamcenter_drop(self):
1192        """
1193        Get the beam center drop (delta y) in y diection due to gravity
1194       
1195        :return delta y: the beam center drop in cm
1196        """
1197        # Check if mass == 0 (X-ray).
1198        if self.mass == 0:
1199            return 0
1200        # Covert unit from A to cm
1201        unit_cm = 1e-08
1202        # Velocity of neutron in horizontal direction (~ actual velocity)
1203        velocity = _PLANK_H / (self.mass * self.wave.wavelength * unit_cm)
1204        # Compute delta y
1205        delta_y = 0.5
1206        delta_y *= _GRAVITY
1207        sampletodetector = self.sample2detector_distance[0] - \
1208                                    self.sample2sample_distance[0]
1209        delta_y *= sampletodetector
1210        delta_y *= (self.source2sample_distance[0] + self.sample2detector_distance[0])
1211        delta_y /= (velocity * velocity)
1212
1213        return delta_y     
1214
1215   
Note: See TracBrowser for help on using the repository browser.