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

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 8c44b2c was 8c44b2c, checked in by Jessica Tumarkin <jtumarki@…>, 13 years ago
  • Property mode set to 100644
File size: 26.9 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                                                             warning=True)
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        return out, cov_x
336   
337    def pr_fit(self, nfunc=5):
338        """
339        This is a direct fit to a given P(r). It assumes that the y data
340        is set to some P(r) distribution that we are trying to reproduce
341        with a set of base functions.
342       
343        This method is provided as a test.
344        """
345        #from scipy import optimize
346       
347        # First, check that the current data is valid
348        if self.is_valid() <= 0:
349            msg = "Invertor.invert: Data arrays are of different length"
350            raise RuntimeError, msg
351       
352        p = numpy.ones(nfunc)
353        t_0 = time.time()
354        out, cov_x, info, mesg, success = optimize.leastsq(self.pr_residuals, p,
355                                                            full_output=1,
356                                                            warning=True)
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            A[i][j] = (Fourier transformed base function for point j)
405           
406        We them choose a number of r-points, n_r, to evaluate the second
407        derivative of P(r) at. This is used as our regularization term.
408        For a vector r of length n_r, the following n_r rows are set to
409            A[i+Npts][j] = (2nd derivative of P(r), d**2(P(r))/d(r)**2,
410            evaluated at r[j])
411           
412        The vector b has its first Npts entries set to
413            b[j] = (I(q) observed for point j)
414           
415        The following n_r entries are set to zero.
416       
417        The result is found by using scipy.linalg.basic.lstsq to invert
418        the matrix and find the coefficients x.
419       
420        :param nfunc: number of base functions to use.
421        :param nr: number of r points to evaluate the 2nd derivative at
422            for the reg. term.
423
424        If the result does not allow us to compute the covariance matrix,
425        a matrix filled with zeros will be returned.
426
427        """
428        # Note: To make sure an array is contiguous:
429        # blah = numpy.ascontiguousarray(blah_original)
430        # ... before passing it to C
431       
432        if self.is_valid() < 0:
433            msg = "Invertor: invalid data; incompatible data lengths."
434            raise RuntimeError, msg
435       
436        self.nfunc = nfunc
437        # a -- An M x N matrix.
438        # b -- An M x nrhs matrix or M vector.
439        npts = len(self.x)
440        nq   = nr
441        sqrt_alpha = math.sqrt(math.fabs(self.alpha))
442        if sqrt_alpha < 0.0:
443            nq = 0
444
445        # If we need to fit the background, add a term
446        if self.has_bck == True:
447            nfunc_0 = nfunc
448            nfunc += 1
449
450        a = numpy.zeros([npts+nq, nfunc])
451        b = numpy.zeros(npts+nq)
452        err = numpy.zeros([nfunc, nfunc])
453       
454        # Construct the a matrix and b vector that represent the problem
455        t_0 = time.time()
456        self._get_matrix(nfunc, nq, a, b)
457             
458        # Perform the inversion (least square fit)
459        c, chi2, rank, n = lstsq(a, b)
460        # Sanity check
461        try:
462            float(chi2)
463        except:
464            chi2 = -1.0
465        self.chi2 = chi2
466               
467        inv_cov = numpy.zeros([nfunc,nfunc])
468        # Get the covariance matrix, defined as inv_cov = a_transposed * a
469        self._get_invcov_matrix(nfunc, nr, a, inv_cov)
470                   
471        # Compute the reg term size for the output
472        sum_sig, sum_reg = self._get_reg_size(nfunc, nr, a)
473                   
474        if math.fabs(self.alpha) > 0:
475            new_alpha = sum_sig/(sum_reg/self.alpha)
476        else:
477            new_alpha = 0.0
478        self.suggested_alpha = new_alpha
479       
480        try:
481            cov = numpy.linalg.pinv(inv_cov)
482            err = math.fabs(chi2/float(npts-nfunc)) * cov
483        except:
484            # We were not able to estimate the errors
485            # Return an empty error matrix
486            pass
487           
488        # Keep a copy of the last output
489        if self.has_bck == False:
490            self.background = 0
491            self.out = c
492            self.cov = err
493        else:
494            self.background = c[0]
495           
496            err_0 = numpy.zeros([nfunc, nfunc])
497            c_0 = numpy.zeros(nfunc)
498           
499            for i in range(nfunc_0):
500                c_0[i] = c[i+1]
501                for j in range(nfunc_0):
502                    err_0[i][j] = err[i+1][j+1]
503                   
504            self.out = c_0
505            self.cov = err_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, message, elapsed =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, cov = 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, cov = 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, cov = 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                    messsage += "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       
619    def to_file(self, path, npts=100):
620        """
621        Save the state to a file that will be readable
622        by SliceView.
623       
624        :param path: path of the file to write
625        :param npts: number of P(r) points to be written
626       
627        """
628        file = open(path, 'w')
629        file.write("#d_max=%g\n" % self.d_max)
630        file.write("#nfunc=%g\n" % self.nfunc)
631        file.write("#alpha=%g\n" % self.alpha)
632        file.write("#chi2=%g\n" % self.chi2)
633        file.write("#elapsed=%g\n" % self.elapsed)
634        file.write("#qmin=%s\n" % str(self.q_min))
635        file.write("#qmax=%s\n" % str(self.q_max))
636        file.write("#slit_height=%g\n" % self.slit_height)
637        file.write("#slit_width=%g\n" % self.slit_width)
638        file.write("#background=%g\n" % self.background)
639        if self.has_bck == True:
640            file.write("#has_bck=1\n")
641        else:
642            file.write("#has_bck=0\n")
643        file.write("#alpha_estimate=%g\n" % self.suggested_alpha)
644        if not self.out == None:
645            if len(self.out) == len(self.cov):
646                for i in range(len(self.out)):
647                    file.write("#C_%i=%s+-%s\n" % (i, str(self.out[i]),
648                                                    str(self.cov[i][i])))
649        file.write("<r>  <Pr>  <dPr>\n")
650        r = numpy.arange(0.0, self.d_max, self.d_max/npts)
651       
652        for r_i in r:
653            (value, err) = self.pr_err(self.out, self.cov, r_i)
654            file.write("%g  %g  %g\n" % (r_i, value, err))
655   
656        file.close()
657     
658       
659    def from_file(self, path):
660        """
661        Load the state of the Invertor from a file,
662        to be able to generate P(r) from a set of
663        parameters.
664       
665        :param path: path of the file to load
666       
667        """
668        #import os
669        #import re
670        if os.path.isfile(path):
671            try:
672                fd = open(path, 'r')
673               
674                buff    = fd.read()
675                lines   = buff.split('\n')
676                for line in lines:
677                    if line.startswith('#d_max='):
678                        toks = line.split('=')
679                        self.d_max = float(toks[1])
680                    elif line.startswith('#nfunc='):
681                        toks = line.split('=')
682                        self.nfunc = int(toks[1])
683                        self.out = numpy.zeros(self.nfunc)
684                        self.cov = numpy.zeros([self.nfunc, self.nfunc])
685                    elif line.startswith('#alpha='):
686                        toks = line.split('=')
687                        self.alpha = float(toks[1])
688                    elif line.startswith('#chi2='):
689                        toks = line.split('=')
690                        self.chi2 = float(toks[1])
691                    elif line.startswith('#elapsed='):
692                        toks = line.split('=')
693                        self.elapsed = float(toks[1])
694                    elif line.startswith('#alpha_estimate='):
695                        toks = line.split('=')
696                        self.suggested_alpha = float(toks[1])
697                    elif line.startswith('#qmin='):
698                        toks = line.split('=')
699                        try:
700                            self.q_min = float(toks[1])
701                        except:
702                            self.q_min = None
703                    elif line.startswith('#qmax='):
704                        toks = line.split('=')
705                        try:
706                            self.q_max = float(toks[1])
707                        except:
708                            self.q_max = None
709                    elif line.startswith('#slit_height='):
710                        toks = line.split('=')
711                        self.slit_height = float(toks[1])
712                    elif line.startswith('#slit_width='):
713                        toks = line.split('=')
714                        self.slit_width = float(toks[1])
715                    elif line.startswith('#background='):
716                        toks = line.split('=')
717                        self.background = float(toks[1])
718                    elif line.startswith('#has_bck='):
719                        toks = line.split('=')
720                        if int(toks[1]) == 1:
721                            self.has_bck = True
722                        else:
723                            self.has_bck = False
724           
725                    # Now read in the parameters
726                    elif line.startswith('#C_'):
727                        toks = line.split('=')
728                        p = re.compile('#C_([0-9]+)')
729                        m = p.search(toks[0])
730                        toks2 = toks[1].split('+-')
731                        i = int(m.group(1))
732                        self.out[i] = float(toks2[0])
733                       
734                        self.cov[i][i] = float(toks2[1])                       
735           
736            except:
737                msg = "Invertor.from_file: corrupted file\n%s" % sys.exc_value
738                raise RuntimeError, msg
739        else:
740            msg = "Invertor.from_file: '%s' is not a file" % str(path) 
741            raise RuntimeError, msg
742       
743         
744if __name__ == "__main__":
745    o = Invertor()
746
747   
748   
749   
750   
Note: See TracBrowser for help on using the repository browser.