root/applications/pmsm/simulator_zdenek/ekf_example/ekf_obj.cpp @ 1241

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

Choleski - covariance 0.01

  • Property svn:eol-style set to native
Line 
1
2#include <estim/kalman.h>
3
4#include "ekf_obj.h"
5#include "../simulator.h"
6
7double minQ(double Q){if (Q>1.0){ return 1.0;} else {return Q;};};
8
9void mat_to_int16(const imat &M, int16 *I){
10        for (int16 i=0;i<M.rows(); i++){
11                for (int16 j=0;j<M.cols(); j++){
12                        *I++ = M(i,j);
13                }
14        }
15}
16void vec_to_int16(const ivec &v, int16 *I){
17        for (int16 i=0;i<v.length(); i++){
18                        *I++ = v(i);
19        }
20}
21
22#ifdef XXX
23///////////////
24void EKFfixed::bayes(const vec &yt, const vec &ut){
25        ekf(yt(0),yt(1),ut(0),ut(1));
26       
27        vec xhat(4);   
28        //UGLY HACK!!! reliance on a predictor!!
29        xhat(0)=zprevod(x_est[0],Qm)*Iref;
30        xhat(1)=zprevod(x_est[1],Qm)*Iref;
31        xhat(2)=zprevod(x_est[2],Qm)*Wref;
32        xhat(3)=zprevod(x_est[3],15)*Thetaref;
33       
34        E.set_mu(xhat);
35       
36        if ( BM::evalll ) {
37/*              //from enorm
38                vec xdif(x,4);//first 4 of x
39                //UGLY HACK!!! reliance on a predictor!!
40/*              xdif(0)=x[0]-zprevod(x_pred[0],Qm)*Iref;
41                xdif(1)=x[1]-zprevod(x_pred[1],Qm)*Iref;
42                xdif(2)=x[2]-zprevod(x_pred[2],Qm)*Wref;
43                xdif(3)=x[3]-zprevod(x_pred[3],15);*/
44               
45//              xdif -=xhat; //(xdif=x-xhat)
46               
47                mat Pfull(4,4);
48                double* Pp=Pfull._data();
49                for(int16 i=0;i<16;i++){*(Pp++) = zprevod(P_est[i],15);}
50               
51                E._R()._M()=Pfull;
52               
53               
54//              BM::ll = -0.5* ( 4 * 1.83787706640935 +log ( det ( Pfull ) ) +xdif* ( inv(Pfull)*xdif ) );*/
55        }
56};
57
58
59void EKFfixed::update_psi(void)
60{
61  int16 t_sin,t_cos,tmp;
62
63  // implementace v PC
64  t_sin=prevod(sin(Thetaref*x_est[3]/32768.),15);
65  t_cos=prevod(cos(Thetaref*x_est[3]/32768.),15);
66
67  PSI[2]=((int32)cB*t_sin)>>15;
68  tmp=((int32)cH*x_est[2])>>15;
69  PSI[3]=((int32)tmp*t_cos)>>15;
70  PSI[6]=-((int32)cB*t_cos)>>15;
71  PSI[7]=((int32)tmp*t_sin)>>15;
72}
73
74
75void EKFfixed::prediction(int16 *ux)
76{
77  int16 t_sin,t_cos, tmp;
78
79  // implementace v PC
80  //t_sin=prevod(sin(Thetaref*x_est[3]/32768.),15);
81  //t_cos=prevod(cos(Thetaref*x_est[3]/32768.),15);
82
83  t_sin=prevod(sin(Thetaref*x_est[3]/Qm),15);
84  t_cos=prevod(cos(Thetaref*x_est[3]/Qm),15);
85 
86  tmp=((int32)cB*x_est[2])>>15;
87  x_pred[0]=((int32)cA*x_est[0]+(int32)tmp*t_sin+(int32)cC*ux[0])>>15;
88  x_pred[1]=((int32)cA*x_est[1]-(int32)tmp*t_cos+(int32)cC*ux[1])>>15;
89  x_pred[2]=x_est[2];
90  x_pred[3]=(((int32)x_est[3]<<15)+(int32)cG*x_est[2])>>15;
91
92  update_psi();
93
94  mmult(PSI,P_est,temp15a,3,3,3);
95//  mtrans(PSI,temp15b,5,5);
96  mmultt(temp15a,PSI,P_pred,3,3,3);
97  maddD(P_pred,Q,3,3);
98}
99
100void EKFfixed::correction(void)
101{
102  int16 Y_error[2];
103  int32 temp30a[4]; /* matrix [2,2] - temporary matrix for inversion */
104
105  choice_P(P_pred,temp15a,3);
106  maddD(temp15a,R,1,1);
107  minv2(temp15a,temp30a);
108        Ry(0,0)=zprevod(temp15a[0],15);
109        Ry(0,1)=zprevod(temp15a[1],15);
110        Ry(1,0)=zprevod(temp15a[2],15);
111        Ry(1,1)=zprevod(temp15a[3],15);
112 
113  mmultDr(P_pred,temp15a,3,3,1,1);
114  mmult1530(temp15a,temp30a,Kalm,3,1,1);
115
116
117  /* estimate the state system */
118  choice_x(x_pred, temp15a);
119  msub(Y_mes,temp15a,Y_error,1,0);
120  mmult(Kalm,Y_error,temp15a,3,1,0);
121  madd(x_pred,temp15a,x_est,3,0);
122
123  /* matrix of covariances - version without MMULTDL() */
124
125/* Version with MMULTDL() */
126  mmultDl(P_pred,temp15a,1,3,3,1);
127
128  mmult(Kalm,temp15a,P_est,3,1,3);
129  msub(P_pred,P_est,P_est,3,3);
130/* END */
131}
132
133
134void EKFfixed::ekf(double ux, double uy, double isxd, double isyd)
135{
136  // vypocet napeti v systemu (x,y)
137  ukalm[0]=prevod(ux/Uref,Qm);
138  ukalm[1]=prevod(uy/Uref,Qm);
139
140  // zadani mereni
141  Y_mes[0]=prevod(isxd/Iref,Qm);
142  Y_mes[1]=prevod(isyd/Iref,Qm);
143
144  ////// vlastni rutina EKF /////////////////////////
145  prediction(ukalm);
146  correction();
147
148  // navrat estimovanych hodnot regulatoru
149  vec& mu = E._mu();
150  (mu)(0)=zprevod(x_est[0],Qm)*Iref;
151  (mu)(1)=zprevod(x_est[1],Qm)*Iref;
152  (mu)(2)=zprevod(x_est[2],Qm)*Wref;
153  (mu)(3)=zprevod(x_est[3],15)*Thetaref;
154}
155
156void EKFfixed::init_ekf(double Tv)
157{
158  // Tuning of matrix Q
159  Q[0]=prevod(.05,15);       // 0.05
160  Q[5]=Q[0];
161  Q[10]=prevod(0.0002,15);      // 1e-3
162  Q[15]=prevod(0.001,15);      // 1e-3
163
164  // Tuning of matrix R
165  R[0]=prevod(0.1,15);         // 0.05
166  R[3]=R[0];
167
168  // Motor model parameters
169  cA=prevod(1-Tv*Rs/Ls,15);
170  cB=prevod(Tv*Wref*Fmag/Iref/Ls,15);
171  cC=prevod(Tv/Ls/Iref*Uref,15);
172//  cD=prevod(1-Tv*Bf/J,15);
173//  cE=prevod(kp*p*p*Tv*Fmag*Iref/J/Wref,15);
174//  cF=prevod(p*Tv*Mref/J/Wref,15);
175  cG=prevod(Tv*Wref*4/Thetaref,15);
176//  cH=prevod(Tv*Wref*Fmag/Iref/Ls*Thetaref,15);
177  cH=prevod(Tv*Wref*Fmag/Iref/Ls,17); //cB in q
178  //  cI=prevod(kp*p*p*Tv*Fmag*Iref/J/Wref*Thetaref);
179
180  /* Init matrix PSI with permanently constant terms */
181  PSI[0]=cA;
182  PSI[5]=PSI[0];
183  PSI[10]=0x7FFF;
184  PSI[14]=cG;
185  PSI[15]=0x7FFF;
186 
187  P_est[0]=0x7FFF;
188  P_est[5]=0x7FFF;
189  P_est[10]=0x7FFF;
190  P_est[15]=0x7FFF;
191}
192#endif
193
194void EKF_UDfix::set_parameters ( const shared_ptr<diffbifn> &pfxu0, const shared_ptr<diffbifn> &phxu0, const mat Q0, const vec R0 ) {
195        pfxu = pfxu0;
196        phxu = phxu0;
197       
198        set_dim ( pfxu0->_dimx() );
199        dimy = phxu0->dimension();
200        dimc = pfxu0->_dimu();
201       
202        vec &_mu = est._mu();
203        // if mu is not set, set it to zeros, just for constant terms of A and C
204        if ( _mu.length() != dimension() ) _mu = zeros ( dimension() );
205        A = zeros ( dimension(), dimension() );
206        C = zeros ( dimy, dimension() );
207       
208        //initialize matrices A C, later, these will be only updated!
209        pfxu->dfdx_cond ( _mu, zeros ( dimc ), A, true );
210        //      pfxu->dfdu_cond ( *_mu,zeros ( dimu ),B,true );
211        phxu->dfdx_cond ( _mu, zeros ( dimc ), C, true );
212        //      phxu->dfdu_cond ( *_mu,zeros ( dimu ),D,true );
213       
214        R = R0;
215        Q = Q0;
216       
217        //
218}
219// aux fnc
220void UDtof(const mat &U, const vec &D, imat &Uf, ivec &Df, const vec &xref){
221        mat P= U*diag(D)*U.T();
222        mat T = diag(1.0/(xref));
223        mat Pf = T*P*T;
224       
225        ldmat Pld(Pf);
226       
227        mat Ut=Pld._L().T()*(1<<15); // U is in q15 -- diagonal is 0!!!
228        Uf=round_i(Ut);
229        Df=round_i(Pld._D()*(1<<15));
230        ivec zer=find(Df==0);
231        for(int16 i=0; i<zer.length(); i++) Df(zer(i))=1;
232}
233
234
235void EKF_UDfix::bayes ( const vec &yt, const vec &cond ) {
236        //preparatory
237        vec &_mu=est._mu();
238        const vec &u=cond;
239        int16 dim = dimension();
240        ///// !!!!!!!!!!!!!!!!
241        U = est._R()._L().T();
242        D =  est._R()._D();
243       
244        ////////////
245       
246        pfxu->dfdx_cond ( _mu, u, A, false ); //update A by a derivative of fx
247        phxu->dfdx_cond ( _mu, u, C, false ); //update A by a derivative of fx
248       
249        mat PhiU = A*U;
250
251        //////
252/*      vec xref(4);
253        xref(0)= 30.0*1.4142;
254        xref(1)= 30.0*1.4142;
255        xref(2)= 6.283185*200.;
256        xref(3) = 3.141593;*/
257        //xref(4) = 34.0;
258       
259       
260        vec Din = D;
261        int16 i,j,k;
262        double sigma;
263        mat G = eye(dim);
264        //////// thorton
265       
266        //move mean;
267        _mu = pfxu->eval(_mu,u);
268       
269        for (i=dim-1; i>=0;i--){
270                sigma = 0.0; 
271                for (j=0; j<dim; j++) {
272                        sigma += PhiU(i,j)*PhiU(i,j) *Din(j); 
273                        sigma += G(i,j)*G(i,j) * Q(j,j);
274                }
275               
276/*              double sigma2= 0.0;
277                for (j=0; j<dim; j++) {
278                        sigma2 += PhiU(i,j)*PhiU(i,j) *Din(j);
279                }
280                sigma2 +=Q(i,i);//*G(i,i)=1.0
281                for (j=i+1; j<dim; j++) {
282                        sigma2 += G(i,j)*G(i,j) * Q(j,j);
283                }*/
284                D(i) = sigma; 
285               
286        /*      UDtof(U,D,Utf,Dtf,xref);
287                cout << "d=sig"<<endl;
288                cout << Dtf << endl;
289        */     
290                for (j=0;j<i;j++){ 
291//                      cout << i << "," << j << endl;
292                        sigma = 0.0; 
293                        for (k=0;k<dim;k++){ 
294                                sigma += PhiU(i,k)*Din(k)*PhiU(j,k); 
295                        }
296                        for (k=0;k<dim;k++){ 
297                                sigma += G(i,k)*Q(k,k)*G(j,k); 
298                        }
299                        //
300                        U(j,i) = sigma/D(i); 
301                       
302/*                      cout << "U=sig/D"<<endl;
303                        UDtof(U,D,Utf,Dtf,xref);
304                        cout << Utf << endl << Dtf << endl;
305                        cout << G << endl << Din << endl<<endl;*/
306                       
307                        for (k=0;k<dim;k++){ 
308                                PhiU(j,k) = PhiU(j,k) - U(j,i)*PhiU(i,k); 
309                        }
310                        for (k=0;k<dim;k++){ 
311                                G(j,k) = G(j,k) - U(j,i)*G(i,k); 
312                        }
313               
314                }
315        }
316
317        // bierman
318       
319        double dz,alpha,gamma,beta,lambda;
320        vec a;
321        vec b;
322        vec yp = phxu->eval(_mu,u);
323        vec xp=_mu; // used in bierman
324       
325       
326        for (int16 iy=0; iy<dimy; iy++){
327                a     = U.T()*C.get_row(iy);    // a is not modified, but
328                b     = elem_mult(D,a);                          // b is modified to become unscaled Kalman gain.
329                dz    = yt(iy) - yp(iy); 
330               
331                alpha = R(iy); 
332                gamma = 1/alpha; 
333                for (j=0;j<dim;j++){
334                        beta   = alpha; 
335                        alpha  = alpha + a(j)*b(j); 
336                        lambda = -a(j)*gamma; 
337                        gamma  = 1.0/alpha; 
338                        D(j) = beta*gamma*D(j); 
339                       
340                        //                      cout << "a: " << alpha << "g: " << gamma << endl;
341                        for (i=0;i<j;i++){
342                                beta   = U(i,j); 
343                                U(i,j) = beta + b(i)*lambda; 
344                                b(i)   = b(i) + b(j)*beta; 
345                        }
346                }
347                double dzs = gamma*dz;  // apply scaling to innovations
348                _mu   = _mu + dzs*b; // multiply by unscaled Kalman gain
349                //cout << "Ub: " << U << endl;
350                //cout << "Db: " << D << endl <<endl;
351               
352        }
353
354       
355        /////
356        est._R().__L()=U.T();
357        est._R().__D()=D;
358       
359        if ( evalll == true ) { //likelihood of observation y
360        }
361}
362
363void EKF_UDfix::from_setting ( const Setting &set ) {
364        BM::from_setting ( set );
365        shared_ptr<diffbifn> IM = UI::build<diffbifn> ( set, "IM", UI::compulsory );
366        shared_ptr<diffbifn> OM = UI::build<diffbifn> ( set, "OM", UI::compulsory );
367       
368        //statistics
369        int16 dim = IM->dimension();
370        vec mu0;
371        if ( !UI::get ( mu0, set, "mu0" ) )
372                mu0 = zeros ( dim );
373       
374        mat P0;
375        vec dP0;
376        if ( UI::get ( dP0, set, "dP0" ) )
377                P0 = diag ( dP0 );
378        else if ( !UI::get ( P0, set, "P0" ) )
379                P0 = eye ( dim );
380       
381        est._mu()=mu0;
382        est._R()=ldmat(P0);
383       
384        //parameters
385        vec dQ, dR;
386        UI::get ( dQ, set, "dQ", UI::compulsory );
387        UI::get ( dR, set, "dR", UI::compulsory );
388        set_parameters ( IM, OM, diag ( dQ ), dR  );
389       
390        UI::get(log_level, set, "log_level", UI::optional);
391}
392
393
394void EKFfixedUD::bayes(const itpp::vec& yt, const itpp::vec& ut)
395{
396ekf(ut[0],ut[1],yt[0],yt[1]);
397}
398
399
400void EKFfixedUD::ekf(double ux, double uy, double isxd, double isyd)
401{
402        // vypocet napeti v systemu (x,y)
403        int16 uf[2];
404        uf[0]=prevod(ux/Uref,15);
405        uf[1]=prevod(uy/Uref,15);
406       
407        int16 Y_mes[2];
408        // zadani mereni
409        Y_mes[0]=prevod(isxd/Iref,15);
410        Y_mes[1]=prevod(isyd/Iref,15);
411               
412        ////// vlastni rutina EKF -- /////////////////////////
413        int16 t_sin,t_cos;
414        int32 tmp;
415       
416        // implementace v PC
417/*      t_sin=prevod(sin(Thetaref*x_est[3]/32768.),15);
418        t_cos=prevod(cos(Thetaref*x_est[3]/32768.),15);*/
419        tmp=prevod(sin(Thetaref*x_est[3]/32768.),15);
420        if (tmp>32767) t_sin =32767; else t_sin=tmp;
421        tmp=prevod(cos(Thetaref*x_est[3]/32768.),15);
422        if (tmp>32767) t_cos =32767; else t_cos=tmp;
423       
424        tmp=((int32)cB*x_est[2])>>15; // q15*q13 -> q13
425        //             q15*q13 +          q13*q15      + q15*q13??
426        x_est[0]=((int32)cA*x_est[0]+(int32)tmp*t_sin+(int32)cC*(uf[0]))>>15; // in q13
427        x_est[1]=((int32)cA*x_est[1]-(int32)tmp*t_cos+(int32)cC*(uf[1]))>>15; // q13
428        x_est[2]=x_est[2];
429        //               q15              + q15*q13
430        x_est[3]=(((int32)x_est[3]<<15)+(int32)cG*x_est[2])>>15;
431       
432        if(x_est[3]>(1<<15)) x_est[3]-=2*(1<<15);
433        if(x_est[3]<-(1<<15)) x_est[3]+=2*(1<<15);
434       
435        //void EKFfixed::update_psi(void)
436        {       
437                PSI[2]=((int32)cB*t_sin)>>15;
438                tmp=((int32)cH*x_est[2])>>15; // ! cH =cB with different scale!!
439                PSI[3]=((int32)tmp*t_cos)>>15;
440                PSI[6]=-((int32)cB*t_cos)>>15;
441                PSI[7]=((int32)tmp*t_sin)>>15;
442        }
443        {
444                int Ai[16];
445                for (int i=0;i<16; i++) Ai[i]=(int)(PSI[i]);
446                ivec Ad(Ai,16);
447                log_level.store(logA,get_from_ivec(Ad));
448        }
449       
450        ///////// end of copy ///////////////
451        mmultAU(PSI,Uf,PSIU,4,4);
452        //thorton(int16 *U, int16 *D, int16 *PSIU, int16 *Q, int16 *G, int16 *Dold, unsigned int16 dimx);
453        thorton_fast(Uf,Df,PSIU,Q,G,Dfold,4);
454       
455         {
456                 int Ui[16];    for (int i=0;i<16; i++) Ui[i]=(int)Uf[i];
457                 
458                ivec Ud(Ui,16);
459                log_level.store(logU,get_from_ivec(Ud));
460               
461                int di[16];     for (int i=0;i<16; i++) di[i]=(int)Df[i];
462                ivec dd(di,4);
463                imat U(Ui,4,4);
464                mat U2=to_mat(U)/32768;
465                mat PP=U2*diag(to_vec(dd))*U2.T();
466                vec pp(PP._data(),16);
467                log_level.store(logP,pp);
468         }
469         {
470                 int Gi[16];    for (int i=0;i<16; i++) Gi[i]=(int)G[i];
471                 ivec Gd(Gi,16);
472                log_level.store(logG,get_from_ivec(Gd));
473        }
474       
475       
476        int16 difz[2];
477        difz[0]=(Y_mes[0]-x_est[0]); // shift to q15!!
478        difz[1]=(Y_mes[1]-x_est[1]);//<<2;
479
480        {
481                int Di[4];               for (int i=0;i<4; i++) Di[i]=(int)Df[i];
482                ivec dd(Di,4);//;dd(0)=Y_mes[0];dd(1)=Y_mes[1]; dd(2)=difz[0]; dd(3)=difz[1];
483                log_level.store(logD,get_from_ivec(dd));
484        }
485       
486        //bierman(int16 *difz, int16 *xp, int16 *U, int16 *D, int16 *R, unsigned int16 dimy, unsigned int16 dimx );
487        int16 dR[2];dR[0]=R[0];dR[1]=R[3];
488        //int16 xb[4]; xb[0]=x_est[0]<<2; xb[1]=x_est[1]<<2;  xb[2]=x_est[2]<<2;  xb[3]=x_est[3]; 
489        bierman_fast(difz,x_est,Uf,Df,dR,2,4);
490        //x_est[0] = xb[0]>>2; x_est[1]=xb[1]>>2; x_est[2]=xb[2]>>2; x_est[3]=xb[3];
491       
492        // navrat estimovanych hodnot regulatoru
493        vec& mu = E._mu();
494        (mu)(0)=zprevod(x_est[0],15)*Iref;
495        (mu)(1)=zprevod(x_est[1],15)*Iref;
496        (mu)(2)=zprevod(x_est[2],15)*Wref;
497        (mu)(3)=zprevod(x_est[3],15)*Thetaref;
498       
499        //x_est[2]=x[2]*32768/Wref;
500        //x_est[3]=x[3]*32768/Thetaref;
501//      mat T=diag(concat(vec_2(Iref,Iref), vec_2(Wref,Thetaref)));
502}
503
504void EKFfixedUD::init_ekf(double Tv)
505{
506        // Tuning of matrix Q
507        Q[0]=prevod(.01,15);       // 0.05
508        Q[5]=Q[0];
509        Q[10]=prevod(0.0005,15);      // 1e-3
510        Q[15]=prevod(0.001,15);      // 1e-3
511
512        Uf[0]=0x7FFF;// !       // 0.05
513        Uf[1]=Uf[2]=Uf[3]=Uf[4]=0;
514        Uf[5]=0x7FFF;//!
515        Uf[6]=Uf[6]=Uf[8]=Uf[9]=0;
516        Uf[10]=0x7FFF;//!      // 1e-3
517        Uf[11]=Uf[12]=Uf[13]=Uf[4]=0;
518        Uf[15]=0x7FFF;      // 1e-3
519       
520        Df[0]=1<<14;
521        Df[1]=1<<14;
522        Df[2]=1<<14;
523        Df[3]=1<<14;
524       
525        // Tuning of matrix R
526        R[0]=prevod(0.05,15);         // 0.05
527        R[3]=R[0];
528       
529        // Motor model parameters
530        cA=prevod(1-Tv*Rs/Ls,15);
531        cB=prevod(Tv*Wref*Fmag/Iref/Ls,15);
532        cC=prevod(Tv/Ls/Iref*Uref,15);
533        //  cD=prevod(1-Tv*Bf/J,15);
534        //  cE=prevod(kp*p*p*Tv*Fmag*Iref/J/Wref,15);
535        //  cF=prevod(p*Tv*Mref/J/Wref,15);
536        cG=prevod(Tv*Wref/Thetaref,15); //in q15!
537        //cH=prevod(Tv*Wref*Fmag/Iref/Ls*Thetaref,15);
538        // cH=prevod(Tv*Wref*Thetaref*Fmag/Iref/Ls,15); // < ===========
539        cH=prevod(Tv*Wref*Thetaref*Fmag/Iref/Ls,15); //
540        //  cI=prevod(kp*p*p*Tv*Fmag*Iref/J/Wref*Thetaref);
541       
542        /* Init matrix PSI with permanently constant terms */
543        PSI[0]=cA;
544        PSI[5]=PSI[0];
545        PSI[10]=0x7FFF;
546        PSI[14]=cG;
547        PSI[15]=0x7FFF;
548       
549}
550
551void EKFfixedCh::bayes(const itpp::vec& yt, const itpp::vec& ut)
552{
553        ekf(ut[0],ut[1],yt[0],yt[1]);
554}
555
556
557void EKFfixedCh::ekf(double ux, double uy, double isxd, double isyd)
558{
559        // vypocet napeti v systemu (x,y)
560        int16 uf[2];
561        uf[0]=prevod(ux/Uref,15);
562        uf[1]=prevod(uy/Uref,15);
563       
564        int16 Y_mes[2];
565        // zadani mereni
566        Y_mes[0]=prevod(isxd/Iref,15);
567        Y_mes[1]=prevod(isyd/Iref,15);
568       
569        ////// vlastni rutina EKF -- /////////////////////////
570        int16 t_sin,t_cos, tmp;
571       
572        // implementace v PC
573        /*      t_sin=prevod(sin(Thetaref*x_est[3]/32768.),15);
574         *      t_cos=prevod(cos(Thetaref*x_est[3]/32768.),15);*/
575        t_sin=prevod(sin(Thetaref*x_est[3]/32768.),15);
576        t_cos=prevod(cos(Thetaref*x_est[3]/32768.),15);
577       
578        tmp=((int32)cB*x_est[2])>>15; // q15*q13 -> q13
579        //             q15*q13 +          q13*q15      + q15*q13??
580        x_est[0]=((int32)cA*x_est[0]+(int32)tmp*t_sin+(int32)cC*(uf[0]))>>15; // in q13
581        x_est[1]=((int32)cA*x_est[1]-(int32)tmp*t_cos+(int32)cC*(uf[1]))>>15; // q13
582        x_est[2]=x_est[2];
583        //               q15              + q15*q13
584        x_est[3]=(((int32)x_est[3]<<15)+(int32)cG*x_est[2])>>15;
585       
586        if(x_est[3]>(1<<15)) x_est[3]-=2*(1<<15);
587        if(x_est[3]<-(1<<15)) x_est[3]+=2*(1<<15);
588       
589        //void EKFfixed::update_psi(void)
590        {       
591                PSI[2]=((int32)cB*t_sin)>>15;
592                tmp=((int32)cH*x_est[2])>>15; // ! cH =cB with different scale!!
593                PSI[3]=((int32)tmp*t_cos)>>15;
594                PSI[6]=-((int32)cB*t_cos)>>15;
595                PSI[7]=((int32)tmp*t_sin)>>15;
596        }
597        {
598/*              ivec Ad(PSI,16);
599                log_level.store(logA,get_from_ivec(Ad));*/
600        }
601       
602        ///////// end of copy ///////////////
603        mmultACh(PSI,Chf,PSICh,4,4);
604        //thorton(int16 *U, int16 *D, int16 *PSIU, int16 *Q, int16 *G, int16 *Dold, unsigned int16 dimx);
605        householder(PSICh,Q,4);
606        //givens(PSICh,Q,4);
607        // COPY
608        for (int16 ii=0; ii<16; ii++){Chf[ii]=PSICh[ii];}
609       
610        {
611               
612                int Chi[16]; for (int i=0;i<16;i++) Chi[i]=Chf[i];
613                ivec Chd(Chi,16);
614                log_level.store(logCh,get_from_ivec(Chd));
615                imat mCh(Chd._data(), 4,4);
616                imat P = mCh*mCh.T();
617                ivec iP(P._data(),16);
618                log_level.store(logP,get_from_ivec(iP));
619        }
620       
621       
622        int16 difz[2];
623        difz[0]=(Y_mes[0]-x_est[0]); // shift to q15!!
624        difz[1]=(Y_mes[1]-x_est[1]);//<<2;
625       
626       
627        //
628        int16 dR[2];dR[0]=R[0];dR[1]=R[3];
629        //int16 xb[4]; xb[0]=x_est[0]<<2; xb[1]=x_est[1]<<2;  xb[2]=x_est[2]<<2;  xb[3]=x_est[3]; 
630
631        carlson(difz,x_est,Chf,dR,2,4);
632        //x_est[0] = xb[0]>>2; x_est[1]=xb[1]>>2; x_est[2]=xb[2]>>2; x_est[3]=xb[3];
633       
634        // navrat estimovanych hodnot regulatoru
635        vec& mu = E._mu();
636        (mu)(0)=zprevod(x_est[0],15)*Iref;
637        (mu)(1)=zprevod(x_est[1],15)*Iref;
638        (mu)(2)=zprevod(x_est[2],15)*Wref;
639        (mu)(3)=zprevod(x_est[3],15)*Thetaref;
640       
641        //x_est[2]=x[2]*32768/Wref;
642        //x_est[3]=x[3]*32768/Thetaref;
643        //      mat T=diag(concat(vec_2(Iref,Iref), vec_2(Wref,Thetaref)));
644}
645
646void EKFfixedCh::init_ekf(double Tv)
647{
648        // Tuning of matrix Q
649        Q[0]=prevod(.01,15);       // 0.05
650        Q[5]=Q[0];
651        Q[10]=prevod(0.01,15);      // 1e-3
652        Q[15]=prevod(0.01,15);      // 1e-3
653       
654        Chf[0]=0x3FFF;// !       // 0.05
655        Chf[1]=Chf[2]=Chf[3]=Chf[4]=0;
656        Chf[5]=0x3FFF;//!
657        Chf[6]=Chf[6]=Chf[8]=Chf[9]=0;
658        Chf[10]=0x3FFF;//!      // 1e-3
659        Chf[11]=Chf[12]=Chf[13]=Chf[4]=0;
660        Chf[15]=0x3FFF;      // 1e-3
661               
662        // Tuning of matrix R
663        R[0]=prevod(0.05,15);         // 0.05
664        R[3]=R[0];
665       
666        // Motor model parameters
667        cA=prevod(1-Tv*Rs/Ls,15);
668        cB=prevod(Tv*Wref*Fmag/Iref/Ls,15);
669        cC=prevod(Tv/Ls/Iref*Uref,15);
670        //  cD=prevod(1-Tv*Bf/J,15);
671        //  cE=prevod(kp*p*p*Tv*Fmag*Iref/J/Wref,15);
672        //  cF=prevod(p*Tv*Mref/J/Wref,15);
673        cG=prevod(Tv*Wref/Thetaref,15); //in q15!
674        //cH=prevod(Tv*Wref*Fmag/Iref/Ls*Thetaref,15);
675        //      cH=prevod(Tv*Wref*Thetaref*Fmag/Iref/Ls,15); //  <=======
676        cH=prevod(Tv*Wref*Fmag/Iref/Ls,15); //
677        //  cI=prevod(kp*p*p*Tv*Fmag*Iref/J/Wref*Thetaref);
678       
679        /* Init matrix PSI with permanently constant terms */
680        PSI[0]=cA;
681        PSI[5]=PSI[0];
682        PSI[10]=0x7FFF;
683        PSI[14]=cG;
684        PSI[15]=0x7FFF;
685       
686}
Note: See TracBrowser for help on using the browser.