root/library/bdm/stat/exp_family.h @ 461

Revision 461, 34.8 kB (checked in by vbarta, 15 years ago)

mpdf (& its dependencies) reformat: now using shared_ptr, moved virtual method bodies to .cpp

  • Property svn:eol-style set to native
Line 
1/*!
2  \file
3  \brief Probability distributions for Exponential Family models.
4  \author Vaclav Smidl.
5
6  -----------------------------------
7  BDM++ - C++ library for Bayesian Decision Making under Uncertainty
8
9  Using IT++ for numerical operations
10  -----------------------------------
11*/
12
13#ifndef EF_H
14#define EF_H
15
16
17#include "../shared_ptr.h"
18#include "../base/bdmbase.h"
19#include "../math/chmat.h"
20
21namespace bdm
22{
23
24
25//! Global Uniform_RNG
26        extern Uniform_RNG UniRNG;
27//! Global Normal_RNG
28        extern Normal_RNG NorRNG;
29//! Global Gamma_RNG
30        extern Gamma_RNG GamRNG;
31
32        /*!
33        * \brief General conjugate exponential family posterior density.
34
35        * More?...
36        */
37
38        class eEF : public epdf
39        {
40                public:
41//      eEF() :epdf() {};
42                        //! default constructor
43                        eEF ( ) :epdf ( ) {};
44                        //! logarithm of the normalizing constant, \f$\mathcal{I}\f$
45                        virtual double lognc() const =0;
46                        //!TODO decide if it is really needed
47                        virtual void dupdate ( mat &v ) {it_error ( "Not implemented" );};
48                        //!Evaluate normalized log-probability
49                        virtual double evallog_nn ( const vec &val ) const{it_error ( "Not implemented" );return 0.0;};
50                        //!Evaluate normalized log-probability
51                        virtual double evallog ( const vec &val ) const {
52                                double tmp;
53                                tmp= evallog_nn ( val )-lognc();
54//                              it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
55                                return tmp;}
56                        //!Evaluate normalized log-probability for many samples
57                        virtual vec evallog ( const mat &Val ) const
58                        {
59                                vec x ( Val.cols() );
60                                for ( int i=0;i<Val.cols();i++ ) {x ( i ) =evallog_nn ( Val.get_col ( i ) ) ;}
61                                return x-lognc();
62                        }
63                        //!Power of the density, used e.g. to flatten the density
64                        virtual void pow ( double p ) {it_error ( "Not implemented" );};
65        };
66
67        /*!
68        * \brief Exponential family model.
69
70        * More?...
71        */
72
73        class mEF : public mpdf
74        {
75
76                public:
77                        //! Default constructor
78                        mEF ( ) :mpdf ( ) {};
79        };
80
81//! Estimator for Exponential family
82        class BMEF : public BM
83        {
84                protected:
85                        //! forgetting factor
86                        double frg;
87                        //! cached value of lognc() in the previous step (used in evaluation of \c ll )
88                        double last_lognc;
89                public:
90                        //! Default constructor (=empty constructor)
91                        BMEF ( double frg0=1.0 ) :BM ( ), frg ( frg0 ) {}
92                        //! Copy constructor
93                        BMEF ( const BMEF &B ) :BM ( B ), frg ( B.frg ), last_lognc ( B.last_lognc ) {}
94                        //!get statistics from another model
95                        virtual void set_statistics ( const BMEF* BM0 ) {it_error ( "Not implemented" );};
96                        //! Weighted update of sufficient statistics (Bayes rule)
97                        virtual void bayes ( const vec &data, const double w ) {};
98                        //original Bayes
99                        void bayes ( const vec &dt );
100                        //!Flatten the posterior according to the given BMEF (of the same type!)
101                        virtual void flatten ( const BMEF * B ) {it_error ( "Not implemented" );}
102                        //!Flatten the posterior as if to keep nu0 data
103//      virtual void flatten ( double nu0 ) {it_error ( "Not implemented" );}
104
105                        BMEF* _copy_ () const {it_error ( "function _copy_ not implemented for this BM" ); return NULL;};
106        };
107
108        template<class sq_T>
109        class mlnorm;
110
111        /*!
112        * \brief Gaussian density with positive definite (decomposed) covariance matrix.
113
114        * More?...
115        */
116        template<class sq_T>
117        class enorm : public eEF
118        {
119                protected:
120                        //! mean value
121                        vec mu;
122                        //! Covariance matrix in decomposed form
123                        sq_T R;
124                public:
125                        //!\name Constructors
126                        //!@{
127
128                        enorm ( ) :eEF ( ), mu ( ),R ( ) {};
129                        enorm ( const vec &mu,const sq_T &R ) {set_parameters ( mu,R );}
130                        void set_parameters ( const vec &mu,const sq_T &R );
131                        void from_setting(const Setting &root);
132                        void validate() {
133                                it_assert(mu.length()==R.rows(),"parameters mismatch");
134                                dim = mu.length();
135                        }
136                        //!@}
137
138                        //! \name Mathematical operations
139                        //!@{
140
141                        //! dupdate in exponential form (not really handy)
142                        void dupdate ( mat &v,double nu=1.0 );
143
144                        vec sample() const;
145
146                        double evallog_nn ( const vec &val ) const;
147                        double lognc () const;
148                        vec mean() const {return mu;}
149                        vec variance() const {return diag ( R.to_mat() );}
150//      mlnorm<sq_T>* condition ( const RV &rvn ) const ; <=========== fails to cmpile. Why?
151                        mpdf* condition ( const RV &rvn ) const ;
152                        enorm<sq_T>* marginal ( const RV &rv ) const;
153//                      epdf* marginal ( const RV &rv ) const;
154                        //!@}
155
156                        //! \name Access to attributes
157                        //!@{
158
159                        vec& _mu() {return mu;}
160                        void set_mu ( const vec mu0 ) { mu=mu0;}
161                        sq_T& _R() {return R;}
162                        const sq_T& _R() const {return R;}
163                        //!@}
164
165        };
166        UIREGISTER(enorm<chmat>);
167        UIREGISTER(enorm<ldmat>);
168        UIREGISTER(enorm<fsqmat>);
169
170
171        /*!
172        * \brief Gauss-inverse-Wishart density stored in LD form
173
174        * For \f$p\f$-variate densities, given rv.count() should be \f$p\times\f$ V.rows().
175        *
176        */
177        class egiw : public eEF
178        {
179                protected:
180                        //! Extended information matrix of sufficient statistics
181                        ldmat V;
182                        //! Number of data records (degrees of freedom) of sufficient statistics
183                        double nu;
184                        //! Dimension of the output
185                        int dimx;
186                        //! Dimension of the regressor
187                        int nPsi;
188                public:
189                        //!\name Constructors
190                        //!@{
191                        egiw() :eEF() {};
192                        egiw ( int dimx0, ldmat V0, double nu0=-1.0 ) :eEF() {set_parameters ( dimx0,V0, nu0 );};
193
194                        void set_parameters ( int dimx0, ldmat V0, double nu0=-1.0 )
195                        {
196                                dimx=dimx0;
197                                nPsi = V0.rows()-dimx;
198                                dim = dimx* ( dimx+nPsi ); // size(R) + size(Theta)
199
200                                V=V0;
201                                if ( nu0<0 )
202                                {
203                                        nu = 0.1 +nPsi +2*dimx +2; // +2 assures finite expected value of R
204                                        // terms before that are sufficient for finite normalization
205                                }
206                                else
207                                {
208                                        nu=nu0;
209                                }
210                        }
211                        //!@}
212
213                        vec sample() const;
214                        vec mean() const;
215                        vec variance() const;
216
217                        //! LS estimate of \f$\theta\f$
218                        vec est_theta() const;
219
220                        //! Covariance of the LS estimate
221                        ldmat est_theta_cov() const;
222
223                        void mean_mat ( mat &M, mat&R ) const;
224                        //! In this instance, val= [theta, r]. For multivariate instances, it is stored columnwise val = [theta_1 theta_2 ... r_1 r_2 ]
225                        double evallog_nn ( const vec &val ) const;
226                        double lognc () const;
227                        void pow ( double p ) {V*=p;nu*=p;};
228
229                        //! \name Access attributes
230                        //!@{
231
232                        ldmat& _V() {return V;}
233                        const ldmat& _V() const {return V;}
234                        double& _nu()  {return nu;}
235                        const double& _nu() const {return nu;}
236                        void from_setting(const Setting &set){
237                                set.lookupValue("nu",nu);
238                                set.lookupValue("dimx",dimx);
239                                mat V;
240                                UI::get(V,set,"V");
241                                set_parameters(dimx, V, nu);
242                                RV* rv=UI::build<RV>(set,"rv");
243                                set_rv(*rv);
244                                delete rv;
245                        }
246                        //!@}
247        };
248        UIREGISTER(egiw);
249
250        /*! \brief Dirichlet posterior density
251
252        Continuous Dirichlet density of \f$n\f$-dimensional variable \f$x\f$
253        \f[
254        f(x|\beta) = \frac{\Gamma[\gamma]}{\prod_{i=1}^{n}\Gamma(\beta_i)} \prod_{i=1}^{n}x_i^{\beta_i-1}
255        \f]
256        where \f$\gamma=\sum_i \beta_i\f$.
257        */
258        class eDirich: public eEF
259        {
260                protected:
261                        //!sufficient statistics
262                        vec beta;
263                public:
264                        //!\name Constructors
265                        //!@{
266
267                        eDirich () : eEF ( ) {};
268                        eDirich ( const eDirich &D0 ) : eEF () {set_parameters ( D0.beta );};
269                        eDirich ( const vec &beta0 ) {set_parameters ( beta0 );};
270                        void set_parameters ( const vec &beta0 )
271                        {
272                                beta= beta0;
273                                dim = beta.length();
274                        }
275                        //!@}
276
277                        vec sample() const {it_error ( "Not implemented" );return vec_1 ( 0.0 );};
278                        vec mean() const {return beta/sum(beta);};
279                        vec variance() const {double gamma =sum(beta); return elem_mult ( beta, ( beta+1 ) ) / ( gamma* ( gamma+1 ) );}
280                        //! In this instance, val is ...
281                        double evallog_nn ( const vec &val ) const
282                        {
283                                double tmp; tmp= ( beta-1 ) *log ( val );
284//                              it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
285                                return tmp;
286                        };
287                        double lognc () const
288                        {
289                                double tmp;
290                                double gam=sum ( beta );
291                                double lgb=0.0;
292                                for ( int i=0;i<beta.length();i++ ) {lgb+=lgamma ( beta ( i ) );}
293                                tmp= lgb-lgamma ( gam );
294//                              it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
295                                return tmp;
296                        };
297                        //!access function
298                        vec& _beta()  {return beta;}
299                        //!Set internal parameters
300        };
301
302//! \brief Estimator for Multinomial density
303        class multiBM : public BMEF
304        {
305                protected:
306                        //! Conjugate prior and posterior
307                        eDirich est;
308                        //! Pointer inside est to sufficient statistics
309                        vec &beta;
310                public:
311                        //!Default constructor
312                        multiBM ( ) : BMEF ( ),est ( ),beta ( est._beta() )
313                        {
314                                if ( beta.length() >0 ) {last_lognc=est.lognc();}
315                                else{last_lognc=0.0;}
316                        }
317                        //!Copy constructor
318                        multiBM ( const multiBM &B ) : BMEF ( B ),est ( B.est ),beta ( est._beta() ) {}
319                        //! Sets sufficient statistics to match that of givefrom mB0
320                        void set_statistics ( const BM* mB0 ) {const multiBM* mB=dynamic_cast<const multiBM*> ( mB0 ); beta=mB->beta;}
321                        void bayes ( const vec &dt )
322                        {
323                                if ( frg<1.0 ) {beta*=frg;last_lognc=est.lognc();}
324                                beta+=dt;
325                                if ( evalll ) {ll=est.lognc()-last_lognc;}
326                        }
327                        double logpred ( const vec &dt ) const
328                        {
329                                eDirich pred ( est );
330                                vec &beta = pred._beta();
331
332                                double lll;
333                                if ( frg<1.0 )
334                                        {beta*=frg;lll=pred.lognc();}
335                                else
336                                        if ( evalll ) {lll=last_lognc;}
337                                        else{lll=pred.lognc();}
338
339                                beta+=dt;
340                                return pred.lognc()-lll;
341                        }
342                        void flatten ( const BMEF* B )
343                        {
344                                const multiBM* E=dynamic_cast<const multiBM*> ( B );
345                                // sum(beta) should be equal to sum(B.beta)
346                                const vec &Eb=E->beta;//const_cast<multiBM*> ( E )->_beta();
347                                beta*= ( sum ( Eb ) /sum ( beta ) );
348                                if ( evalll ) {last_lognc=est.lognc();}
349                        }
350                        const epdf& posterior() const {return est;};
351                        const eDirich* _e() const {return &est;};
352                        void set_parameters ( const vec &beta0 )
353                        {
354                                est.set_parameters ( beta0 );
355                                if ( evalll ) {last_lognc=est.lognc();}
356                        }
357        };
358
359        /*!
360         \brief Gamma posterior density
361
362         Multivariate Gamma density as product of independent univariate densities.
363         \f[
364         f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
365         \f]
366        */
367
368        class egamma : public eEF
369        {
370                protected:
371                        //! Vector \f$\alpha\f$
372                        vec alpha;
373                        //! Vector \f$\beta\f$
374                        vec beta;
375                public :
376                        //! \name Constructors
377                        //!@{
378                        egamma ( ) :eEF ( ), alpha ( 0 ), beta ( 0 ) {};
379                        egamma ( const vec &a, const vec &b ) {set_parameters ( a, b );};
380                        void set_parameters ( const vec &a, const vec &b ) {alpha=a,beta=b;dim = alpha.length();};
381                        //!@}
382
383                        vec sample() const;
384                        //! TODO: is it used anywhere?
385//      mat sample ( int N ) const;
386                        double evallog ( const vec &val ) const;
387                        double lognc () const;
388                        //! Returns poiter to alpha and beta. Potentially dengerous: use with care!
389                        vec& _alpha() {return alpha;}
390                        vec& _beta() {return beta;}
391                        vec mean() const {return elem_div ( alpha,beta );}
392                        vec variance() const {return elem_div ( alpha,elem_mult ( beta,beta ) ); }
393                       
394                        //! Load from structure with elements:
395                        //!  \code
396                        //! { alpha = [...];         // vector of alpha
397                        //!   beta = [...];          // vector of beta
398                        //!   rv = {class="RV",...}  // description
399                        //! }
400                        //! \endcode
401                        //!@}
402                        void from_setting(const Setting &set){
403                                epdf::from_setting(set); // reads rv
404                                UI::get(alpha,set,"alpha");
405                                UI::get(beta,set,"beta");
406                                validate();
407                        }
408                        void validate(){
409                                it_assert(alpha.length() ==beta.length(), "parameters do not match");
410                                dim =alpha.length();
411                        }
412        };
413UIREGISTER(egamma);
414        /*!
415         \brief Inverse-Gamma posterior density
416
417         Multivariate inverse-Gamma density as product of independent univariate densities.
418         \f[
419         f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
420         \f]
421
422        Vector \f$\beta\f$ has different meaning (in fact it is 1/beta as used in definition of iG)
423
424         Inverse Gamma can be converted to Gamma using \f[
425         x\sim iG(a,b) => 1/x\sim G(a,1/b)
426        \f]
427        This relation is used in sampling.
428         */
429
430        class eigamma : public egamma
431        {
432                protected:
433                public :
434                        //! \name Constructors
435                        //! All constructors are inherited
436                        //!@{
437                        //!@}
438
439                        vec sample() const {return 1.0/egamma::sample();};
440                        //! Returns poiter to alpha and beta. Potentially dangerous: use with care!
441                        vec mean() const {return elem_div ( beta,alpha-1 );}
442                        vec variance() const {vec mea=mean(); return elem_div ( elem_mult ( mea,mea ),alpha-2 );}
443        };
444        /*
445        //! Weighted mixture of epdfs with external owned components.
446        class emix : public epdf {
447        protected:
448                int n;
449                vec &w;
450                Array<epdf*> Coms;
451        public:
452        //! Default constructor
453                emix ( const RV &rv, vec &w0): epdf(rv), n(w0.length()), w(w0), Coms(n) {};
454                void set_parameters( int &i, double wi, epdf* ep){w(i)=wi;Coms(i)=ep;}
455                vec mean(){vec pom; for(int i=0;i<n;i++){pom+=Coms(i)->mean()*w(i);} return pom;};
456                vec sample() {it_error ( "Not implemented" );return 0;}
457        };
458        */
459
460//!  Uniform distributed density on a rectangular support
461
462        class euni: public epdf
463        {
464                protected:
465//! lower bound on support
466                        vec low;
467//! upper bound on support
468                        vec high;
469//! internal
470                        vec distance;
471//! normalizing coefficients
472                        double nk;
473//! cache of log( \c nk )
474                        double lnk;
475                public:
476                        //! \name Constructors
477                        //!@{
478                        euni ( ) :epdf ( ) {}
479                        euni ( const vec &low0, const vec &high0 ) {set_parameters ( low0,high0 );}
480                        void set_parameters ( const vec &low0, const vec &high0 )
481                        {
482                                distance = high0-low0;
483                                it_assert_debug ( min ( distance ) >0.0,"bad support" );
484                                low = low0;
485                                high = high0;
486                                nk = prod ( 1.0/distance );
487                                lnk = log ( nk );
488                                dim = low.length();
489                        }
490                        //!@}
491
492                        double eval ( const vec &val ) const  {return nk;}
493                        double evallog ( const vec &val ) const  {
494                                if (any(val<low) && any(val>high)) {return inf;}
495                                else return lnk;
496                        }
497                        vec sample() const
498                        {
499                                vec smp ( dim );
500#pragma omp critical
501                                UniRNG.sample_vector ( dim ,smp );
502                                return low+elem_mult ( distance,smp );
503                        }
504                        //! set values of \c low and \c high
505                        vec mean() const {return ( high-low ) /2.0;}
506                        vec variance() const {return ( pow ( high,2 ) +pow ( low,2 ) +elem_mult ( high,low ) ) /3.0;}
507                        //! Load from structure with elements:
508                        //!  \code
509                        //! { high = [...];          // vector of upper bounds
510                        //!   low = [...];           // vector of lower bounds
511                        //!   rv = {class="RV",...}  // description of RV
512                        //! }
513                        //! \endcode
514                        //!@}
515                        void from_setting(const Setting &set){
516                                epdf::from_setting(set); // reads rv and rvc
517                                UI::get(high,set,"high");
518                                UI::get(low,set,"low");
519                        }
520        };
521
522
523        /*!
524         \brief Normal distributed linear function with linear function of mean value;
525
526         Mean value \f$mu=A*rvc+mu_0\f$.
527        */
528        template<class sq_T>
529        class mlnorm : public mEF
530        {
531                protected:
532                        //! Internal epdf that arise by conditioning on \c rvc
533                        shared_ptr<enorm<sq_T> > iepdf;
534                        mat A;
535                        vec mu_const;
536                        vec& _mu; //cached epdf.mu;
537                public:
538                        //! \name Constructors
539                        //!@{
540                        mlnorm():iepdf(new enorm<sq_T>()), _mu(iepdf->_mu()) { set_ep(iepdf); };
541                        mlnorm (const mat &A, const vec &mu0, const sq_T &R ) :iepdf(new enorm<sq_T>()), _mu(iepdf->_mu())
542                        {
543                                set_ep(iepdf); set_parameters ( A,mu0,R );
544                        }
545
546                        //! Set \c A and \c R
547                        void set_parameters ( const  mat &A, const vec &mu0, const sq_T &R );
548                        //!@}
549                        //! Set value of \c rvc . Result of this operation is stored in \c epdf use function \c _ep to access it.
550                        void condition ( const vec &cond );
551
552                        //!access function
553                        vec& _mu_const() {return mu_const;}
554                        //!access function
555                        mat& _A() {return A;}
556                        //!access function
557                        mat _R() { return iepdf->_R().to_mat(); }
558
559                        template<class sq_M>
560                        friend std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_M> &ml );
561                       
562                        void from_setting(const Setting &set){
563                                mpdf::from_setting(set);
564                                UI::get(A,set,"A");
565                                UI::get(mu_const,set,"const");
566                                mat R0;
567                                UI::get(R0,set,"R");
568                                set_parameters(A,mu_const,R0);
569                        };
570        };
571        UIREGISTER(mlnorm<ldmat>);
572        UIREGISTER(mlnorm<fsqmat>);
573        UIREGISTER(mlnorm<chmat>);
574
575//! Mpdf with general function for mean value
576        template<class sq_T>
577        class mgnorm : public mEF
578        {
579                protected:
580                        //! Internal epdf that arise by conditioning on \c rvc
581                        shared_ptr<enorm<sq_T> > iepdf;
582                        vec &mu;
583                        fnc* g;
584                public:
585                        //!default constructor
586                        mgnorm():iepdf(new enorm<sq_T>()), mu(iepdf->_mu()) { set_ep(iepdf); }
587                        //!set mean function
588                        void set_parameters ( fnc* g0, const sq_T &R0 ) {g=g0; iepdf->set_parameters ( zeros ( g->dimension() ), R0 );}
589                        void condition ( const vec &cond ) {mu=g->eval ( cond );};
590
591
592                        /*! UI for mgnorm
593
594                        The mgnorm is constructed from a structure with fields:
595                        \code
596                        system = {
597                                type = "mgnorm";
598                                // function for mean value evolution
599                                g = {type="fnc"; ... }
600
601                                // variance
602                                R = [1, 0,
603                                         0, 1];
604                                // --OR --
605                                dR = [1, 1];
606
607                                // == OPTIONAL ==
608
609                                // description of y variables
610                                y = {type="rv"; names=["y", "u"];};
611                                // description of u variable
612                                u = {type="rv"; names=[];}
613                        };
614                        \endcode
615
616                        Result if
617                        */
618
619                        void from_setting( const Setting &set ) 
620                        {       
621                                fnc* g = UI::build<fnc>( set, "g" );
622
623                                mat R;
624                                if ( set.exists( "dR" ) )
625                                {
626                                        vec dR;
627                                        UI::get( dR, set, "dR" );
628                                        R=diag(dR);
629                                }
630                                else 
631                                        UI::get( R, set, "R");                                 
632               
633                                set_parameters(g,R);
634                        }
635
636                        /*void mgnorm::to_setting( Setting &set ) const
637                        {       
638                                Transport::to_setting( set );
639
640                                Setting &kilometers_setting = set.add("kilometers", Setting::TypeInt );
641                                kilometers_setting = kilometers;
642
643                                UI::save( passengers, set, "passengers" );
644                        }*/
645
646        };
647
648        UIREGISTER(mgnorm<chmat>);
649
650
651        /*! (Approximate) Student t density with linear function of mean value
652
653        The internal epdf of this class is of the type of a Gaussian (enorm).
654        However, each conditioning is trying to assure the best possible approximation by taking into account the zeta function. See [] for reference.
655
656        Perhaps a moment-matching technique?
657        */
658        class mlstudent : public mlnorm<ldmat>
659        {
660                protected:
661                        ldmat Lambda;
662                        ldmat &_R;
663                        ldmat Re;
664                public:
665                        mlstudent ( ) :mlnorm<ldmat> (),
666                                        Lambda (),      _R ( iepdf->_R() ) {}
667                        void set_parameters ( const mat &A0, const vec &mu0, const ldmat &R0, const ldmat& Lambda0 )
668                        {
669                                it_assert_debug ( A0.rows() ==mu0.length(),"" );
670                                it_assert_debug ( R0.rows() ==A0.rows(),"" );
671
672                                iepdf->set_parameters ( mu0,Lambda ); //
673                                A = A0;
674                                mu_const = mu0;
675                                Re=R0;
676                                Lambda = Lambda0;
677                        }
678                        void condition ( const vec &cond )
679                        {
680                                _mu = A*cond + mu_const;
681                                double zeta;
682                                //ugly hack!
683                                if ( ( cond.length() +1 ) ==Lambda.rows() )
684                                {
685                                        zeta = Lambda.invqform ( concat ( cond, vec_1 ( 1.0 ) ) );
686                                }
687                                else
688                                {
689                                        zeta = Lambda.invqform ( cond );
690                                }
691                                _R = Re;
692                                _R*= ( 1+zeta );// / ( nu ); << nu is in Re!!!!!!
693                        };
694
695        };
696        /*!
697        \brief  Gamma random walk
698
699        Mean value, \f$\mu\f$, of this density is given by \c rvc .
700        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
701        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
702
703        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
704        */
705        class mgamma : public mEF
706        {
707                protected:
708                        //! Internal epdf that arise by conditioning on \c rvc
709                        shared_ptr<egamma> iepdf;
710
711                        //! Constant \f$k\f$
712                        double k;
713
714                        //! cache of iepdf.beta
715                        vec &_beta;
716
717                public:
718                        //! Constructor
719                        mgamma():iepdf(new egamma()), k(0),
720                            _beta(iepdf->_beta()) {
721                            set_ep(iepdf);
722                        }
723
724                        //! Set value of \c k
725                        void set_parameters(double k, const vec &beta0);
726
727                        void condition ( const vec &val ) {_beta=k/val;};
728                        //! Load from structure with elements:
729                        //!  \code
730                        //! { alpha = [...];         // vector of alpha
731                        //!   k = 1.1;               // multiplicative constant k
732                        //!   rv = {class="RV",...}  // description of RV
733                        //!   rvc = {class="RV",...} // description of RV in condition
734                        //! }
735                        //! \endcode
736                        //!@}
737                        void from_setting(const Setting &set){
738                                mpdf::from_setting(set); // reads rv and rvc
739                                vec betatmp; // ugly but necessary
740                                UI::get(betatmp,set,"beta");
741                                set.lookupValue("k",k);
742                                set_parameters(k,betatmp);
743                        }
744        };
745        UIREGISTER(mgamma);
746       
747        /*!
748        \brief  Inverse-Gamma random walk
749
750        Mean value, \f$ \mu \f$, of this density is given by \c rvc .
751        Standard deviation of the random walk is proportional to one \f$ k \f$-th the mean.
752        This is achieved by setting \f$ \alpha=\mu/k^2+2 \f$ and \f$ \beta=\mu(\alpha-1)\f$.
753
754        The standard deviation of the walk is then: \f$ \mu/\sqrt(k)\f$.
755         */
756        class migamma : public mEF
757        {
758                protected:
759                        //! Internal epdf that arise by conditioning on \c rvc
760                        shared_ptr<eigamma> iepdf;
761
762                        //! Constant \f$k\f$
763                        double k;
764
765                        //! cache of iepdf.alpha
766                        vec &_alpha;
767
768                        //! cache of iepdf.beta
769                        vec &_beta;
770
771                public:
772                        //! \name Constructors
773                        //!@{
774                        migamma():iepdf(new eigamma()),
775                            k(0),
776                            _alpha(iepdf->_alpha()),
777                            _beta(iepdf->_beta()) {
778                            set_ep(iepdf);
779                        }
780
781                        migamma(const migamma &m):iepdf(m.iepdf),
782                            k(0),
783                            _alpha(iepdf->_alpha()),
784                            _beta(iepdf->_beta()) {
785                            set_ep(iepdf);
786                        }
787                        //!@}
788
789                        //! Set value of \c k
790                        void set_parameters ( int len, double k0 )
791                        {
792                                k=k0;
793                                iepdf->set_parameters ( ( 1.0/ ( k*k ) +2.0 ) *ones ( len ) /*alpha*/, ones ( len ) /*beta*/ );
794                                dimc = dimension();
795                        };
796                        void condition ( const vec &val )
797                        {
798                                _beta=elem_mult ( val, ( _alpha-1.0 ) );
799                        };
800        };
801
802
803        /*!
804        \brief  Gamma random walk around a fixed point
805
806        Mean value, \f$\mu\f$, of this density is given by a geometric combination of \c rvc and given fixed point, \f$p\f$. \f$l\f$ is the coefficient of the geometric combimation
807        \f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
808
809        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
810        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
811
812        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
813        */
814        class mgamma_fix : public mgamma
815        {
816                protected:
817                        //! parameter l
818                        double l;
819                        //! reference vector
820                        vec refl;
821                public:
822                        //! Constructor
823                        mgamma_fix ( ) : mgamma ( ),refl () {};
824                        //! Set value of \c k
825                        void set_parameters ( double k0 , vec ref0, double l0 )
826                        {
827                                mgamma::set_parameters ( k0, ref0 );
828                                refl=pow ( ref0,1.0-l0 );l=l0;
829                                dimc=dimension();
830                        };
831
832                        void condition ( const vec &val ) {vec mean=elem_mult ( refl,pow ( val,l ) ); _beta=k/mean;};
833        };
834
835
836        /*!
837        \brief  Inverse-Gamma random walk around a fixed point
838
839        Mean value, \f$\mu\f$, of this density is given by a geometric combination of \c rvc and given fixed point, \f$p\f$. \f$l\f$ is the coefficient of the geometric combimation
840        \f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
841
842        ==== Check == vv =
843        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
844        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
845
846        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
847         */
848        class migamma_ref : public migamma
849        {
850                protected:
851                        //! parameter l
852                        double l;
853                        //! reference vector
854                        vec refl;
855                public:
856                        //! Constructor
857                        migamma_ref ( ) : migamma (),refl ( ) {};
858                        //! Set value of \c k
859                        void set_parameters ( double k0 , vec ref0, double l0 )
860                        {
861                                migamma::set_parameters ( ref0.length(), k0 );
862                                refl=pow ( ref0,1.0-l0 );
863                                l=l0;
864                                dimc = dimension();
865                        };
866
867                        void condition ( const vec &val )
868                        {
869                                vec mean=elem_mult ( refl,pow ( val,l ) );
870                                migamma::condition ( mean );
871                        };
872
873                        /*! UI for migamma_ref
874
875                        The migamma_ref is constructed from a structure with fields:
876                        \code
877                        system = {
878                                type = "migamma_ref";
879                                ref = [1e-5; 1e-5; 1e-2 1e-3];            // reference vector
880                                l = 0.999;                                // constant l
881                                k = 0.1;                                  // constant k
882                               
883                                // == OPTIONAL ==
884                                // description of y variables
885                                y = {type="rv"; names=["y", "u"];};
886                                // description of u variable
887                                u = {type="rv"; names=[];}
888                        };
889                        \endcode
890
891                        Result if
892                         */
893                        void from_setting( const Setting &set );
894
895                        // TODO dodelat void to_setting( Setting &set ) const;
896        };
897
898
899        UIREGISTER(migamma_ref);
900
901        /*! Log-Normal probability density
902         only allow diagonal covariances!
903
904        Density of the form \f$ \log(x)\sim \mathcal{N}(\mu,\sigma^2) \f$ , i.e.
905        \f[
906        x \sim \frac{1}{x\sigma\sqrt{2\pi}}\exp{-\frac{1}{2\sigma^2}(\log(x)-\mu)}
907        \f]
908
909        */
910        class elognorm: public enorm<ldmat>
911        {
912                public:
913                        vec sample() const {return exp ( enorm<ldmat>::sample() );};
914                        vec mean() const {vec var=enorm<ldmat>::variance();return exp ( mu - 0.5*var );};
915
916        };
917
918        /*!
919        \brief  Log-Normal random walk
920
921        Mean value, \f$\mu\f$, is...
922
923        ==== Check == vv =
924        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
925        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
926
927        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
928         */
929        class mlognorm : public mpdf
930        {
931                protected:
932                        shared_ptr<elognorm> eno;
933
934                        //! parameter 1/2*sigma^2
935                        double sig2;
936
937                        //! access
938                        vec &mu;
939                public:
940                        //! Constructor
941                        mlognorm():eno(new elognorm()),
942                            sig2(0),
943                            mu(eno->_mu()) {
944                            set_ep(eno);
945                        }
946
947                        //! Set value of \c k
948                        void set_parameters ( int size, double k )
949                        {
950                                sig2 = 0.5*log ( k*k+1 );
951                                eno->set_parameters ( zeros ( size ),2*sig2*eye ( size ) );
952
953                                dimc = size;
954                        };
955
956                        void condition ( const vec &val )
957                        {
958                                mu=log ( val )-sig2;//elem_mult ( refl,pow ( val,l ) );
959                        };
960
961                        /*! UI for mlognorm
962
963                        The mlognorm is constructed from a structure with fields:
964                        \code
965                        system = {
966                                type = "mlognorm";
967                                k = 0.1;                                  // constant k
968                                mu0 = [1., 1.];
969                               
970                                // == OPTIONAL ==
971                                // description of y variables
972                                y = {type="rv"; names=["y", "u"];};
973                                // description of u variable
974                                u = {type="rv"; names=[];}
975                        };
976                        \endcode
977
978                         */
979                        void from_setting( const Setting &set );
980
981                        // TODO dodelat void to_setting( Setting &set ) const;
982
983        };
984       
985        UIREGISTER(mlognorm);
986
987        /*! inverse Wishart density defined on Choleski decomposition
988
989        */
990        class eWishartCh : public epdf
991        {
992                protected:
993                        //! Upper-Triagle of Choleski decomposition of \f$ \Psi \f$
994                        chmat Y;
995                        //! dimension of matrix  \f$ \Psi \f$
996                        int p;
997                        //! degrees of freedom  \f$ \nu \f$
998                        double delta;
999                public:
1000                        void set_parameters ( const mat &Y0, const double delta0 ) {Y=chmat ( Y0 );delta=delta0; p=Y.rows(); dim = p*p; }
1001                        mat sample_mat() const
1002                        {
1003                                mat X=zeros ( p,p );
1004
1005                                //sample diagonal
1006                                for ( int i=0;i<p;i++ )
1007                                {
1008                                        GamRNG.setup ( 0.5* ( delta-i ) , 0.5 ); // no +1 !! index if from 0
1009#pragma omp critical
1010                                        X ( i,i ) =sqrt ( GamRNG() );
1011                                }
1012                                //do the rest
1013                                for ( int i=0;i<p;i++ )
1014                                {
1015                                        for ( int j=i+1;j<p;j++ )
1016                                        {
1017#pragma omp critical
1018                                                X ( i,j ) =NorRNG.sample();
1019                                        }
1020                                }
1021                                return X*Y._Ch();// return upper triangular part of the decomposition
1022                        }
1023                        vec sample () const
1024                        {
1025                                return vec ( sample_mat()._data(),p*p );
1026                        }
1027                        //! fast access function y0 will be copied into Y.Ch.
1028                        void setY ( const mat &Ch0 ) {copy_vector ( dim,Ch0._data(), Y._Ch()._data() );}
1029                        //! fast access function y0 will be copied into Y.Ch.
1030                        void _setY ( const vec &ch0 ) {copy_vector ( dim, ch0._data(), Y._Ch()._data() ); }
1031                        //! access function
1032                        const chmat& getY()const {return Y;}
1033        };
1034
1035        class eiWishartCh: public epdf
1036        {
1037                protected:
1038                        eWishartCh W;
1039                        int p;
1040                        double delta;
1041                public:
1042                        void set_parameters ( const mat &Y0, const double delta0) {
1043                                delta = delta0;
1044                                W.set_parameters ( inv ( Y0 ),delta0 ); 
1045                                dim = W.dimension(); p=Y0.rows();
1046                        }
1047                        vec sample() const {mat iCh; iCh=inv ( W.sample_mat() ); return vec ( iCh._data(),dim );}
1048                        void _setY ( const vec &y0 )
1049                        {
1050                                mat Ch ( p,p );
1051                                mat iCh ( p,p );
1052                                copy_vector ( dim, y0._data(), Ch._data() );
1053                               
1054                                iCh=inv ( Ch );
1055                                W.setY ( iCh );
1056                        }
1057                        virtual double evallog ( const vec &val ) const {
1058                                chmat X(p);
1059                                const chmat& Y=W.getY();
1060                                 
1061                                copy_vector(p*p,val._data(),X._Ch()._data());
1062                                chmat iX(p);X.inv(iX);
1063                                // compute 
1064//                              \frac{ |\Psi|^{m/2}|X|^{-(m+p+1)/2}e^{-tr(\Psi X^{-1})/2} }{ 2^{mp/2}\Gamma_p(m/2)},
1065                                mat M=Y.to_mat()*iX.to_mat();
1066                               
1067                                double log1 = 0.5*p*(2*Y.logdet())-0.5*(delta+p+1)*(2*X.logdet())-0.5*trace(M); 
1068                                //Fixme! Multivariate gamma omitted!! it is ok for sampling, but not otherwise!!
1069                               
1070/*                              if (0) {
1071                                        mat XX=X.to_mat();
1072                                        mat YY=Y.to_mat();
1073                                       
1074                                        double log2 = 0.5*p*log(det(YY))-0.5*(delta+p+1)*log(det(XX))-0.5*trace(YY*inv(XX));
1075                                        cout << log1 << "," << log2 << endl;
1076                                }*/
1077                                return log1;                           
1078                        };
1079                       
1080        };
1081
1082        class rwiWishartCh : public mpdf
1083        {
1084                protected:
1085                        shared_ptr<eiWishartCh> eiW;
1086                        //!square root of \f$ \nu-p-1 \f$ - needed for computation of \f$ \Psi \f$ from conditions
1087                        double sqd;
1088                        //reference point for diagonal
1089                        vec refl;
1090                        double l;
1091                        int p;
1092
1093                public:
1094                        rwiWishartCh():eiW(new eiWishartCh()),
1095                            sqd(0), l(0), p(0) {
1096                            set_ep(eiW);
1097                        }
1098
1099                        void set_parameters ( int p0, double k, vec ref0, double l0  )
1100                        {
1101                                p=p0;
1102                                double delta = 2/(k*k)+p+3;
1103                                sqd=sqrt ( delta-p-1 );
1104                                l=l0;
1105                                refl=pow(ref0,1-l);
1106                               
1107                                eiW->set_parameters ( eye ( p ),delta );
1108                                dimc=eiW->dimension();
1109                        }
1110                        void condition ( const vec &c ) {
1111                                vec z=c;
1112                                int ri=0;
1113                                for(int i=0;i<p*p;i+=(p+1)){//trace diagonal element
1114                                        z(i) = pow(z(i),l)*refl(ri);
1115                                        ri++;
1116                                }
1117
1118                                eiW->_setY ( sqd*z );
1119                        }
1120        };
1121
1122//! Switch between various resampling methods.
1123        enum RESAMPLING_METHOD { MULTINOMIAL = 0, STRATIFIED = 1, SYSTEMATIC = 3 };
1124        /*!
1125        \brief Weighted empirical density
1126
1127        Used e.g. in particle filters.
1128        */
1129        class eEmp: public epdf
1130        {
1131                protected :
1132                        //! Number of particles
1133                        int n;
1134                        //! Sample weights \f$w\f$
1135                        vec w;
1136                        //! Samples \f$x^{(i)}, i=1..n\f$
1137                        Array<vec> samples;
1138                public:
1139                        //! \name Constructors
1140                        //!@{
1141                        eEmp ( ) :epdf ( ),w ( ),samples ( ) {};
1142                        //! copy constructor
1143                        eEmp ( const eEmp &e ) : epdf ( e ), w ( e.w ), samples ( e.samples ) {};
1144                        //!@}
1145
1146                        //! Set samples and weights
1147                        void set_statistics ( const vec &w0, const epdf* pdf0 );
1148                        //! Set samples and weights
1149                        void set_statistics ( const epdf* pdf0 , int n ) {set_statistics ( ones ( n ) /n,pdf0 );};
1150                        //! Set sample
1151                        void set_samples ( const epdf* pdf0 );
1152                        //! Set sample
1153                        void set_parameters ( int n0, bool copy=true ) {n=n0; w.set_size ( n0,copy );samples.set_size ( n0,copy );};
1154                        //! Potentially dangerous, use with care.
1155                        vec& _w()  {return w;};
1156                        //! Potentially dangerous, use with care.
1157                        const vec& _w() const {return w;};
1158                        //! access function
1159                        Array<vec>& _samples() {return samples;};
1160                        //! access function
1161                        const Array<vec>& _samples() const {return samples;};
1162                        //! Function performs resampling, i.e. removal of low-weight samples and duplication of high-weight samples such that the new samples represent the same density.
1163                        ivec resample ( RESAMPLING_METHOD method=SYSTEMATIC );
1164                        //! inherited operation : NOT implemneted
1165                        vec sample() const {it_error ( "Not implemented" );return 0;}
1166                        //! inherited operation : NOT implemneted
1167                        double evallog ( const vec &val ) const {it_error ( "Not implemented" );return 0.0;}
1168                        vec mean() const
1169                        {
1170                                vec pom=zeros ( dim );
1171                                for ( int i=0;i<n;i++ ) {pom+=samples ( i ) *w ( i );}
1172                                return pom;
1173                        }
1174                        vec variance() const
1175                        {
1176                                vec pom=zeros ( dim );
1177                                for ( int i=0;i<n;i++ ) {pom+=pow ( samples ( i ),2 ) *w ( i );}
1178                                return pom-pow ( mean(),2 );
1179                        }
1180                        //! For this class, qbounds are minimum and maximum value of the population!
1181                        void qbounds ( vec &lb, vec &ub, double perc=0.95 ) const
1182                        {
1183                                // lb in inf so than it will be pushed below;
1184                                lb.set_size ( dim );
1185                                ub.set_size ( dim );
1186                                lb = std::numeric_limits<double>::infinity();
1187                                ub = -std::numeric_limits<double>::infinity();
1188                                int j;
1189                                for ( int i=0;i<n;i++ )
1190                                {
1191                                        for ( j=0;j<dim; j++ )
1192                                        {
1193                                                if ( samples ( i ) ( j ) <lb ( j ) ) {lb ( j ) =samples ( i ) ( j );}
1194                                                if ( samples ( i ) ( j ) >ub ( j ) ) {ub ( j ) =samples ( i ) ( j );}
1195                                        }
1196                                }
1197                        }
1198        };
1199
1200
1201////////////////////////
1202
1203        template<class sq_T>
1204        void enorm<sq_T>::set_parameters ( const vec &mu0, const sq_T &R0 )
1205        {
1206//Fixme test dimensions of mu0 and R0;
1207                mu = mu0;
1208                R = R0;
1209                validate();
1210        };
1211
1212        template<class sq_T>
1213        void enorm<sq_T>::from_setting(const Setting &set){
1214                epdf::from_setting(set); //reads rv
1215               
1216                UI::get(mu,set,"mu");
1217                mat Rtmp;// necessary for conversion
1218                UI::get(Rtmp,set,"R");
1219                R=Rtmp; // conversion
1220                validate();
1221        }
1222
1223        template<class sq_T>
1224        void enorm<sq_T>::dupdate ( mat &v, double nu )
1225        {
1226                //
1227        };
1228
1229// template<class sq_T>
1230// void enorm<sq_T>::tupdate ( double phi, mat &vbar, double nubar ) {
1231//      //
1232// };
1233
1234        template<class sq_T>
1235        vec enorm<sq_T>::sample() const
1236        {
1237                vec x ( dim );
1238#pragma omp critical
1239                NorRNG.sample_vector ( dim,x );
1240                vec smp = R.sqrt_mult ( x );
1241
1242                smp += mu;
1243                return smp;
1244        };
1245
1246// template<class sq_T>
1247// double enorm<sq_T>::eval ( const vec &val ) const {
1248//      double pdfl,e;
1249//      pdfl = evallog ( val );
1250//      e = exp ( pdfl );
1251//      return e;
1252// };
1253
1254        template<class sq_T>
1255        double enorm<sq_T>::evallog_nn ( const vec &val ) const
1256        {
1257                // 1.83787706640935 = log(2pi)
1258                double tmp=-0.5* ( R.invqform ( mu-val ) );// - lognc();
1259                return  tmp;
1260        };
1261
1262        template<class sq_T>
1263        inline double enorm<sq_T>::lognc () const
1264        {
1265                // 1.83787706640935 = log(2pi)
1266                double tmp=0.5* ( R.cols() * 1.83787706640935 +R.logdet() );
1267                return tmp;
1268        };
1269
1270        template<class sq_T>
1271        void mlnorm<sq_T>::set_parameters ( const mat &A0, const vec &mu0, const sq_T &R0 )
1272        {
1273                it_assert_debug ( A0.rows() ==mu0.length(),"" );
1274                it_assert_debug ( A0.rows() ==R0.rows(),"" );
1275
1276                iepdf->set_parameters(zeros(A0.rows()), R0);
1277                A = A0;
1278                mu_const = mu0;
1279                dimc = A0.cols();
1280        }
1281
1282// template<class sq_T>
1283// vec mlnorm<sq_T>::samplecond (const  vec &cond, double &lik ) {
1284//      this->condition ( cond );
1285//      vec smp = epdf.sample();
1286//      lik = epdf.eval ( smp );
1287//      return smp;
1288// }
1289
1290// template<class sq_T>
1291// mat mlnorm<sq_T>::samplecond (const vec &cond, vec &lik, int n ) {
1292//      int i;
1293//      int dim = rv.count();
1294//      mat Smp ( dim,n );
1295//      vec smp ( dim );
1296//      this->condition ( cond );
1297//
1298//      for ( i=0; i<n; i++ ) {
1299//              smp = epdf.sample();
1300//              lik ( i ) = epdf.eval ( smp );
1301//              Smp.set_col ( i ,smp );
1302//      }
1303//
1304//      return Smp;
1305// }
1306
1307        template<class sq_T>
1308        void mlnorm<sq_T>::condition ( const vec &cond )
1309        {
1310                _mu = A*cond + mu_const;
1311//R is already assigned;
1312        }
1313
1314        template<class sq_T>
1315        enorm<sq_T>* enorm<sq_T>::marginal ( const RV &rvn ) const
1316        {
1317                it_assert_debug ( isnamed(), "rv description is not assigned" );
1318                ivec irvn = rvn.dataind ( rv );
1319
1320                sq_T Rn ( R,irvn ); //select rows and columns of R
1321
1322                enorm<sq_T>* tmp = new enorm<sq_T>;
1323                tmp->set_rv ( rvn );
1324                tmp->set_parameters ( mu ( irvn ), Rn );
1325                return tmp;
1326        }
1327
1328        template<class sq_T>
1329        mpdf* enorm<sq_T>::condition ( const RV &rvn ) const
1330        {
1331
1332                it_assert_debug ( isnamed(),"rvs are not assigned" );
1333
1334                RV rvc = rv.subt ( rvn );
1335                it_assert_debug ( ( rvc._dsize() +rvn._dsize() ==rv._dsize() ),"wrong rvn" );
1336                //Permutation vector of the new R
1337                ivec irvn = rvn.dataind ( rv );
1338                ivec irvc = rvc.dataind ( rv );
1339                ivec perm=concat ( irvn , irvc );
1340                sq_T Rn ( R,perm );
1341
1342                //fixme - could this be done in general for all sq_T?
1343                mat S=Rn.to_mat();
1344                //fixme
1345                int n=rvn._dsize()-1;
1346                int end=R.rows()-1;
1347                mat S11 = S.get ( 0,n, 0, n );
1348                mat S12 = S.get ( 0, n , rvn._dsize(), end );
1349                mat S22 = S.get ( rvn._dsize(), end, rvn._dsize(), end );
1350
1351                vec mu1 = mu ( irvn );
1352                vec mu2 = mu ( irvc );
1353                mat A=S12*inv ( S22 );
1354                sq_T R_n ( S11 - A *S12.T() );
1355
1356                mlnorm<sq_T>* tmp=new mlnorm<sq_T> ( );
1357                tmp->set_rv ( rvn ); tmp->set_rvc ( rvc );
1358                tmp->set_parameters ( A,mu1-A*mu2,R_n );
1359                return tmp;
1360        }
1361
1362///////////
1363
1364        template<class sq_T>
1365        std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_T> &ml )
1366        {
1367                os << "A:"<< ml.A<<endl;
1368                os << "mu:"<< ml.mu_const<<endl;
1369                os << "R:" << ml.iepdf->_R().to_mat() <<endl;
1370                return os;
1371        };
1372
1373}
1374#endif //EF_H
Note: See TracBrowser for help on using the browser.