source: sasview/pr_inversion/src/sans/pr/invertor.py @ f425805

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

Re #5 improving tests

  • Property mode set to 100644
File size: 27.0 KB
Line 
1"""
2Module to perform P(r) inversion.
3The module contains the Invertor class.
4"""
5
6import numpy
7import sys
8import math
9import time
10import copy
11import os
12import re
13from numpy.linalg import lstsq
14from scipy import optimize
15from sans.pr.core.pr_inversion import Cinvertor
16
17def help():
18    """
19    Provide general online help text
20    Future work: extend this function to allow topic selection
21    """
22    info_txt  = "The inversion approach is based on Moore, J. Appl. Cryst. " 
23    info_txt += "(1980) 13, 168-175.\n\n"
24    info_txt += "P(r) is set to be equal to an expansion of base functions "
25    info_txt += "of the type "
26    info_txt += "phi_n(r) = 2*r*sin(pi*n*r/D_max). The coefficient of each "
27    info_txt += "base functions "
28    info_txt += "in the expansion is found by performing a least square fit "
29    info_txt += "with the "
30    info_txt += "following fit function:\n\n"
31    info_txt += "chi**2 = sum_i[ I_meas(q_i) - I_th(q_i) ]**2/error**2 +"
32    info_txt += "Reg_term\n\n"
33    info_txt += "where I_meas(q) is the measured scattering intensity and "
34    info_txt += "I_th(q) is "
35    info_txt += "the prediction from the Fourier transform of the P(r) "
36    info_txt += "expansion. "
37    info_txt += "The Reg_term term is a regularization term set to the second"
38    info_txt += " derivative "
39    info_txt += "d**2P(r)/dr**2 integrated over r. It is used to produce "
40    info_txt += "a smooth P(r) output.\n\n"
41    info_txt += "The following are user inputs:\n\n"
42    info_txt += "   - Number of terms: the number of base functions in the P(r)"
43    info_txt += " expansion.\n\n"
44    info_txt += "   - Regularization constant: a multiplicative constant "
45    info_txt += "to set the size of "
46    info_txt += "the regularization term.\n\n"
47    info_txt += "   - Maximum distance: the maximum distance between any "
48    info_txt += "two points in the system.\n"
49     
50    return info_txt
51   
52
53class Invertor(Cinvertor):
54    """
55    Invertor class to perform P(r) inversion
56   
57    The problem is solved by posing the problem as  Ax = b,
58    where x is the set of coefficients we are looking for.
59   
60    Npts is the number of points.
61   
62    In the following i refers to the ith base function coefficient.
63    The matrix has its entries j in its first Npts rows set to
64        A[j][i] = (Fourier transformed base function for point j)
65       
66    We them choose a number of r-points, n_r, to evaluate the second
67    derivative of P(r) at. This is used as our regularization term.
68    For a vector r of length n_r, the following n_r rows are set to
69        A[j+Npts][i] = (2nd derivative of P(r), d**2(P(r))/d(r)**2,
70        evaluated at r[j])
71       
72    The vector b has its first Npts entries set to
73        b[j] = (I(q) observed for point j)
74       
75    The following n_r entries are set to zero.
76   
77    The result is found by using scipy.linalg.basic.lstsq to invert
78    the matrix and find the coefficients x.
79   
80    Methods inherited from Cinvertor:
81    - get_peaks(pars): returns the number of P(r) peaks
82    - oscillations(pars): returns the oscillation parameters for the output P(r)
83    - get_positive(pars): returns the fraction of P(r) that is above zero
84    - get_pos_err(pars): returns the fraction of P(r) that is 1-sigma above zero
85    """
86    ## Chisqr of the last computation
87    chi2  = 0
88    ## Time elapsed for last computation
89    elapsed = 0
90    ## Alpha to get the reg term the same size as the signal
91    suggested_alpha = 0
92    ## Last number of base functions used
93    nfunc = 10
94    ## Last output values
95    out = None
96    ## Last errors on output values
97    cov = None
98    ## Background value
99    background = 0
100    ## Information dictionary for application use
101    info = {}
102   
103   
104    def __init__(self):
105        Cinvertor.__init__(self)
106       
107    def __setstate__(self, state):
108        """
109        restore the state of invertor for pickle
110        """
111        (self.__dict__,self.alpha, self.d_max,
112         self.q_min, self.q_max, 
113         self.x, self.y,
114         self.err, self.has_bck,
115         self.slit_height, self.slit_width) = state
116   
117    def __reduce_ex__(self, proto):
118        """
119        Overwrite the __reduce_ex__
120        """
121
122        state = (self.__dict__, 
123                 self.alpha, self.d_max,
124                 self.q_min, self.q_max,
125                 self.x, self.y,
126                 self.err, self.has_bck, 
127                 self.slit_height, self.slit_width,
128                 )
129        return (Invertor,tuple(), state, None, None)
130   
131    def __setattr__(self, name, value):
132        """
133        Set the value of an attribute.
134        Access the parent class methods for
135        x, y, err, d_max, q_min, q_max and alpha
136        """
137        if   name == 'x':
138            if 0.0 in value:
139                msg = "Invertor: one of your q-values is zero. "
140                msg += "Delete that entry before proceeding"
141                raise ValueError, msg
142            return self.set_x(value)
143        elif name == 'y':
144            return self.set_y(value)
145        elif name == 'err':
146            value2 = abs(value)
147            return self.set_err(value2)
148        elif name == 'd_max':
149            return self.set_dmax(value)
150        elif name == 'q_min':
151            if value == None:
152                return self.set_qmin(-1.0)
153            return self.set_qmin(value)
154        elif name == 'q_max':
155            if value == None:
156                return self.set_qmax(-1.0)
157            return self.set_qmax(value)
158        elif name == 'alpha':
159            return self.set_alpha(value)
160        elif name == 'slit_height':
161            return self.set_slit_height(value)
162        elif name == 'slit_width':
163            return self.set_slit_width(value)
164        elif name == 'has_bck':
165            if value == True:
166                return self.set_has_bck(1)
167            elif value == False:
168                return self.set_has_bck(0)
169            else:
170                raise ValueError, "Invertor: has_bck can only be True or False"
171           
172        return Cinvertor.__setattr__(self, name, value)
173   
174    def __getattr__(self, name):
175        """
176        Return the value of an attribute
177        """
178        #import numpy
179        if name == 'x':
180            out = numpy.ones(self.get_nx())
181            self.get_x(out)
182            return out
183        elif name == 'y':
184            out = numpy.ones(self.get_ny())
185            self.get_y(out)
186            return out
187        elif name == 'err':
188            out = numpy.ones(self.get_nerr())
189            self.get_err(out)
190            return out
191        elif name == 'd_max':
192            return self.get_dmax()
193        elif name == 'q_min':
194            qmin = self.get_qmin()
195            if qmin < 0:
196                return None
197            return qmin
198        elif name == 'q_max':
199            qmax = self.get_qmax()
200            if qmax < 0:
201                return None
202            return qmax
203        elif name == 'alpha':
204            return self.get_alpha()
205        elif name == 'slit_height':
206            return self.get_slit_height()
207        elif name == 'slit_width':
208            return self.get_slit_width()
209        elif name == 'has_bck':
210            value = self.get_has_bck()
211            if value == 1:
212                return True
213            else:
214                return False
215        elif name in self.__dict__:
216            return self.__dict__[name]
217        return None
218   
219    def clone(self):
220        """
221        Return a clone of this instance
222        """
223        #import copy
224       
225        invertor = Invertor()
226        invertor.chi2    = self.chi2
227        invertor.elapsed = self.elapsed
228        invertor.nfunc   = self.nfunc
229        invertor.alpha   = self.alpha
230        invertor.d_max   = self.d_max
231        invertor.q_min   = self.q_min
232        invertor.q_max   = self.q_max
233       
234        invertor.x = self.x
235        invertor.y = self.y
236        invertor.err = self.err
237        invertor.has_bck = self.has_bck
238        invertor.slit_height = self.slit_height
239        invertor.slit_width  = self.slit_width
240       
241        invertor.info = copy.deepcopy(self.info)
242       
243        return invertor
244   
245    def invert(self, nfunc=10, nr=20):
246        """
247        Perform inversion to P(r)
248       
249        The problem is solved by posing the problem as  Ax = b,
250        where x is the set of coefficients we are looking for.
251       
252        Npts is the number of points.
253       
254        In the following i refers to the ith base function coefficient.
255        The matrix has its entries j in its first Npts rows set to
256            A[i][j] = (Fourier transformed base function for point j)
257           
258        We them choose a number of r-points, n_r, to evaluate the second
259        derivative of P(r) at. This is used as our regularization term.
260        For a vector r of length n_r, the following n_r rows are set to
261            A[i+Npts][j] = (2nd derivative of P(r), d**2(P(r))/d(r)**2, evaluated at r[j])
262           
263        The vector b has its first Npts entries set to
264            b[j] = (I(q) observed for point j)
265           
266        The following n_r entries are set to zero.
267       
268        The result is found by using scipy.linalg.basic.lstsq to invert
269        the matrix and find the coefficients x.
270       
271        :param nfunc: number of base functions to use.
272        :param nr: number of r points to evaluate the 2nd derivative at for the reg. term.
273        :return: c_out, c_cov - the coefficients with covariance matrix
274       
275        """
276        # Reset the background value before proceeding
277        self.background = 0.0
278        return self.lstsq(nfunc, nr=nr)
279   
280    def iq(self, out, q):
281        """
282        Function to call to evaluate the scattering intensity
283       
284        :param args: c-parameters, and q
285        :return: I(q)
286       
287        """
288        return Cinvertor.iq(self, out, q) + self.background
289   
290    def invert_optimize(self, nfunc=10, nr=20):
291        """
292        Slower version of the P(r) inversion that uses scipy.optimize.leastsq.
293       
294        This probably produce more reliable results, but is much slower.
295        The minimization function is set to
296        sum_i[ (I_obs(q_i) - I_theo(q_i))/err**2 ] + alpha * reg_term,
297        where the reg_term is given by Svergun: it is the integral of
298        the square of the first derivative
299        of P(r), d(P(r))/dr, integrated over the full range of r.
300       
301        :param nfunc: number of base functions to use.
302        :param nr: number of r points to evaluate the 2nd derivative at
303            for the reg. term.
304       
305        :return: c_out, c_cov - the coefficients with covariance matrix
306       
307        """
308       
309        #from scipy import optimize
310        #import time
311       
312        self.nfunc = nfunc
313        # First, check that the current data is valid
314        if self.is_valid() <= 0:
315            msg = "Invertor.invert: Data array are of different length"
316            raise RuntimeError, msg
317       
318        p = numpy.ones(nfunc)
319        t_0 = time.time()
320        out, cov_x, _, _, _ = optimize.leastsq(self.residuals,
321                                                            p, full_output=1)
322       
323        # Compute chi^2
324        res = self.residuals(out)
325        chisqr = 0
326        for i in range(len(res)):
327            chisqr += res[i]
328       
329        self.chi2 = chisqr
330
331        # Store computation time
332        self.elapsed = time.time() - t_0
333       
334        if cov_x is None:
335            cov_x = numpy.ones([nfunc,nfunc])
336            cov_x *= math.fabs(chisqr)
337        return out, cov_x
338   
339    def pr_fit(self, nfunc=5):
340        """
341        This is a direct fit to a given P(r). It assumes that the y data
342        is set to some P(r) distribution that we are trying to reproduce
343        with a set of base functions.
344       
345        This method is provided as a test.
346        """
347        #from scipy import optimize
348       
349        # First, check that the current data is valid
350        if self.is_valid() <= 0:
351            msg = "Invertor.invert: Data arrays are of different length"
352            raise RuntimeError, msg
353       
354        p = numpy.ones(nfunc)
355        t_0 = time.time()
356        out, cov_x, info, mesg, success = optimize.leastsq(self.pr_residuals, p,
357                                                            full_output=1)
358       
359        # Compute chi^2
360        res = self.pr_residuals(out)
361        chisqr = 0
362        for i in range(len(res)):
363            chisqr += res[i]
364       
365        self.chisqr = chisqr
366       
367        # Store computation time
368        self.elapsed = time.time() - t_0
369
370        return out, cov_x
371   
372    def pr_err(self, c, c_cov, r):
373        """   
374        Returns the value of P(r) for a given r, and base function
375        coefficients, with error.
376       
377        :param c: base function coefficients
378        :param c_cov: covariance matrice of the base function coefficients
379        :param r: r-value to evaluate P(r) at
380       
381        :return: P(r)
382       
383        """
384        return self.get_pr_err(c, c_cov, r)
385       
386    def _accept_q(self, q):
387        """
388        Check q-value against user-defined range
389        """
390        if not self.q_min == None and q < self.q_min:
391            return False
392        if not self.q_max == None and q > self.q_max:
393            return False
394        return True
395       
396    def lstsq(self, nfunc=5, nr=20):
397        """
398        The problem is solved by posing the problem as  Ax = b,
399        where x is the set of coefficients we are looking for.
400       
401        Npts is the number of points.
402       
403        In the following i refers to the ith base function coefficient.
404        The matrix has its entries j in its first Npts rows set to
405            A[i][j] = (Fourier transformed base function for point j)
406           
407        We them choose a number of r-points, n_r, to evaluate the second
408        derivative of P(r) at. This is used as our regularization term.
409        For a vector r of length n_r, the following n_r rows are set to
410            A[i+Npts][j] = (2nd derivative of P(r), d**2(P(r))/d(r)**2,
411            evaluated at r[j])
412           
413        The vector b has its first Npts entries set to
414            b[j] = (I(q) observed for point j)
415           
416        The following n_r entries are set to zero.
417       
418        The result is found by using scipy.linalg.basic.lstsq to invert
419        the matrix and find the coefficients x.
420       
421        :param nfunc: number of base functions to use.
422        :param nr: number of r points to evaluate the 2nd derivative at
423            for the reg. term.
424
425        If the result does not allow us to compute the covariance matrix,
426        a matrix filled with zeros will be returned.
427
428        """
429        # Note: To make sure an array is contiguous:
430        # blah = numpy.ascontiguousarray(blah_original)
431        # ... before passing it to C
432       
433        if self.is_valid() < 0:
434            msg = "Invertor: invalid data; incompatible data lengths."
435            raise RuntimeError, msg
436       
437        self.nfunc = nfunc
438        # a -- An M x N matrix.
439        # b -- An M x nrhs matrix or M vector.
440        npts = len(self.x)
441        nq   = nr
442        sqrt_alpha = math.sqrt(math.fabs(self.alpha))
443        if sqrt_alpha < 0.0:
444            nq = 0
445
446        # If we need to fit the background, add a term
447        if self.has_bck == True:
448            nfunc_0 = nfunc
449            nfunc += 1
450
451        a = numpy.zeros([npts+nq, nfunc])
452        b = numpy.zeros(npts+nq)
453        err = numpy.zeros([nfunc, nfunc])
454       
455        # Construct the a matrix and b vector that represent the problem
456        t_0 = time.time()
457        try:
458            self._get_matrix(nfunc, nq, a, b)
459        except:
460            raise RuntimeError, "Invertor: could not invert I(Q)\n  %s" % sys.exc_value
461             
462        # Perform the inversion (least square fit)
463        c, chi2, rank, n = lstsq(a, b)
464        # Sanity check
465        try:
466            float(chi2)
467        except:
468            chi2 = -1.0
469        self.chi2 = chi2
470               
471        inv_cov = numpy.zeros([nfunc,nfunc])
472        # Get the covariance matrix, defined as inv_cov = a_transposed * a
473        self._get_invcov_matrix(nfunc, nr, a, inv_cov)
474                   
475        # Compute the reg term size for the output
476        sum_sig, sum_reg = self._get_reg_size(nfunc, nr, a)
477                   
478        if math.fabs(self.alpha) > 0:
479            new_alpha = sum_sig/(sum_reg/self.alpha)
480        else:
481            new_alpha = 0.0
482        self.suggested_alpha = new_alpha
483       
484        try:
485            cov = numpy.linalg.pinv(inv_cov)
486            err = math.fabs(chi2/float(npts-nfunc)) * cov
487        except:
488            # We were not able to estimate the errors
489            # Return an empty error matrix
490            pass
491           
492        # Keep a copy of the last output
493        if self.has_bck == False:
494            self.background = 0
495            self.out = c
496            self.cov = err
497        else:
498            self.background = c[0]
499           
500            err_0 = numpy.zeros([nfunc, nfunc])
501            c_0 = numpy.zeros(nfunc)
502           
503            for i in range(nfunc_0):
504                c_0[i] = c[i+1]
505                for j in range(nfunc_0):
506                    err_0[i][j] = err[i+1][j+1]
507                   
508            self.out = c_0
509            self.cov = err_0
510           
511        return self.out, self.cov
512       
513    def estimate_numterms(self, isquit_func=None):
514        """
515        Returns a reasonable guess for the
516        number of terms
517       
518        :param isquit_func: reference to thread function to call to
519                            check whether the computation needs to
520                            be stopped.
521       
522        :return: number of terms, alpha, message
523       
524        """
525        from num_term import Num_terms
526        estimator = Num_terms(self.clone())
527        try:
528            return estimator.num_terms(isquit_func)
529        except:
530            # If we fail, estimate alpha and return the default
531            # number of terms
532            best_alpha, message, elapsed =self.estimate_alpha(self.nfunc)
533            return self.nfunc, best_alpha, "Could not estimate number of terms"
534                   
535    def estimate_alpha(self, nfunc):
536        """
537        Returns a reasonable guess for the
538        regularization constant alpha
539       
540        :param nfunc: number of terms to use in the expansion.
541       
542        :return: alpha, message, elapsed
543       
544        where alpha is the estimate for alpha,
545        message is a message for the user,
546        elapsed is the computation time
547        """
548        #import time
549        try:           
550            pr = self.clone()
551           
552            # T_0 for computation time
553            starttime = time.time()
554            elapsed = 0
555           
556            # If the current alpha is zero, try
557            # another value
558            if pr.alpha <= 0:
559                pr.alpha = 0.0001
560                 
561            # Perform inversion to find the largest alpha
562            out, cov = pr.invert(nfunc)
563            elapsed = time.time() - starttime
564            initial_alpha = pr.alpha
565            initial_peaks = pr.get_peaks(out)
566   
567            # Try the inversion with the estimated alpha
568            pr.alpha = pr.suggested_alpha
569            out, cov = pr.invert(nfunc)
570   
571            npeaks = pr.get_peaks(out)
572            # if more than one peak to start with
573            # just return the estimate
574            if npeaks > 1:
575                #message = "Your P(r) is not smooth,
576                #please check your inversion parameters"
577                message = None
578                return pr.suggested_alpha, message, elapsed
579            else:
580               
581                # Look at smaller values
582                # We assume that for the suggested alpha, we have 1 peak
583                # if not, send a message to change parameters
584                alpha = pr.suggested_alpha
585                best_alpha = pr.suggested_alpha
586                found = False
587                for i in range(10):
588                    pr.alpha = (0.33)**(i+1) * alpha
589                    out, cov = pr.invert(nfunc)
590                   
591                    peaks = pr.get_peaks(out)
592                    if peaks > 1:
593                        found = True
594                        break
595                    best_alpha = pr.alpha
596                   
597                # If we didn't find a turning point for alpha and
598                # the initial alpha already had only one peak,
599                # just return that
600                if not found and initial_peaks == 1 and \
601                    initial_alpha<best_alpha:
602                    best_alpha = initial_alpha
603                   
604                # Check whether the size makes sense
605                message = ''
606               
607                if not found:
608                    message = None
609                elif best_alpha >= 0.5 * pr.suggested_alpha:
610                    # best alpha is too big, return a
611                    # reasonable value
612                    message  = "The estimated alpha for your system is too "
613                    messsage += "large. "
614                    message += "Try increasing your maximum distance."
615               
616                return best_alpha, message, elapsed
617   
618        except:
619            message = "Invertor.estimate_alpha: %s" % sys.exc_value
620            return 0, message, elapsed
621   
622       
623    def to_file(self, path, npts=100):
624        """
625        Save the state to a file that will be readable
626        by SliceView.
627       
628        :param path: path of the file to write
629        :param npts: number of P(r) points to be written
630       
631        """
632        file = open(path, 'w')
633        file.write("#d_max=%g\n" % self.d_max)
634        file.write("#nfunc=%g\n" % self.nfunc)
635        file.write("#alpha=%g\n" % self.alpha)
636        file.write("#chi2=%g\n" % self.chi2)
637        file.write("#elapsed=%g\n" % self.elapsed)
638        file.write("#qmin=%s\n" % str(self.q_min))
639        file.write("#qmax=%s\n" % str(self.q_max))
640        file.write("#slit_height=%g\n" % self.slit_height)
641        file.write("#slit_width=%g\n" % self.slit_width)
642        file.write("#background=%g\n" % self.background)
643        if self.has_bck == True:
644            file.write("#has_bck=1\n")
645        else:
646            file.write("#has_bck=0\n")
647        file.write("#alpha_estimate=%g\n" % self.suggested_alpha)
648        if not self.out == None:
649            if len(self.out) == len(self.cov):
650                for i in range(len(self.out)):
651                    file.write("#C_%i=%s+-%s\n" % (i, str(self.out[i]),
652                                                    str(self.cov[i][i])))
653        file.write("<r>  <Pr>  <dPr>\n")
654        r = numpy.arange(0.0, self.d_max, self.d_max/npts)
655       
656        for r_i in r:
657            (value, err) = self.pr_err(self.out, self.cov, r_i)
658            file.write("%g  %g  %g\n" % (r_i, value, err))
659   
660        file.close()
661     
662       
663    def from_file(self, path):
664        """
665        Load the state of the Invertor from a file,
666        to be able to generate P(r) from a set of
667        parameters.
668       
669        :param path: path of the file to load
670       
671        """
672        #import os
673        #import re
674        if os.path.isfile(path):
675            try:
676                fd = open(path, 'r')
677               
678                buff    = fd.read()
679                lines   = buff.split('\n')
680                for line in lines:
681                    if line.startswith('#d_max='):
682                        toks = line.split('=')
683                        self.d_max = float(toks[1])
684                    elif line.startswith('#nfunc='):
685                        toks = line.split('=')
686                        self.nfunc = int(toks[1])
687                        self.out = numpy.zeros(self.nfunc)
688                        self.cov = numpy.zeros([self.nfunc, self.nfunc])
689                    elif line.startswith('#alpha='):
690                        toks = line.split('=')
691                        self.alpha = float(toks[1])
692                    elif line.startswith('#chi2='):
693                        toks = line.split('=')
694                        self.chi2 = float(toks[1])
695                    elif line.startswith('#elapsed='):
696                        toks = line.split('=')
697                        self.elapsed = float(toks[1])
698                    elif line.startswith('#alpha_estimate='):
699                        toks = line.split('=')
700                        self.suggested_alpha = float(toks[1])
701                    elif line.startswith('#qmin='):
702                        toks = line.split('=')
703                        try:
704                            self.q_min = float(toks[1])
705                        except:
706                            self.q_min = None
707                    elif line.startswith('#qmax='):
708                        toks = line.split('=')
709                        try:
710                            self.q_max = float(toks[1])
711                        except:
712                            self.q_max = None
713                    elif line.startswith('#slit_height='):
714                        toks = line.split('=')
715                        self.slit_height = float(toks[1])
716                    elif line.startswith('#slit_width='):
717                        toks = line.split('=')
718                        self.slit_width = float(toks[1])
719                    elif line.startswith('#background='):
720                        toks = line.split('=')
721                        self.background = float(toks[1])
722                    elif line.startswith('#has_bck='):
723                        toks = line.split('=')
724                        if int(toks[1]) == 1:
725                            self.has_bck = True
726                        else:
727                            self.has_bck = False
728           
729                    # Now read in the parameters
730                    elif line.startswith('#C_'):
731                        toks = line.split('=')
732                        p = re.compile('#C_([0-9]+)')
733                        m = p.search(toks[0])
734                        toks2 = toks[1].split('+-')
735                        i = int(m.group(1))
736                        self.out[i] = float(toks2[0])
737                       
738                        self.cov[i][i] = float(toks2[1])                       
739           
740            except:
741                msg = "Invertor.from_file: corrupted file\n%s" % sys.exc_value
742                raise RuntimeError, msg
743        else:
744            msg = "Invertor.from_file: '%s' is not a file" % str(path) 
745            raise RuntimeError, msg
746       
747         
748if __name__ == "__main__":
749    o = Invertor()
750
751   
752   
753   
754   
Note: See TracBrowser for help on using the repository browser.