source: sasview/src/sans/pr/invertor.py @ 5777106

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

Moving things around. Will definitely not build.

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