source: sasview/src/sans/pr/invertor.py @ 51f14603

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 51f14603 was 51f14603, checked in by Peter Parker, 10 years ago

Refs #202 - Fix Sphinx build errors (not including park-1.2.1/). Most warnings remain.

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