root/library/bdm/estim/arx.cpp @ 1038

Revision 1038, 12.0 kB (checked in by smidl, 14 years ago)

small steps in partforg...

  • Property svn:eol-style set to native
Line 
1#include "arx.h"
2namespace bdm {
3
4void ARX::bayes_weighted ( const vec &yt, const vec &cond, const double w ) {
5        bdm_assert_debug ( yt.length() == dimy, "BM::bayes yt is of size "+num2str(yt.length())+" expected dimension is "+num2str(dimy) );
6        bdm_assert_debug ( cond.length() == rgrlen , "BM::bayes cond is of size "+num2str(cond.length())+" expected dimension is "+num2str(rgrlen) );
7       
8        BMEF::bayes_weighted(yt,cond,w); //potential discount scheduling
9       
10        double lnc;
11        //cache
12        ldmat &V = est._V();
13        double &nu = est._nu();
14
15        dyad.set_subvector ( 0, yt );
16        if (cond.length()>0)
17                dyad.set_subvector ( dimy, cond );
18        // possible "1" is there from the beginning
19
20        if ( frg < 1.0 ) {
21                est.pow ( frg ); // multiply V and nu
22
23
24                //stabilize
25                ldmat V0 = alter_est._V(); //$ copy
26                double &nu0 = alter_est._nu();
27
28                V0 *= ( 1 - frg );
29                V += V0; //stabilization
30                nu += ( 1 - frg ) * nu0;
31
32                // recompute loglikelihood of new "prior"
33                if ( evalll ) {
34                        last_lognc = est.lognc();
35                }
36        }
37        V.opupdt ( dyad, w );
38        nu += w;
39
40        // log(sqrt(2*pi)) = 0.91893853320467
41        if ( evalll ) {
42                lnc = est.lognc();
43                ll = lnc - last_lognc - 0.91893853320467;
44                last_lognc = lnc;
45        }
46}
47
48double ARX::logpred ( const vec &yt, const vec &cond ) const {
49        egiw pred ( est );
50        ldmat &V = pred._V();
51        double &nu = pred._nu();
52
53        double lll;
54        vec dyad_p = dyad;
55        dyad_p.set_subvector ( 0, yt );
56        dyad_p.set_subvector(dimy,cond);
57       
58        if ( frg < 1.0 ) {
59                pred.pow ( frg );
60                lll = pred.lognc();
61        } else//should be save: last_lognc is changed only by bayes;
62                if ( evalll ) {
63                        lll = last_lognc;
64                } else {
65                        lll = pred.lognc();
66                }
67
68        V.opupdt ( dyad_p, 1.0 );
69        nu += 1.0;
70        // log(sqrt(2*pi)) = 0.91893853320467
71        return pred.lognc() - lll - 0.91893853320467;
72}
73
74void ARX::flatten ( const BMEF* B , double weight =1.0) {
75        const ARX* A = dynamic_cast<const ARX*> ( B );
76        // nu should be equal to B.nu
77        est.pow ( A->posterior()._nu() / posterior()._nu() *weight);
78        if ( evalll ) {
79                last_lognc = est.lognc();
80        }
81}
82
83ARX* ARX::_copy ( ) const {
84        ARX* Tmp = new ARX ( *this );
85        return Tmp;
86}
87
88void ARX::set_statistics ( const BMEF* B0 ) {
89        const ARX* A0 = dynamic_cast<const ARX*> ( B0 );
90
91        bdm_assert_debug ( dimension() == A0->dimension(), "Statistics of different dimensions" );
92        set_statistics ( A0->dimensiony(), A0->posterior()._V(), A0->posterior()._nu() );
93}
94
95enorm<ldmat>* ARX::epredictor ( const vec &cond ) const {
96        bdm_assert_debug ( cond.length() == rgrlen , "ARX::epredictor cond is of size "+num2str(cond.length())+" expected dimension is "+num2str(rgrlen) );
97       
98        mat mu ( dimy, posterior()._V().rows() - dimy );
99        mat R ( dimy, dimy );
100
101        vec ext_rgr;
102        if (have_constant){
103                ext_rgr = concat(cond,vec_1(1.0));
104        } else {
105                ext_rgr = cond;
106        }
107       
108        enorm<ldmat>* tmp;
109        tmp = new enorm<ldmat> ( );
110        //TODO: too hackish
111        if ( yrv._dsize() > 0 ) {
112        }
113
114        est.mean_mat ( mu, R ); //mu =
115        //correction for student-t  -- TODO check if correct!!
116        //R*=nu/(nu-2);
117        if (mu.cols()>0) {// nonempty egiw
118                mat p_mu = mu.T() * ext_rgr;    //the result is one column
119                tmp->set_parameters ( p_mu.get_col ( 0 ), ldmat ( R ) );
120        } else {
121                tmp->set_parameters ( zeros( R.rows() ), ldmat ( R ) );
122        }
123        if (dimy==yrv._dsize())
124                tmp->set_rv(yrv);
125        return tmp;
126}
127
128mlstudent* ARX::predictor_student ( ) const {
129        const ldmat &V = posterior()._V();
130
131        mat mu ( dimy, V.rows() - dimy );
132        mat R ( dimy, dimy );
133        mlstudent* tmp;
134        tmp = new mlstudent ( );
135
136        est.mean_mat ( mu, R ); //
137        mu = mu.T();
138
139        int end = V._L().rows() - 1;
140        ldmat Lam ( V._L() ( dimy, end, dimy, end ), V._D() ( dimy, end ) );  //exp val of R
141
142
143        if ( have_constant ) { // no constant term
144                //Assume the constant term is the last one:
145                if ( mu.cols() > 1 ) {
146                        tmp->set_parameters ( mu.get_cols ( 0, mu.cols() - 2 ), mu.get_col ( mu.cols() - 1 ), ldmat ( R ), Lam );
147                } else {
148                        tmp->set_parameters ( mat ( dimy, dimc ), mu.get_col ( mu.cols() - 1 ), ldmat ( R ), Lam );
149                }
150        } else {
151                // no constant term
152                tmp->set_parameters ( mu, zeros ( dimy ), ldmat ( R ), Lam );
153        }
154        return tmp;
155}
156
157
158
159/*! \brief Return the best structure
160@param Eg a copy of GiW density that is being examined
161@param Eg0 a copy of prior GiW density before estimation
162@param Egll likelihood of the current Eg
163@param indices current indices
164\return best likelihood in the structure below the given one
165*/
166double egiw_bestbelow ( egiw Eg, egiw Eg0, double Egll, ivec &indices ) { //parameter Eg is a copy!
167        ldmat Vo = Eg._V(); //copy
168        ldmat Vo0 = Eg._V(); //copy
169        ldmat& Vp = Eg._V(); // pointer into Eg
170        ldmat& Vp0 = Eg._V(); // pointer into Eg
171        int end = Vp.rows() - 1;
172        int i;
173        mat Li;
174        mat Li0;
175        double maxll = Egll;
176        double tmpll = Egll;
177        double belll = Egll;
178
179        ivec tmpindices;
180        ivec maxindices = indices;
181
182
183        cout << "bb:(" << indices << ") ll=" << Egll << endl;
184
185        //try to remove only one rv
186        for ( i = 0; i < end; i++ ) {
187                //copy original
188                Li = Vo._L();
189                Li0 = Vo0._L();
190                //remove stuff
191                Li.del_col ( i + 1 );
192                Li0.del_col ( i + 1 );
193                Vp.ldform ( Li, Vo._D() );
194                Vp0.ldform ( Li0, Vo0._D() );
195                tmpll = Eg.lognc() - Eg0.lognc(); // likelihood is difference of norm. coefs.
196
197                cout << "i=(" << i << ") ll=" << tmpll << endl;
198
199                //
200                if ( tmpll > Egll ) { //increase of the likelihood
201                        tmpindices = indices;
202                        tmpindices.del ( i );
203                        //search for a better match in this substructure
204                        belll = egiw_bestbelow ( Eg, Eg0, tmpll, tmpindices );
205                        if ( belll > maxll ) { //better match found
206                                maxll = belll;
207                                maxindices = tmpindices;
208                        }
209                }
210        }
211        indices = maxindices;
212        return maxll;
213}
214
215ivec ARX::structure_est ( const egiw &est0 ) {
216        ivec ind = linspace ( 1, est.dimension() - 1 );
217        egiw_bestbelow ( est, est0, est.lognc() - est0.lognc(), ind );
218        return ind;
219}
220
221
222
223ivec ARX::structure_est_LT ( const egiw &est0 ) {
224        //some stuff with beliefs etc.
225        ivec belief = vec_1 ( 2 );        // default belief
226        int nbest = 1;           // nbest: how many regressors are returned
227        int nrep = 5;         // nrep: number of random repetions of structure estimation
228        double lambda   = 0.9;
229        int k = 2;
230       
231        Array<str_aux> o2;
232       
233        ivec ind = bdm::straux1(est._V(),est._nu(), est0._V(), est0._nu(), belief, nbest, nrep, lambda, k, o2);
234       
235        return ind;
236}
237
238void ARX::from_setting ( const Setting &set ) {
239        BMEF::from_setting(set);
240       
241        UI::get (rgr, set, "rgr", UI::compulsory );
242       
243        dimy = yrv._dsize();
244        bdm_assert(dimy>0,"ARX::yrv should not be empty");
245        rgrlen = rgr._dsize();
246
247        int constant;
248        if ( !UI::get ( constant, set, "constant", UI::optional ) ) {
249                have_constant = true;
250        } else {
251                have_constant = constant > 0;
252        }
253        dimc = rgrlen;
254        rvc = rgr;
255
256        //init
257        shared_ptr<egiw> pri = UI::build<egiw> ( set, "prior", UI::optional );
258        if (pri){
259                set_prior(pri.get());
260        } else {
261                shared_ptr<egiw> post = UI::build<egiw> ( set, "posterior", UI::optional );
262                set_prior(post.get());
263        }
264               
265       
266        shared_ptr<egiw> alt = UI::build<egiw> ( set, "alternative", UI::optional );
267        if ( alt ) {
268                bdm_assert ( alt->_dimx() == dimy, "alternative is not compatible" );
269                bdm_assert ( alt->_V().rows() == dimy + rgrlen + int(have_constant==true), "alternative is not compatible" );
270                alter_est.set_parameters ( alt->_dimx(), alt->_V(), alt->_nu() );
271                alter_est.validate();
272        } 
273        // frg handled by BMEF
274
275}
276
277void ARX::set_prior(const epdf *pri){
278        const egiw * eg=dynamic_cast<const egiw*>(pri);
279        if ( eg ) {
280                bdm_assert ( eg->_dimx() == dimy, "prior is not compatible" );
281                bdm_assert ( eg->_V().rows() == dimy + rgrlen + int(have_constant==true), "prior is not compatible" );
282                est.set_parameters ( eg->_dimx(), eg->_V(), eg->_nu() );
283                est.validate();
284        } else {
285                est.set_parameters ( dimy, zeros ( dimy + rgrlen +int(have_constant==true)) );
286                set_prior_default ( est );
287        }
288        //check alternative
289        if (alter_est.dimension()!=dimension()){
290                alter_est = est;
291        }
292}
293
294void ARXpartialforg::bayes ( const vec &val, const vec &cond ) {
295    #define LOG2  0.69314718055995
296    vec frg = cond.right(cond.length() - rgrlen);
297    vec cond_rgr = cond.left(rgrlen);       // regression vector
298
299    int dimV = est._V().cols();
300    int nparams = dimV - 1;                 // number of parameters
301    int nalternatives = 1 << nparams;    // number of alternatives
302
303    // Permutation matrix
304    mat perm_matrix = ones(nalternatives, nparams);
305    int i, j, period, idx_from, idx_to, start, end;
306    for(i = 0; i < nparams; i++) {
307                idx_from = 1 << i;
308                period   = ( idx_from << 1 );
309        idx_to   = period - 1;
310                j = 0;
311        start = idx_from;
312        end = idx_to;
313                while ( start < nalternatives ) {
314            perm_matrix.set_submatrix(start, end, i, i, 0);
315            j++;
316            start = idx_from + period * j;
317            end = idx_to + period * j;
318        }
319    }
320
321    // Array of egiws for approximation
322    Array<egiw*> egiw_array(nalternatives + 1);
323    // No. of conditioning rows in LD
324    int nalternatives_cond, position;
325
326    for(int i = 0; i < nalternatives; i++) {
327        // vector defining alternatives
328        vec vec_alt = perm_matrix.get_row(i);
329
330        // full alternative or filtered
331        if( sum(vec_alt) == vec_alt.length() )  {
332            egiw_array(i) = &alter_est;
333            continue;
334        } else if( sum(vec_alt) == 0 ) {
335            egiw_array(i) = &est;
336            continue;
337        }
338
339        nalternatives_cond = (int) sum(vec_alt) + 1;
340        ivec vec_perm(0);                   // permutation vector
341
342        for(int j = 0; j < nparams; j++) {
343            position = dimV - j - 2;
344            if ( vec_alt(position) == 0 ) {
345                vec_perm.ins(j, position + 1);
346            }
347            else {
348                vec_perm.ins(0, position + 1);
349            }
350        }
351        vec_perm.ins(0, 0);
352
353        ldmat filt (est._V(), vec_perm);
354        ldmat alt (alter_est._V(), vec_perm);
355
356        mat tmpL(dimV, dimV);
357        tmpL.set_rows( 0, alt._L().get_rows(0, nalternatives_cond - 1) );
358        tmpL.set_rows( nalternatives_cond, filt._L().get_rows(nalternatives_cond, nparams) );
359
360        vec tmpD(dimV);
361        tmpD.set_subvector( 0, alt._D()(0, nalternatives_cond - 1) );
362        tmpD.set_subvector( nalternatives_cond, filt._D()(nalternatives_cond, nparams) );
363
364        ldmat tmpLD (tmpL, tmpD);
365
366        vec_perm = sort_index(vec_perm);
367        ldmat newLD (tmpLD, vec_perm);
368
369        egiw_array(i) = new egiw(1, newLD, alter_est._nu());
370    }
371
372    // Approximation
373    double sumVecCommon;                // frequently used term
374    vec vecNu(nalternatives);           // vector of nus of components
375    vec vecD(nalternatives);            // vector of LS reminders
376    vec vecCommon(nalternatives);       // vector of common parts
377    mat matVecsTheta;                           // matrix whose rows are theta vects.
378
379    for (i = 0; i < nalternatives; i++) {
380        vecNu.shift_left( egiw_array(i)->_nu() );
381        vecD.shift_left( egiw_array(i)->_V()._D()(0) );
382        matVecsTheta.append_row( egiw_array(i)->est_theta()  );
383    }
384
385    vecCommon = elem_mult ( frg, elem_div(vecNu, vecD) );
386    sumVecCommon = sum(vecCommon);
387
388    // approximation of est. regr. coefficients
389    vec aprEstTheta(nparams);  aprEstTheta.zeros();
390    for (i = 0; i < nalternatives; i++) {
391        aprEstTheta +=  matVecsTheta.get_row(i) * vecCommon(i);
392    }
393    aprEstTheta /= sumVecCommon;
394   
395    // approximation of degr. of freedom
396    double aprNu;
397    double A = log( sumVecCommon );             // Term 'A' in equation
398    for ( int i = 0; i < nalternatives; i++ ) {
399        A += frg(i) * ( log( vecD(i) ) - psi( 0.5 * vecNu(i) ) );
400    }
401    aprNu = ( 1 + sqrt(1 + 4 * (A - LOG2)/3 ) ) / ( 2 * (A - LOG2) );
402
403    // approximation of LS reminder D(0,0)
404    double aprD = aprNu / sumVecCommon;
405
406    // Aproximation of covariance of LS est.
407    mat aprC = zeros(nparams, nparams);
408    for ( int i = 0; i < nalternatives; i++ ) {
409        aprC += egiw_array(i)->est_theta_cov().to_mat() * frg(i); 
410        vec tmp = ( matVecsTheta.get_row(i) - aprEstTheta );
411        aprC += vecCommon(i) * outer_product( tmp, tmp);
412    }
413
414    // Construct GiW pdf
415    ldmat aprCinv ( inv(aprC) );
416    vec D = concat( aprD, aprCinv._D() );
417    mat L = eye(dimV);
418    L.set_submatrix(1, 0, aprCinv._L() * aprEstTheta);
419    L.set_submatrix(1, 1, aprCinv._L());
420    ldmat aprLD (L, D);
421    est = egiw(1, aprLD, aprNu);
422
423        if ( evalll ) {
424                last_lognc = est.lognc();
425        }
426               
427    // update
428    ARX::bayes ( val, cond_rgr );
429}
430
431}
Note: See TracBrowser for help on using the browser.