source: sasview/Invariant/invariant.py @ aafa962

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

invariant: minor code refactor

  • Property mode set to 100644
File size: 37.5 KB
Line 
1"""
2    This module implements invariant and its related computations.
3    @author: Gervaise B. Alina/UTK
4"""
5#TODO: Need to decide if/how to use smearing for extrapolation
6import math 
7import numpy
8
9from DataLoader.data_info import Data1D as LoaderData1D
10from DataLoader.qsmearing import smear_selection
11
12# The minimum q-value to be used when extrapolating
13Q_MINIMUM  = 1e-5
14
15# The maximum q-value to be used when extrapolating
16Q_MAXIMUM  = 10
17
18# Number of steps in the extrapolation
19INTEGRATION_NSTEPS = 1000
20
21class Transform(object):
22    """
23        Define interface that need to compute a function or an inverse
24        function given some x, y
25    """
26    def linearize_data(self, x, y, dy=None):
27        """
28            Linearize data so that a linear fit can be performed.
29            Filter out the data that can't be transformed.
30            @param x: array of q values
31            @param y: array of I(q) values
32            @param dy: array of dI(q) values
33        """
34        return NotImplemented
35   
36    def linearize_q_value(self, value):
37        """
38            Transform the input q-value for linearization
39        """
40        return NotImplemented
41
42    def extract_model_parameters(self, a, b):
43        """
44            set private member
45        """
46        return NotImplemented
47     
48    def evaluate_model(self, x):
49        """
50            Returns an array f(x) values where f is the Transform function.
51        """
52        return NotImplemented
53   
54class Guinier(Transform):
55    """
56        class of type Transform that performs operations related to guinier
57        function
58    """
59    def __init__(self, scale=1, radius=60):
60        Transform.__init__(self)
61        self.scale = scale
62        self.radius = radius
63   
64    def linearize_data(self, x, y, dy=None):
65        """
66            Linearize data so that a linear fit can be performed.
67            Filter out the data that can't be transformed.
68            @param x: array of q values
69            @param y: array of I(q) values
70            @param dy: array of dI(q) values
71        """
72        # Check that the vector lengths are equal
73        assert(len(x)==len(y))
74        if dy is not None:
75            assert(len(x)==len(dy))
76        else:
77            dy = numpy.array([math.sqrt(math.fabs(j)) for j in y])
78       
79        data_points = zip(x,y,dy)
80       
81        # Transform the data
82        output_points = [(self.linearize_q_value(p[0]),
83                          math.log(p[1]),
84                          p[2]/p[1]) for p in data_points if p[0]>0 and p[1]>0]
85       
86        x_out, y_out, dy_out = zip(*output_points)
87
88        return numpy.asarray(x_out), numpy.asarray(y_out), numpy.asarray(dy_out)
89   
90    def linearize_q_value(self, value):
91        """
92            Transform the input q-value for linearization
93            @param value: q-value
94            @return: q*q
95        """
96        return value * value
97   
98    def extract_model_parameters(self, a, b):
99        """
100           assign new value to the scale and the radius
101        """
102        b = math.sqrt(-3 * b)
103        a = math.exp(a)
104        self.scale = a
105        self.radius = b
106        return a, b
107       
108    def evaluate_model(self, x):
109        """
110            return F(x)= scale* e-((radius*x)**2/3)
111        """
112        return self._guinier(x)
113             
114    def _guinier(self, x):
115        """
116            Retrive the guinier function after apply an inverse guinier function
117            to x
118            Compute a F(x) = scale* e-((radius*x)**2/3).
119            @param x: a vector of q values
120            @param scale: the scale value
121            @param radius: the guinier radius value
122            @return F(x)
123        """   
124        # transform the radius of coming from the inverse guinier function to a
125        # a radius of a guinier function
126        if self.radius <= 0:
127            raise ValueError("Rg expected positive value, but got %s"%self.radius) 
128        value = numpy.array([math.exp(-((self.radius * i)**2/3)) for i in x ]) 
129        return self.scale * value
130
131class PowerLaw(Transform):
132    """
133        class of type transform that perform operation related to power_law
134        function
135    """
136    def __init__(self, scale=1, power=4):
137        Transform.__init__(self)
138        self.scale = scale
139        self.power = power
140       
141    def linearize_data(self, x, y, dy=None):
142        """
143            Linearize data so that a linear fit can be performed.
144            Filter out the data that can't be transformed.
145            @param x: array of q values
146            @param y: array of I(q) values
147            @param dy: array of dI(q) values
148        """
149        # Check that the vector lengths are equal
150        assert(len(x)==len(y))
151        if dy is not None:
152            assert(len(x)==len(dy))
153        else:
154            dy = numpy.array([math.sqrt(math.fabs(j)) for j in y])
155       
156        data_points = zip(x,y,dy)
157       
158        # Transform the data
159        output_points = [(self.linearize_q_value(p[0]),
160                          math.log(p[1]),
161                          p[2]/p[1]) for p in data_points if p[0]>0 and p[1]>0]
162       
163        x_out, y_out, dy_out = zip(*output_points)
164
165        return numpy.asarray(x_out), numpy.asarray(y_out), numpy.asarray(dy_out)
166       
167    def linearize_q_value(self, value):
168        """
169            Transform the input q-value for linearization
170            @param value: q-value
171            @return: log(q)
172        """
173        return math.log(value)
174   
175    def extract_model_parameters(self, a, b):
176        """
177            Assign new value to the scale and the power
178        """
179        b = -1 * b
180        a = math.exp(a)
181        self.power = b
182        self.scale = a
183        return a, b
184       
185    def evaluate_model(self, x):
186        """
187            given a scale and a radius transform x, y using a power_law
188            function
189        """
190        return self._power_law(x)
191       
192    def _power_law(self, x):
193        """
194            F(x) = scale* (x)^(-power)
195                when power= 4. the model is porod
196                else power_law
197            The model has three parameters:
198            @param x: a vector of q values
199            @param power: power of the function
200            @param scale : scale factor value
201            @param F(x)
202        """
203        if self.power <= 0:
204            raise ValueError("Power_law function expected positive power, but got %s"%self.power)
205        if self.scale <= 0:
206            raise ValueError("scale expected positive value, but got %s"%self.scale) 
207       
208        value = numpy.array([ math.pow(i, -self.power) for i in x ]) 
209        return self.scale * value
210
211class Extrapolator:
212    """
213        Extrapolate I(q) distribution using a given model
214    """
215    def __init__(self, data):
216        """
217            Determine a and b given a linear equation y = ax + b
218            @param Data: data containing x and y  such as  y = ax + b
219        """
220        self.data  = data
221       
222        # Set qmin as the lowest non-zero value
223        self.qmin = Q_MINIMUM
224        for q_value in self.data.x:
225            if q_value > 0: 
226                self.qmin = q_value
227                break
228        self.qmax = max(self.data.x)
229       
230        #get the smear object of data
231        self.smearer = smear_selection(self.data)
232        # Set the q-range information to allow smearing
233        self.set_fit_range()
234     
235    def set_fit_range(self, qmin=None, qmax=None):
236        """ to set the fit range"""
237        if qmin is not None:
238            self.qmin = qmin
239        if qmax is not None:
240            self.qmax = qmax
241           
242        # Determine the range needed in unsmeared-Q to cover
243        # the smeared Q range
244        self._qmin_unsmeared = self.qmin
245        self._qmax_unsmeared = self.qmax   
246       
247        if self.smearer is not None:
248            self._first_unsmeared_bin, self._last_unsmeared_bin = self.smearer.get_bin_range(self.qmin, self.qmax)
249            self._qmin_unsmeared = self.data.x[self._first_unsmeared_bin]
250            self._qmax_unsmeared = self.data.x[self._last_unsmeared_bin]
251           
252        # Identify the bin range for the unsmeared and smeared spaces
253        self.idx = (self.data.x >= self.qmin) & (self.data.x <= self.qmax)
254        self.idx_unsmeared = (self.data.x >= self._qmin_unsmeared) \
255                                & (self.data.x <= self._qmax_unsmeared)
256 
257    def fit(self, power=None):
258        """
259           Fit data for y = ax + b  return a and b
260           @param power = a fixed, otherwise None
261        """
262        fx = numpy.zeros(len(self.data.x))
263
264         # Uncertainty
265        if type(self.data.dy)==numpy.ndarray and len(self.data.dy)==len(self.data.x):
266            sigma = self.data.dy
267        else:
268            sigma = numpy.ones(len(self.data.x))
269           
270        # Compute theory data f(x)
271        fx[self.idx_unsmeared] = self.data.y[self.idx_unsmeared]
272        ## Smear theory data
273        if self.smearer is not None:
274            fx = self.smearer(fx, self._first_unsmeared_bin,self._last_unsmeared_bin)
275       
276        fx[self.idx_unsmeared] = fx[self.idx_unsmeared]/sigma[self.idx_unsmeared]
277       
278        ##power is given only for function = power_law   
279        if power != None:
280            sigma2 = sigma * sigma
281            a = -(power)
282            b = (numpy.sum(fx[self.idx]/sigma[self.idx]) - a*numpy.sum(self.data.x[self.idx]/sigma2[self.idx]))/numpy.sum(numpy.ones(len(sigma2[self.idx]))/sigma2[self.idx])
283            return a, b
284        else:
285            A = numpy.vstack([ self.data.x[self.idx]/sigma[self.idx],
286                               numpy.ones(len(self.data.x[self.idx]))/sigma[self.idx]]).T
287           
288            a, b = numpy.linalg.lstsq(A, fx[self.idx])[0]
289           
290            return a, b
291       
292
293class InvariantCalculator(object):
294    """
295        Compute invariant if data is given.
296        Can provide volume fraction and surface area if the user provides
297        Porod constant  and contrast values.
298        @precondition:  the user must send a data of type DataLoader.Data1D
299                        the user provide background and scale values.
300                       
301        @note: Some computations depends on each others.
302    """
303   
304   
305    def __init__(self, data, background=0, scale=1 ):
306        """
307            Initialize variables
308            @param data: data must be of type DataLoader.Data1D
309            @param background: Background value. The data will be corrected before processing
310            @param scale: Scaling factor for I(q). The data will be corrected before processing
311        """
312        # Background and scale should be private data member if the only way to
313        # change them are by instantiating a new object.
314        self._background = background
315        self._scale = scale
316       
317        # The data should be private
318        self._data = self._get_data(data)
319     
320        # Since there are multiple variants of Q*, you should force the
321        # user to use the get method and keep Q* a private data member
322        self._qstar = None
323       
324        # You should keep the error on Q* so you can reuse it without
325        # recomputing the whole thing.
326        self._qstar_err = 0
327       
328        # Extrapolation parameters
329        self._low_extrapolation_npts = 4
330        self._low_extrapolation_function = Guinier()
331        self._low_extrapolation_power = None
332   
333        self._high_extrapolation_npts = 4
334        self._high_extrapolation_function = PowerLaw()
335        self._high_extrapolation_power = None
336       
337    def _get_data(self, data):
338        """
339            @note this function must be call before computing any type
340             of invariant
341            @return data= self._scale *data - self._background
342        """
343        if not issubclass(data.__class__, LoaderData1D):
344            #Process only data that inherited from DataLoader.Data_info.Data1D
345            raise ValueError,"Data must be of type DataLoader.Data1D"
346        new_data = (self._scale * data) - self._background   
347        new_data.dxl = data.dxl
348        new_data.dxw = data.dxw
349        return  new_data
350     
351    def _fit(self, model, qmin=Q_MINIMUM, qmax=Q_MAXIMUM, power=None):
352        """
353            fit data with function using
354            data= self._get_data()
355            fx= Functor(data , function)
356            y = data.y
357            slope, constant = linalg.lstsq(y,fx)
358            @param qmin: data first q value to consider during the fit
359            @param qmax: data last q value to consider during the fit
360            @param power : power value to consider for power-law
361            @param function: the function to use during the fit
362            @return a: the scale of the function
363            @return b: the other parameter of the function for guinier will be radius
364                    for power_law will be the power value
365        """
366        # Linearize the data in preparation for fitting
367        fit_x, fit_y, fit_dy = model.linearize_data(self._data.x, self._data.y, self._data.dy)
368     
369        qmin = model.linearize_q_value(qmin)
370        qmax = model.linearize_q_value(qmax)
371       
372        # Create a new Data1D object for processing
373        fit_data = LoaderData1D(x=fit_x, y=fit_y, dy=fit_dy)
374        fit_data.dxl = self._data.dxl
375        fit_data.dxw = self._data.dxw   
376        extrapolator = Extrapolator(data=fit_data)
377        extrapolator.set_fit_range(qmin=qmin, qmax=qmax)
378        b, a = extrapolator.fit(power=power) 
379       
380        return model.extract_model_parameters(a=a, b=b)
381   
382    def _get_qstar(self, data):
383        """
384            Compute invariant for data
385            @param data: data to use to compute invariant.
386             
387        """
388        if data is None:
389            return 0
390        if data.is_slit_smeared():
391            return self._get_qstar_smear(data)
392        else:
393            return self._get_qstar_unsmear(data)
394       
395    def _get_qstar_unsmear(self, data):
396        """
397            Compute invariant for pinhole data.
398            This invariant is given by:
399       
400                q_star = x0**2 *y0 *dx0 +x1**2 *y1 *dx1
401                            + ..+ xn**2 *yn *dxn
402                           
403            where n >= len(data.x)-1
404            dxi = 1/2*(xi+1 - xi) + (xi - xi-1)
405            dx0 = (x1 - x0)/2
406            dxn = (xn - xn-1)/2
407            @param data: the data to use to compute invariant.
408            @return q_star: invariant value for pinhole data. q_star > 0
409        """
410        if len(data.x) <= 1 or len(data.y) <= 1 or len(data.x)!= len(data.y):
411            msg =  "Length x and y must be equal"
412            msg += " and greater than 1; got x=%s, y=%s"%(len(data.x), len(data.y))
413            raise ValueError, msg
414        else:
415            n = len(data.x)- 1
416            #compute the first delta q
417            dx0 = (data.x[1] - data.x[0])/2
418            #compute the last delta q
419            dxn = (data.x[n] - data.x[n-1])/2
420            sum = 0
421            sum += data.x[0] * data.x[0] * data.y[0] * dx0
422            sum += data.x[n] * data.x[n] * data.y[n] * dxn
423           
424            if len(data.x) == 2:
425                return sum
426            else:
427                #iterate between for element different from the first and the last
428                for i in xrange(1, n-1):
429                    dxi = (data.x[i+1] - data.x[i-1])/2
430                    sum += data.x[i] * data.x[i] * data.y[i] * dxi
431                return sum
432           
433    def _get_qstar_smear(self, data):
434        """
435            Compute invariant for slit-smeared data.
436            This invariant is given by:
437                q_star = x0*dxl *y0*dx0 + x1*dxl *y1 *dx1
438                            + ..+ xn*dxl *yn *dxn
439            where n >= len(data.x)-1
440            dxi = 1/2*(xi+1 - xi) + (xi - xi-1)
441            dx0 = (x1 - x0)/2
442            dxn = (xn - xn-1)/2
443            dxl: slit smear value
444           
445            @return q_star: invariant value for slit smeared data.
446        """
447        if not data.is_slit_smeared():
448            msg = "_get_qstar_smear need slit smear data "
449            msg += "Hint :dxl= %s , dxw= %s"%(str(data.dxl), str(data.dxw))
450            raise ValueError, msg
451
452        if len(data.x) <= 1 or len(data.y) <= 1 or len(data.x) != len(data.y)\
453              or len(data.x)!= len(data.dxl):
454            msg = "x, dxl, and y must be have the same length and greater than 1"
455            msg +=" len x=%s, y=%s, dxl=%s"%(len(data.x),len(data.y),len(data.dxl))
456            raise ValueError, msg
457        else:
458            n = len(data.x)-1
459            #compute the first delta
460            dx0 = (data.x[1] - data.x[0])/2
461            #compute the last delta
462            dxn = (data.x[n] - data.x[n-1])/2
463            sum = 0
464            sum += data.x[0] * data.dxl[0] * data.y[0] * dx0
465            sum += data.x[n] * data.dxl[n] * data.y[n] * dxn
466           
467            if len(data.x)==2:
468                return sum
469            else:
470                #iterate between for element different from the first and the last
471                for i in xrange(1, n-1):
472                    dxi = (data.x[i+1] - data.x[i-1])/2
473                    sum += data.x[i] * data.dxl[i] * data.y[i] * dxi
474                return sum
475       
476    def _get_qstar_uncertainty(self, data=None):
477        """
478            Compute uncertainty of invariant value
479            Implementation:
480                if data is None:
481                    data = self.data
482           
483                if data.slit smear:
484                        return self._get_qstar_smear_uncertainty(data)
485                    else:
486                        return self._get_qstar_unsmear_uncertainty(data)
487         
488            @param: data use to compute the invariant which allow uncertainty
489            computation.
490            @return: uncertainty
491        """
492        if data is None:
493            data = self._data
494        if data.is_slit_smeared():
495            return self._get_qstar_smear_uncertainty(data)
496        else:
497            return self._get_qstar_unsmear_uncertainty(data)
498       
499    def _get_qstar_unsmear_uncertainty(self, data):
500        """
501            Compute invariant uncertainty with with pinhole data.
502            This uncertainty is given as follow:
503               dq_star = math.sqrt[(x0**2*(dy0)*dx0)**2 +
504                    (x1**2 *(dy1)*dx1)**2 + ..+ (xn**2 *(dyn)*dxn)**2 ]
505            where n >= len(data.x)-1
506            dxi = 1/2*(xi+1 - xi) + (xi - xi-1)
507            dx0 = (x1 - x0)/2
508            dxn = (xn - xn-1)/2
509            dyn: error on dy
510           
511            @param data:
512            note: if data doesn't contain dy assume dy= math.sqrt(data.y)
513        """
514        if len(data.x) <= 1 or len(data.y) <= 1 or \
515            len(self.data.x) != len(self.data.y):
516            msg = "Length of data.x and data.y must be equal"
517            msg += " and greater than 1; got x=%s, y=%s"%(len(data.x),
518                                                         len(data.y))
519            raise ValueError, msg
520        else:
521            #Create error for data without dy error
522            if (data.dy is None) or (not data.dy):
523                dy = math.sqrt(y) 
524            else:
525                dy = data.dy
526               
527            n = len(data.x) - 1
528            #compute the first delta
529            dx0 = (data.x[1] - data.x[0])/2
530            #compute the last delta
531            dxn= (data.x[n] - data.x[n-1])/2
532            sum = 0
533            sum += (data.x[0] * data.x[0] * dy[0] * dx0)**2
534            sum += (data.x[n] * data.x[n] * dy[n] * dxn)**2
535            if len(data.x) == 2:
536                return math.sqrt(sum)
537            else:
538                #iterate between for element different from the first and the last
539                for i in xrange(1, n-1):
540                    dxi = (data.x[i+1] - data.x[i-1])/2
541                    sum += (data.x[i] * data.x[i] * dy[i] * dxi)**2
542                return math.sqrt(sum)
543       
544    def _get_qstar_smear_uncertainty(self, data):
545        """
546            Compute invariant uncertainty with slit smeared data.
547            This uncertainty is given as follow:
548                dq_star = x0*dxl *dy0 *dx0 + x1*dxl *dy1 *dx1
549                            + ..+ xn*dxl *dyn *dxn
550            where n >= len(data.x)-1
551            dxi = 1/2*(xi+1 - xi) + (xi - xi-1)
552            dx0 = (x1 - x0)/2
553            dxn = (xn - xn-1)/2
554            dxl: slit smearing value
555            dyn : error on dy
556            @param data: data of type Data1D where the scale is applied
557                        and the background is subtracted.
558         
559            note: if data doesn't contain dy assume dy= math.sqrt(data.y)
560        """
561        if not data.is_slit_smeared():
562            msg = "_get_qstar_smear_uncertainty need slit smear data "
563            msg += "Hint :dxl= %s , dxw= %s"%(str(data.dxl), str(data.dxw))
564            raise ValueError, msg
565
566        if len(data.x) <= 1 or len(data.y) <= 1 or len(data.x) != len(data.y)\
567                or len(data.x) != len(data.dxl):
568            msg = "x, dxl, and y must be have the same length and greater than 1"
569            raise ValueError, msg
570        else:
571            #Create error for data without dy error
572            if (data.dy is None) or (not data.dy):
573                dy = math.sqrt(y)
574            else:
575                dy = data.dy
576               
577            n = len(data.x) - 1
578            #compute the first delta
579            dx0 = (data.x[1] - data.x[0])/2
580            #compute the last delta
581            dxn = (data.x[n] - data.x[n-1])/2
582            sum = 0
583            sum += (data.x[0] * data.dxl[0] * dy[0] * dx0)**2
584            sum += (data.x[n] * data.dxl[n] * dy[n] * dxn)**2
585           
586            if len(data.x) == 2:
587                return math.sqrt(sum)
588            else:
589                #iterate between for element different from the first and the last
590                for i in xrange(1, n-1):
591                    dxi = (data.x[i+1] - data.x[i-1])/2
592                    sum += (data.x[i] * data.dxl[i] * dy[i] * dxi)**2
593                return math.sqrt(sum)
594   
595    def get_qstar_low(self):
596        """
597            Compute the invariant for extrapolated data at low q range.
598           
599            Implementation:
600                data = self.get_extra_data_low()
601                return self._get_qstar()
602               
603            @return q_star: the invariant for data extrapolated at low q.
604        """
605        data = self.get_extra_data_low()
606        return self._get_qstar(data=data)
607       
608    def get_qstar_high(self):
609        """
610            Compute the invariant for extrapolated data at high q range.
611           
612            Implementation:
613                data = self.get_extra_data_high()
614                return self._get_qstar()
615               
616            @return q_star: the invariant for data extrapolated at high q.
617        """
618        data = self.get_extra_data_high()
619        return self._get_qstar(data=data)
620       
621    def _get_extrapolated_data(self, model, npts=INTEGRATION_NSTEPS,
622                              q_start=Q_MINIMUM, q_end=Q_MAXIMUM):
623        """
624            @return extrapolate data create from data
625        """
626        #create new Data1D to compute the invariant
627        q = numpy.linspace(start=q_start,
628                               stop=q_end,
629                               num=npts,
630                               endpoint=True)
631        iq = model.evaluate_model(q)
632         
633        # Determine whether we are extrapolating to high or low q-values
634        # If the data is slit smeared, get the index of the slit dimension array entry
635        # that we will use to smear the extrapolated data.
636        dxl = None
637        dxw = None
638
639        if self._data.is_slit_smeared():
640            if q_start<min(self._data.x):
641                smear_index = 0
642            elif q_end>max(self._data.x):
643                smear_index = len(self._data.x)-1
644            else:
645                raise RuntimeError, "Extrapolation can only be evaluated for points outside the data Q range"
646           
647            if self._data.dxl is not None :
648                dxl = numpy.ones(INTEGRATION_NSTEPS)
649                dxl = dxl * self._data.dxl[smear_index]
650               
651            if self._data.dxw is not None :
652                dxw = numpy.ones(INTEGRATION_NSTEPS)
653                dxw = dxw * self._data.dxw[smear_index]
654             
655        result_data = LoaderData1D(x=q, y=iq)
656        result_data.dxl = dxl
657        result_data.dxw = dxw
658        return result_data
659       
660    def get_extra_data_low(self):
661        """
662            This method creates a new data set from the invariant calculator.
663           
664            It will use the extrapolation parameters kept as private data members.
665           
666            self._low_extrapolation_npts is the number of data points to use in to fit.
667            self._low_extrapolation_function will be used as the fit function.
668           
669           
670           
671            It takes npts first points of data, fits them with a given model
672            then uses the new parameters resulting from the fit to create a new data set.
673           
674            The new data first point is Q_MINIMUM.
675           
676            The last point of the new data is the first point of the original data.
677            the number of q points of this data is INTEGRATION_NSTEPS.
678           
679            @return: a new data of type Data1D
680        """
681       
682        # Data boundaries for fitting
683        qmin = self._data.x[0]
684        qmax = self._data.x[self._low_extrapolation_npts - 1]
685       
686        # Extrapolate the low-Q data
687        a, b = self._fit(model=self._low_extrapolation_function,
688                         qmin=qmin,
689                         qmax=qmax,
690                         power=self._low_extrapolation_power)
691        #q_start point
692        q_start = Q_MINIMUM
693        if Q_MINIMUM >= qmin:
694            q_start = qmin/10
695       
696        data_min = self._get_extrapolated_data(model=self._low_extrapolation_function,
697                                               npts=INTEGRATION_NSTEPS,
698                                               q_start=q_start, q_end=qmin)
699        return data_min
700         
701    def get_extra_data_high(self):
702        """
703            This method creates a new data from the invariant calculator.
704           
705            It takes npts last points of data, fits them with a given model
706            (for this function only power_law will be use), then uses
707            the new parameters resulting from the fit to create a new data set.
708            The first point is the last point of data.
709            The last point of the new data is Q_MAXIMUM.
710            The number of q points of this data is INTEGRATION_NSTEPS.
711
712           
713            @return: a new data of type Data1D
714        """
715        # Data boundaries for fitting
716        x_len = len(self._data.x) - 1
717        qmin = self._data.x[x_len - (self._high_extrapolation_npts - 1)]
718        qmax = self._data.x[x_len]
719        q_end = Q_MAXIMUM
720       
721        # fit the data with a model to get the appropriate parameters
722        a, b = self._fit(model=self._high_extrapolation_function,
723                         qmin=qmin,
724                         qmax=qmax,
725                         power=self._high_extrapolation_power)
726       
727        #create new Data1D to compute the invariant
728        data_max = self._get_extrapolated_data(model=self._high_extrapolation_function,
729                                               npts=INTEGRATION_NSTEPS,
730                                               q_start=qmax, q_end=q_end)
731        return data_max
732     
733    def set_extrapolation(self, range, npts=4, function=None, power=None):
734        """
735            Set the extrapolation parameters for the high or low Q-range.
736            Note that this does not turn extrapolation on or off.
737            @param range: a keyword set the type of extrapolation . type string
738            @param npts: the numbers of q points of data to consider for extrapolation
739            @param function: a keyword to select the function to use for extrapolation.
740                of type string.
741            @param power: an power to apply power_low function
742               
743        """
744        range = range.lower()
745        if range not in ['high', 'low']:
746            raise ValueError, "Extrapolation range should be 'high' or 'low'"
747        function = function.lower()
748        if function not in ['power_law', 'guinier']:
749            raise ValueError, "Extrapolation function should be 'guinier' or 'power_law'"
750       
751        if range == 'high':
752            if function != 'power_law':
753                raise ValueError, "Extrapolation only allows a power law at high Q"
754            self._high_extrapolation_npts  = npts
755            self._high_extrapolation_power = power
756        else:
757            if function == 'power_law':
758                self._low_extrapolation_function = PowerLaw()
759            else:
760                self._low_extrapolation_function = Guinier()
761            self._low_extrapolation_npts  = npts
762            self._low_extrapolation_power = power
763       
764    def get_qstar(self, extrapolation=None):
765        """
766            Compute the invariant of the local copy of data.
767            Implementation:
768                if slit smear:
769                    qstar_0 = self._get_qstar_smear()
770                else:
771                    qstar_0 = self._get_qstar_unsmear()
772                if extrapolation is None:
773                    return qstar_0   
774                if extrapolation==low:
775                    return qstar_0 + self.get_qstar_low()
776                elif extrapolation==high:
777                    return qstar_0 + self.get_qstar_high()
778                elif extrapolation==both:
779                    return qstar_0 + self.get_qstar_low() + self.get_qstar_high()
780           
781            @param extrapolation: string to apply optional extrapolation   
782            @return q_star: invariant of the data within data's q range
783           
784            @warning: When using setting data to Data1D , the user is responsible of
785            checking that the scale and the background are properly apply to the data
786           
787            @warning: if error occur self.get_qstar_low() or self.get_qstar_low()
788            their returned value will be ignored
789        """
790        qstar_0 = self._get_qstar(data=self._data)
791       
792        if extrapolation is None:
793            self._qstar = qstar_0
794            return self._qstar
795        # Compute invariant plus invaraint of extrapolated data
796        extrapolation = extrapolation.lower()   
797        if extrapolation == "low":
798            self._qstar = qstar_0 + self.get_qstar_low()
799            return self._qstar
800        elif extrapolation == "high":
801            self._qstar = qstar_0 + self.get_qstar_high()
802            return self._qstar
803        elif extrapolation == "both":
804            self._qstar = qstar_0 + self.get_qstar_low() + self.get_qstar_high()
805            return self._qstar
806       
807    def get_surface(self,contrast, porod_const):
808        """
809            Compute the surface of the data.
810           
811            Implementation:
812                V= self.get_volume_fraction(contrast)
813       
814              Compute the surface given by:
815                surface = (2*pi *V(1- V)*porod_const)/ q_star
816               
817            @param contrast: contrast value to compute the volume
818            @param porod_const: Porod constant to compute the surface
819            @return: specific surface
820        """
821        if contrast is None or porod_const is None:
822            return 
823        #Check whether we have Q star
824        if self._qstar is None:
825            self._qstar = self.get_star()
826        if self._qstar == 0:
827            raise RuntimeError("Cannot compute surface, invariant value is zero")
828        # Compute the volume
829        volume = self.get_volume_fraction(contrast)
830        return 2 * math.pi * volume *(1 - volume) * float(porod_const)/self._qstar
831       
832    def get_volume_fraction(self, contrast):
833        """
834            Compute volume fraction is deduced as follow:
835           
836            q_star = 2*(pi*contrast)**2* volume( 1- volume)
837            for k = 10^(-8)*q_star/(2*(pi*|contrast|)**2)
838            we get 2 values of volume:
839                 with   1 - 4 * k >= 0
840                 volume1 = (1- sqrt(1- 4*k))/2
841                 volume2 = (1+ sqrt(1- 4*k))/2
842           
843            q_star: the invariant value included extrapolation is applied
844                         unit  1/A^(3)*1/cm
845                    q_star = self.get_qstar()
846                   
847            the result returned will be 0<= volume <= 1
848           
849            @param contrast: contrast value provides by the user of type float.
850                     contrast unit is 1/A^(2)= 10^(16)cm^(2)
851            @return: volume fraction
852            @note: volume fraction must have no unit
853        """
854        if contrast is None:
855            return 
856        if contrast < 0:
857            raise ValueError, "contrast must be greater than zero" 
858       
859        if self._qstar is None:
860            self._qstar = self.get_qstar()
861       
862        if self._qstar < 0:
863            raise RuntimeError, "invariant must be greater than zero"
864       
865        # Compute intermediate constant
866        k =  1.e-8 * self._qstar/(2 * (math.pi * math.fabs(float(contrast)))**2)
867        #Check discriminant value
868        discrim = 1 - 4 * k
869       
870        # Compute volume fraction
871        if discrim < 0:
872            raise RuntimeError, "could not compute the volume fraction: negative discriminant"
873        elif discrim == 0:
874            return 1/2
875        else:
876            volume1 = 0.5 * (1 - math.sqrt(discrim))
877            volume2 = 0.5 * (1 + math.sqrt(discrim))
878           
879            if 0 <= volume1 and volume1 <= 1:
880                return volume1
881            elif 0 <= volume2 and volume2 <= 1: 
882                return volume2
883            raise RuntimeError, "could not compute the volume fraction: inconsistent results"
884   
885    def get_qstar_with_error(self, extrapolation=None):
886        """
887            Compute the invariant uncertainty.
888            This uncertainty computation depends on whether or not the data is
889            smeared.
890            @return: invariant,  the invariant uncertainty
891                return self._get_qstar(), self._get_qstar_smear_uncertainty()   
892        """
893        if self._qstar is None:
894            self._qstar = self.get_qstar(extrapolation=extrapolation)
895        if self._qstar_err is None:
896            self._qstar_err = self._get_qstar_smear_uncertainty()
897           
898        return self._qstar, self._qstar_err
899   
900    def get_volume_fraction_with_error(self, contrast):
901        """
902            Compute uncertainty on volume value as well as the volume fraction
903            This uncertainty is given by the following equation:
904            dV = 0.5 * (4*k* dq_star) /(2* math.sqrt(1-k* q_star))
905                                 
906            for k = 10^(-8)*q_star/(2*(pi*|contrast|)**2)
907           
908            q_star: the invariant value including extrapolated value if existing
909            dq_star: the invariant uncertainty
910            dV: the volume uncertainty
911            @param contrast: contrast value
912            @return: V, dV = self.get_volume_fraction_with_error(contrast), dV
913        """
914        if contrast is None:
915            return 
916        self._qstar, self._qstar_err = self.get_qstar_with_error()
917       
918        volume = self.get_volume_fraction(contrast)
919        if self._qstar < 0:
920            raise ValueError, "invariant must be greater than zero"
921       
922        k =  1.e-8 * self._qstar /(2 * (math.pi* math.fabs(float(contrast)))**2)
923        #check value inside the sqrt function
924        value = 1 - k * self._qstar
925        if (value) <= 0:
926            raise ValueError, "Cannot compute incertainty on volume"
927        # Compute uncertainty
928        uncertainty = (0.5 * 4 * k * self._qstar_err)/(2 * math.sqrt(1 - k * self._qstar))
929       
930        return volume, math.fabs(uncertainty)
931   
932    def get_surface_with_error(self, contrast, porod_const):
933        """
934            Compute uncertainty of the surface value as well as thesurface value
935            this uncertainty is given as follow:
936           
937            dS = porod_const *2*pi[( dV -2*V*dV)/q_star
938                 + dq_star(v-v**2)
939                 
940            q_star: the invariant value including extrapolated value if existing
941            dq_star: the invariant uncertainty
942            V: the volume fraction value
943            dV: the volume uncertainty
944           
945            @param contrast: contrast value
946            @param porod_const: porod constant value
947            @return S, dS: the surface, with its uncertainty
948        """
949        if contrast is None or porod_const is None:
950            return 
951        v, dv = self.get_volume_fraction_with_error(contrast)
952        self._qstar, self._qstar_err = self.get_qstar_with_error()
953        if self._qstar <= 0:
954            raise ValueError, "invariant must be greater than zero"
955        ds = porod_const * 2 * math.pi * (( dv - 2 * v * dv)/ self._qstar\
956                 + self._qstar_err * ( v - v**2))
957        s = self.get_surface(contrast=contrast, porod_const=porod_const)
958        return s, ds
Note: See TracBrowser for help on using the repository browser.