root/library/bdm/estim/kalman.cpp

Revision 1467, 17.1 kB (checked in by smidl, 12 years ago)

vykreslovani v kalmanovi

  • Property svn:eol-style set to native
Line 
1#include "kalman.h"
2
3namespace bdm {
4
5using std::endl;
6
7
8
9void KalmanFull::bayes ( const vec &yt, const vec &cond ) {
10    bdm_assert_debug ( yt.length() == ( dimy ), "KalmanFull::bayes wrong size of dt, " +
11                       num2str(yt.length()) + ", expected size is " + num2str(dimy) );
12    bdm_assert_debug ( cond.length() == ( dimc ), "KalmanFull::bayes wrong size of cond, " +
13                       num2str(cond.length()) + ", expected size is " + num2str(dimc) );
14
15    const vec &u = cond; // in this case cond=ut
16    const vec &y = yt;
17
18    vec& mu = est._mu();
19    mat &P = est._R();
20    mat& _Ry = fy._R();
21    vec& yp = fy._mu();
22    //Time update
23    mu = A * mu + B * u;
24    P  = A * P * A.transpose() + ( mat ) Q;
25
26    //Data update
27    _Ry = C * P * C.transpose() + ( mat ) R;
28    _K = P * C.transpose() * inv ( _Ry );
29    P -= _K * C * P; // P = P -KCP;
30    yp = C * mu + D * u;
31    mu += _K * ( y - yp );
32
33    if ( evalll ) {
34        ll = fy.evallog ( y );
35    }
36};
37
38
39
40/////////////////////////////// EKFS
41EKFfull::EKFfull ( ) : KalmanFull () {};
42
43void EKFfull::set_parameters ( const shared_ptr<diffbifn> &pfxu0, const shared_ptr<diffbifn> &phxu0, const mat Q0, const mat R0 ) {
44    pfxu = pfxu0;
45    phxu = phxu0;
46
47    set_dim ( pfxu0->_dimx() );
48    dimy = phxu0->dimension();
49    dimc = pfxu0->_dimu();
50    est.set_parameters ( zeros ( dimension() ), eye ( dimension() ) );
51
52    A.set_size ( dimension(), dimension() );
53    C.set_size ( dimy, dimension() );
54    //initialize matrices A C, later, these will be only updated!
55    pfxu->dfdx_cond ( est._mu(), zeros ( dimc ), A, true );
56    B.clear();
57    phxu->dfdx_cond ( est._mu(), zeros ( dimc ), C, true );
58    D.clear();
59
60    R = R0;
61    Q = Q0;
62}
63
64void EKFfull::bayes ( const vec &yt, const vec &cond ) {
65    bdm_assert_debug ( yt.length() == ( dimy ), "EKFull::bayes wrong size of dt" );
66    bdm_assert_debug ( cond.length() == ( dimc ), "EKFull::bayes wrong size of dt" );
67
68    const vec &u = cond;
69    const vec &y = yt; //lazy to change it
70    vec &mu = est._mu();
71    mat &P = est._R();
72    mat& _Ry = fy._R();
73    vec& yp = fy._mu();
74
75    pfxu->dfdx_cond ( mu, zeros ( dimc ), A, true );
76    phxu->dfdx_cond ( mu, zeros ( dimc ), C, true );
77
78    //Time update
79    mu = pfxu->eval ( mu, u );// A*mu + B*u;
80    P  = A * P * A.transpose() + ( mat ) Q;
81
82    //Data update
83    _Ry = C * P * C.transpose() + ( mat ) R;
84    _K = P * C.transpose() * inv ( _Ry );
85    P -= _K * C * P; // P = P -KCP;
86    yp = phxu->eval ( mu, u );
87    mu += _K * ( y - yp );
88
89    if ( BM::evalll ) {
90        ll = fy.evallog ( y );
91    }
92};
93
94
95
96void KalmanCh::set_parameters ( const mat &A0, const mat &B0, const mat &C0, const mat &D0, const chmat &Q0, const chmat &R0 ) {
97
98    ( ( StateSpace<chmat>* ) this )->set_parameters ( A0, B0, C0, D0, Q0, R0 );
99
100    _K = zeros ( dimension(), dimy );
101}
102
103void KalmanCh::initialize() {
104    preA = zeros ( dimy + dimension() + dimension(), dimy + dimension() );
105//      preA.clear();
106    preA.set_submatrix ( 0, 0, R._Ch() );
107    preA.set_submatrix ( dimy + dimension(), dimy, Q._Ch() );
108}
109
110void KalmanCh::bayes ( const vec &yt, const vec &cond ) {
111    bdm_assert_debug ( yt.length() == dimy, "yt mismatch" );
112    bdm_assert_debug ( cond.length() == dimc, "yt mismatch" );
113
114    const vec &u = cond;
115    const vec &y = yt;
116    vec pom ( dimy );
117
118    chmat &_P = est._R();
119    vec &_mu = est._mu();
120    mat _K ( dimension(), dimy );
121    chmat &_Ry = fy._R();
122    vec &_yp = fy._mu();
123    //TODO get rid of Q in qr()!
124    //  mat Q;
125
126    //R and Q are already set in set_parameters()
127    preA.set_submatrix ( dimy, 0, ( _P._Ch() ) *C.T() );
128    //Fixme can be more efficient if .T() is not used
129    preA.set_submatrix ( dimy, dimy, ( _P._Ch() ) *A.T() );
130
131    if ( !qr ( preA, postA ) ) {
132        bdm_warning ( "QR in KalmanCh unstable!" );
133    }
134
135    ( _Ry._Ch() ) = postA ( 0, dimy - 1, 0, dimy - 1 );
136    _K = inv ( A ) * ( postA ( 0, dimy - 1 , dimy, dimy + dimension() - 1 ) ).T();
137    ( _P._Ch() ) = postA ( dimy, dimy + dimension() - 1, dimy, dimy + dimension() - 1 );
138
139    _mu = A * ( _mu ) + B * u;
140    _yp = C * _mu - D * u;
141
142    backward_substitution ( _Ry._Ch(), ( y - _yp ), pom );
143    _mu += ( _K ) * pom;
144
145    /*          cout << "P:" <<_P.to_mat() <<endl;
146                cout << "Ry:" <<_Ry.to_mat() <<endl;
147                cout << "_K:" <<_K <<endl;*/
148
149    if ( evalll == true ) { //likelihood of observation y
150        ll = fy.evallog ( y );
151    }
152}
153
154void StateCanonical::connect_mlnorm ( const mlnorm<fsqmat> &ml ) {
155    //get ids of yrv
156    const RV &yrv = ml._rv();
157    //need to determine u_t - it is all in _rvc that is not in ml._rv()
158    RV rgr0 = ml._rvc().remove_time();
159    RV urv = rgr0.subt ( yrv );
160
161    //We can do only 1d now... :(
162    bdm_assert ( yrv._dsize() == 1, "Only for SISO so far..." );
163
164    // create names for
165    RV xrv; //empty
166    RV Crv; //empty
167    int td = ml._rvc().mint();
168    // assuming strictly proper function!!!
169    for ( int t = -1; t >= td; t-- ) {
170        xrv.add ( yrv.copy_t ( t ) );
171        Crv.add ( urv.copy_t ( t ) );
172    }
173
174    // get mapp
175    th2A.set_connection ( xrv, ml._rvc() );
176    th2C.set_connection ( Crv, ml._rvc() );
177    th2D.set_connection ( urv, ml._rvc() );
178
179    //set matrix sizes
180    this->A = zeros ( xrv._dsize(), xrv._dsize() );
181    for ( int j = 1; j < xrv._dsize(); j++ ) {
182        A ( j, j - 1 ) = 1.0;    // off diagonal
183    }
184    this->B = zeros ( xrv._dsize(), 1 );
185    this->B ( 0 ) = 1.0;
186    this->C = zeros ( 1, xrv._dsize() );
187    this->D = zeros ( 1, urv._dsize() );
188    this->Q = zeros ( xrv._dsize(), xrv._dsize() );
189    // R is set by update
190
191    //set cache
192    this->A1row = zeros ( xrv._dsize() );
193    this->C1row = zeros ( xrv._dsize() );
194    this->D1row = zeros ( urv._dsize() );
195
196    update_from ( ml );
197    validate();
198};
199
200void StateCanonical::update_from ( const mlnorm<fsqmat> &ml ) {
201
202    vec theta = ml._A().get_row ( 0 ); // this
203
204    th2A.filldown ( theta, A1row );
205    th2C.filldown ( theta, C1row );
206    th2D.filldown ( theta, D1row );
207
208    R = ml._R();
209
210    A.set_row ( 0, A1row );
211    C.set_row ( 0, C1row + D1row ( 0 ) *A1row );
212    D.set_row ( 0, D1row );
213}
214
215void StateFromARX::connect_mlnorm ( const mlnorm<chmat> &ml, RV &xrv, RV &urv ) {
216    //get ids of yrv
217    const RV &yrv = ml._rv();
218    //need to determine u_t - it is all in _rvc that is not in ml._rv()
219    const RV &rgr = ml._rvc();
220    RV rgr0 = rgr.remove_time();
221    urv = rgr0.subt ( yrv );
222
223    // create names for state variables
224    xrv = yrv;
225
226    int y_multiplicity = -rgr.mint ( yrv );
227    int y_block_size = yrv.length() * ( y_multiplicity ); // current yt + all delayed yts
228    for ( int m = 0; m < y_multiplicity - 1; m++ ) { // ========= -1 is important see arx2statespace_notes
229        xrv.add ( yrv.copy_t ( -m - 1 ) ); //add delayed yt
230    }
231    //! temporary RV for connection to ml.rvc, since notation of xrv and ml.rvc does not match
232    RV xrv_ml = xrv.copy_t ( -1 );
233
234    // add regressors
235    ivec u_block_sizes ( urv.length() ); // no_blocks = yt + unique rgr
236    for ( int r = 0; r < urv.length(); r++ ) {
237        RV R = urv.subselect ( vec_1 ( r ) ); //non-delayed element of rgr
238        int r_size = urv.size ( r );
239        int r_multiplicity = -rgr.mint ( R );
240        u_block_sizes ( r ) = r_size * r_multiplicity ;
241        for ( int m = 0; m < r_multiplicity; m++ ) {
242            xrv.add ( R.copy_t ( -m - 1 ) ); //add delayed yt
243            xrv_ml.add ( R.copy_t ( -m - 1 ) ); //add delayed yt
244        }
245    }
246    // add constant
247    if ( any ( ml._mu_const() != 0.0 ) ) {
248        have_constant = true;
249        xrv.add ( RV ( "bdm_reserved_constant_one", 1 ) );
250    } else {
251        have_constant = false;
252    }
253
254
255    // get mapp
256    th2A.set_connection ( xrv_ml, ml._rvc() );
257    th2B.set_connection ( urv, ml._rvc() );
258
259    //set matrix sizes
260    this->A = zeros ( xrv._dsize(), xrv._dsize() );
261    //create y block
262    diagonal_part ( this->A, yrv._dsize(), 0, y_block_size - yrv._dsize() );
263
264    this->B = zeros ( xrv._dsize(), urv._dsize() );
265    //add diagonals for rgr
266    int active_x = y_block_size;
267    int active_Bcol = 0;
268    for ( int r = 0; r < urv.length(); r++ ) {
269        if (u_block_sizes(r)>0) {
270            diagonal_part ( this->A, active_x + urv.size ( r ), active_x, u_block_sizes ( r ) - urv.size ( r ) );
271            this->B.set_submatrix ( active_x, active_Bcol, eye ( urv.size ( r ) ) );
272            active_Bcol+=u_block_sizes(r);
273        }
274        active_x += u_block_sizes ( r );
275    }
276    this->C = zeros ( yrv._dsize(), xrv._dsize() );
277    this->C.set_submatrix ( 0, 0, eye ( yrv._dsize() ) );
278    this->D = zeros ( yrv._dsize(), urv._dsize() );
279    this->R.setCh ( zeros ( yrv._dsize(), yrv._dsize() ) );
280    this->Q.setCh ( zeros ( xrv._dsize(), xrv._dsize() ) );
281    // Q is set by update
282
283    update_from ( ml );
284    validate();
285}
286
287void StateFromARX::update_from ( const mlnorm<chmat> &ml ) {
288    vec Arow = zeros ( A.cols() );
289    vec Brow = zeros ( B.cols() );
290    //  ROW- WISE EVALUATION =====
291    for ( int i = 0; i < ml._rv()._dsize(); i++ ) {
292
293        vec theta = ml._A().get_row ( i );
294
295        th2A.filldown ( theta, Arow );
296        if ( have_constant ) {
297            // constant is always at the end no need for datalink
298            Arow ( A.cols() - 1 ) = ml._mu_const() ( i );
299        }
300        this->A.set_row ( i, Arow );
301
302        th2B.filldown ( theta, Brow );
303        this->B.set_row ( i, Brow );
304    }
305    this->Q._Ch().set_submatrix ( 0, 0, ml.__R()._Ch() );
306
307}
308
309
310void EKFCh::set_parameters ( const shared_ptr<diffbifn> &pfxu0, const shared_ptr<diffbifn> &phxu0, const chmat Q0, const chmat R0 ) {
311    pfxu = pfxu0;
312    phxu = phxu0;
313
314    set_dim ( pfxu0->_dimx() );
315    dimy = phxu0->dimension();
316    dimc = pfxu0->_dimu();
317
318    vec &_mu = est._mu();
319    // if mu is not set, set it to zeros, just for constant terms of A and C
320    if ( _mu.length() != dimension() ) _mu = zeros ( dimension() );
321    A = zeros ( dimension(), dimension() );
322    C = zeros ( dimy, dimension() );
323    preA = zeros ( dimy + dimension() + dimension(), dimy + dimension() );
324
325    //initialize matrices A C, later, these will be only updated!
326    pfxu->dfdx_cond ( _mu, zeros ( dimc ), A, true );
327//      pfxu->dfdu_cond ( *_mu,zeros ( dimu ),B,true );
328    B.clear();
329    phxu->dfdx_cond ( _mu, zeros ( dimc ), C, true );
330//      phxu->dfdu_cond ( *_mu,zeros ( dimu ),D,true );
331    D.clear();
332
333    R = R0;
334    Q = Q0;
335
336    // Cholesky special!
337    preA.clear();
338    preA.set_submatrix ( 0, 0, R._Ch() );
339    preA.set_submatrix ( dimy + dimension(), dimy, Q._Ch() );
340}
341
342
343void EKFCh::bayes ( const vec &yt, const vec &cond ) {
344
345    vec pom ( dimy );
346    const vec &u = cond;
347    const vec &y = yt;
348    vec &_mu = est._mu();
349    chmat &_P = est._R();
350    chmat &_Ry = fy._R();
351    vec &_yp = fy._mu();
352
353    pfxu->dfdx_cond ( _mu, u, A, false ); //update A by a derivative of fx
354    phxu->dfdx_cond ( _mu, u, C, false ); //update A by a derivative of fx
355
356    //R and Q are already set in set_parameters()
357    preA.set_submatrix ( dimy, 0, ( _P._Ch() ) *C.T() );
358    //Fixme can be more efficient if .T() is not used
359    preA.set_submatrix ( dimy, dimy, ( _P._Ch() ) *A.T() );
360
361//      mat Sttm = _P->to_mat();
362//      cout << preA <<endl;
363//      cout << "_mu:" << _mu <<endl;
364
365    if ( !qr ( preA, postA ) ) {
366        bdm_warning ( "QR in EKFCh unstable!" );
367    }
368
369
370    ( _Ry._Ch() ) = postA ( 0, dimy - 1, 0, dimy - 1 );
371    _K = inv ( A ) * ( postA ( 0, dimy - 1 , dimy, dimy + dimension() - 1 ) ).T();
372    ( _P._Ch() ) = postA ( dimy, dimy + dimension() - 1, dimy, dimy + dimension() - 1 );
373
374//      mat iRY = inv(_Ry->to_mat());
375//      mat Stt = Sttm - Sttm * C.T() * iRY * C * Sttm;
376//      mat _K2 = Stt*C.T()*inv(R.to_mat());
377
378    // prediction
379    _mu = pfxu->eval ( _mu , u );
380    _yp = phxu->eval ( _mu, u );
381
382    /*  vec mu2 = *_mu + ( _K2 ) * ( y-*_yp );*/
383
384    //correction //= initial value is already prediction!
385    backward_substitution ( _Ry._Ch(), ( y - _yp ), pom );
386    _mu += ( _K ) * pom ;
387
388    /*  _K = (_P->to_mat())*C.transpose() * ( _iRy->to_mat() );
389        *_mu = pfxu->eval ( *_mu ,u ) + ( _K )* ( y-*_yp );*/
390
391//              cout << "P:" <<_P.to_mat() <<endl;
392//              cout << "Ry:" <<_Ry.to_mat() <<endl;
393//      cout << "_mu:" <<_mu <<endl;
394//      cout << "dt:" <<dt <<endl;
395
396    if ( evalll == true ) { //likelihood of observation y
397        ll = fy.evallog ( y );
398    }
399}
400
401void EKFCh::from_setting ( const Setting &set ) {
402    BM::from_setting ( set );
403    shared_ptr<diffbifn> IM = UI::build<diffbifn> ( set, "IM", UI::compulsory );
404    shared_ptr<diffbifn> OM = UI::build<diffbifn> ( set, "OM", UI::compulsory );
405
406    //statistics
407    int dim = IM->dimension();
408    vec mu0;
409    if ( !UI::get ( mu0, set, "mu0" ) )
410        mu0 = zeros ( dim );
411
412    mat P0;
413    vec dP0;
414    if ( UI::get ( dP0, set, "dP0" ) )
415        P0 = diag ( dP0 );
416    else if ( !UI::get ( P0, set, "P0" ) )
417        P0 = eye ( dim );
418
419    set_statistics ( mu0, P0 );
420
421    //parameters
422    vec dQ, dR;
423    UI::get ( dQ, set, "dQ", UI::compulsory );
424    UI::get ( dR, set, "dR", UI::compulsory );
425    set_parameters ( IM, OM, diag ( dQ ), diag ( dR ) );
426        UI::get ( log_level, set, "log_level", UI::optional );
427}
428
429void MultiModel::from_setting ( const Setting &set ) {
430    Array<EKFCh*> A;
431    UI::get ( A, set, "models", UI::compulsory );
432
433    set_parameters ( A );
434    set_yrv ( A ( 0 )->_yrv() );
435    //set_rv(A(0)->_rv());
436}
437
438void EKF_UD::set_parameters ( const shared_ptr<diffbifn> &pfxu0, const shared_ptr<diffbifn> &phxu0, const mat Q0, const vec R0 ) {
439        pfxu = pfxu0;
440        phxu = phxu0;
441       
442        set_dim ( pfxu0->_dimx() );
443        dimy = phxu0->dimension();
444        dimc = pfxu0->_dimu();
445       
446        vec &_mu = est._mu();
447        // if mu is not set, set it to zeros, just for constant terms of A and C
448        if ( _mu.length() != dimension() ) _mu = zeros ( dimension() );
449        A = zeros ( dimension(), dimension() );
450        C = zeros ( dimy, dimension() );
451       
452        //initialize matrices A C, later, these will be only updated!
453        pfxu->dfdx_cond ( _mu, zeros ( dimc ), A, true );
454        //      pfxu->dfdu_cond ( *_mu,zeros ( dimu ),B,true );
455        phxu->dfdx_cond ( _mu, zeros ( dimc ), C, true );
456        //      phxu->dfdu_cond ( *_mu,zeros ( dimu ),D,true );
457       
458        R = R0;
459        Q = Q0;
460       
461        //
462}
463
464
465void EKF_UD::bayes ( const vec &yt, const vec &cond ) {
466        //preparatory
467        vec &_mu=est._mu();
468        const vec &u=cond;
469        int dim = dimension();
470       
471        U = est._R()._L().T();
472        D =  est._R()._D();
473       
474        ////////////
475       
476        pfxu->dfdx_cond ( _mu, u, A, false ); //update A by a derivative of fx
477        phxu->dfdx_cond ( _mu, u, C, false ); //update A by a derivative of fx
478
479        mat PhiU = A*U;
480       
481        vec Din = D;
482        int i,j,k;
483        double sigma;
484        mat G = eye(dim);
485        //////// thorton
486       
487        //move mean;
488        _mu = pfxu->eval(_mu,u);
489       
490        for (i=dim-1; i>=0;i--){
491                sigma = 0.0; 
492                for (j=0; j<dim; j++) {
493                        sigma += PhiU(i,j)*PhiU(i,j) *Din(j); 
494                        sigma += G(i,j)*G(i,j) * Q(j,j); 
495                }
496                D(i) = sigma; 
497                //if (D(i)>(1<<17)) D(i)=(1<<17);
498                /////////////// !!!!!!!!!!!!!
499               
500                for (j=0;j<i;j++){ 
501                        sigma = 0.0; 
502                        for (k=0;k<dim;k++){ 
503                                sigma += PhiU(i,k)*Din(k)*PhiU(j,k); 
504                        }
505                        for (k=0;k<dim;k++){ 
506                                sigma += G(i,k)*Q(k,k)*G(j,k); 
507                        }
508                        //
509                        U(j,i) = sigma/D(i); 
510                        for (k=0;k<dim;k++){ 
511                                PhiU(j,k) = PhiU(j,k) - U(j,i)*PhiU(i,k); 
512                        }
513                        for (k=0;k<dim;k++){ 
514                                G(j,k) = G(j,k) - U(j,i)*G(i,k); 
515                        }
516                }
517        }
518
519        if ( log_level[logU] ){
520                //vec tmp=vec(U._data(),dimension()*dimension());
521                log_level.store(logU, vec(U._data(),dimension()*dimension()));
522               
523                log_level.store(logD, D);
524               
525//              log_level.store(logA, vec(A._data(), dimension()*dimension()));
526
527//              log_level.store(logC, (1<<15)*vec(C._data(), dimension()*dimensiony()));
528               
529        }
530        if ( log_level[logG] ){
531                log_level.store(logG,vec(G._data(),dimension()*dimension()));
532        }
533        //cout << "Ut: " << U << endl;
534        //cout << "Dt: " << D << endl;
535        // bierman
536       
537        double dz,alpha,gamma,beta,lambda;
538        vec a;
539        vec b;
540        vec yp = phxu->eval(_mu,u);
541        for (int iy=0; iy<dimy; iy++){
542                a     = U.T()*C.get_row(iy);    // a is not modified, but
543                b     = elem_mult(D,a);                          // b is modified to become unscaled Kalman gain.
544                dz    = yt(iy) - yp(iy); 
545
546                alpha = R(iy); 
547                gamma = 1/alpha; 
548                for (j=0;j<dim;j++){
549                        beta   = alpha; 
550                        alpha  = alpha + a(j)*b(j); 
551                        lambda = -a(j)*gamma; 
552                        gamma  = 1.0/alpha; 
553                        D(j) = beta*gamma*D(j); 
554
555//                      cout << "a: " << alpha << "g: " << gamma << endl;
556                        for (i=0;i<j;i++){
557                                beta   = U(i,j); 
558                                U(i,j) = beta + b(i)*lambda; 
559                                b(i)   = b(i) + b(j)*beta; 
560                        }
561                }
562                double dzs = gamma*dz;  // apply scaling to innovations
563                _mu   = _mu + dzs*b; // multiply by unscaled Kalman gain
564               
565                //cout << "Ub: " << U << endl;
566                //cout << "Db: " << D << endl <<endl;
567               
568        }
569       
570        // artificial round of U
571//      U = round(U*(1<<15))/(1<<15);
572               
573/////
574est._R().__L()=U.T();
575est._R().__D()=D;
576
577        if ( evalll == true ) { //likelihood of observation y
578        }
579}
580
581void EKF_UD::from_setting ( const Setting &set ) {
582        BM::from_setting ( set );
583        shared_ptr<diffbifn> IM = UI::build<diffbifn> ( set, "IM", UI::compulsory );
584        shared_ptr<diffbifn> OM = UI::build<diffbifn> ( set, "OM", UI::compulsory );
585       
586        //statistics
587        int dim = IM->dimension();
588        vec mu0;
589        if ( !UI::get ( mu0, set, "mu0" ) )
590                mu0 = zeros ( dim );
591       
592        mat P0;
593        vec dP0;
594        if ( UI::get ( dP0, set, "dP0" ) )
595                P0 = diag ( dP0 );
596        else if ( !UI::get ( P0, set, "P0" ) )
597                P0 = eye ( dim );
598       
599        est._mu()=mu0;
600        est._R()=ldmat(P0);
601       
602        //parameters
603        vec dQ, dR;
604        UI::get ( dQ, set, "dQ", UI::compulsory );
605        UI::get ( dR, set, "dR", UI::compulsory );
606        set_parameters ( IM, OM, diag ( dQ ), dR  );
607       
608        UI::get(log_level, set, "log_level", UI::optional);
609}
610
611}
Note: See TracBrowser for help on using the browser.