root/library/bdm/estim/kalman.cpp @ 1187

Revision 1173, 17.1 kB (checked in by smidl, 14 years ago)

tests in kalman

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