root/library/bdm/itpp_ext.cpp @ 996

Revision 996, 11.6 kB (checked in by smidl, 14 years ago)

Scheduling of forgetting factor

  • Property svn:eol-style set to native
RevLine 
[6]1//
2// C++ Implementation: itpp_ext
3//
[145]4// Description:
[6]5//
6//
7// Author: smidl <smidl@utia.cas.cz>, (C) 2008
8//
9// Copyright: See COPYING file that comes with this distribution
10//
11//
12
[32]13#include "itpp_ext.h"
[6]14
[343]15#ifndef M_PI
16#define M_PI        3.14159265358979323846
17#endif
[37]18// from  algebra/lapack.h
19
[737]20extern "C" {  /* QR factorization of a general matrix A  */
21        void dgeqrf_ ( int *m, int *n, double *a, int *lda, double *tau, double *work,
22        int *lwork, int *info );
[37]23};
24
[737]25namespace itpp {
26vec empty_vec = vec ( 0 );
27
28Array<int> to_Arr ( const ivec &indices ) {
29        Array<int> a ( indices.size() );
30        for ( int i = 0; i < a.size(); i++ ) {
31                a ( i ) = indices ( i );
[145]32        }
[477]33        return a;
34}
[6]35
[737]36ivec linspace ( int from, int to ) {
[477]37        int n = to - from + 1;
38        int i;
[737]39        it_assert_debug ( n > 0, "wrong linspace" );
40        ivec iv ( n );
41        for ( i = 0; i < n; i++ ) iv ( i ) = from + i;
[477]42        return iv;
43};
[86]44
[737]45void set_subvector ( vec &ov, const ivec &iv, const vec &v ) {
46        it_assert_debug ( ( iv.length() <= v.length() ),
[477]47                          "Vec<>::set_subvector(ivec, vec<Num_T>): Indexing out "
[737]48                          "of range of v" );
49        for ( int i = 0; i < iv.length(); i++ ) {
50                it_assert_debug ( iv ( i ) < ov.length(),
51                                  "Vec<>::set_subvector(ivec, vec<Num_T>): Indexing out "
52                                  "of range of v" );
53                ov ( iv ( i ) ) = v ( i );
[145]54        }
[477]55}
[145]56
[996]57//! return vector of elements at positions given by indexlist
[737]58vec get_vec ( const vec &v, const ivec &indexlist ) {
[477]59        int size = indexlist.size();
[737]60        vec temp ( size );
61        for ( int i = 0; i < size; ++i ) {
62                temp ( i ) = v._data() [indexlist ( i ) ];
[180]63        }
[477]64        return temp;
65}
[180]66
[276]67// Gamma
[404]68#define log std::log
69#define exp std::exp
70#define sqrt std::sqrt
71#define R_FINITE std::isfinite
[11]72
[737]73bvec operator> ( const vec &t1, const vec &t2 ) {
74        it_assert_debug ( t1.length() == t2.length(), "Vec<>::operator>(): different size of vectors" );
75        bvec temp ( t1.length() );
76        for ( int i = 0; i < t1.length(); i++ )
77                temp ( i ) = ( t1[i] > t2[i] );
[477]78        return temp;
79}
[32]80
[737]81bvec operator< ( const vec &t1, const vec &t2 ) {
82        it_assert_debug ( t1.length() == t2.length(), "Vec<>::operator>(): different size of vectors" );
83        bvec temp ( t1.length() );
84        for ( int i = 0; i < t1.length(); i++ )
85                temp ( i ) = ( t1[i] < t2[i] );
[477]86        return temp;
87}
[145]88
89
[737]90bvec operator& ( const bvec &a, const bvec &b ) {
91        it_assert_debug ( b.size() == a.size(), "operator&(): Vectors of different lengths" );
[477]92
[737]93        bvec temp ( a.size() );
94        for ( int i = 0; i < a.size(); i++ ) {
95                temp ( i ) = a ( i ) & b ( i );
[145]96        }
[477]97        return temp;
98}
[145]99
[737]100bvec operator| ( const bvec &a, const bvec &b ) {
101        it_assert_debug ( b.size() == a.size(), "operator|(): Vectors of different lengths" );
[145]102
[737]103        bvec temp ( a.size() );
104        for ( int i = 0; i < a.size(); i++ ) {
105                temp ( i ) = a ( i ) | b ( i );
[145]106        }
[477]107        return temp;
108}
[32]109
[712]110//! poor man's operator vec(bvec) - copied for svn version of itpp
[737]111ivec get_from_bvec ( const ivec &v, const bvec &binlist ) {
112        int size = binlist.size();
113        it_assert_debug ( v.size() == size, "Vec<>::operator()(bvec &): "
114                          "Wrong size of binlist vector" );
115        ivec temp ( size );
116        int j = 0;
117        for ( int i = 0; i < size; ++i )
118                if ( binlist ( i ) == bin ( 1 ) )
119                        temp ( j++ ) = v ( i );
120        temp.set_size ( j, true );
121        return temp;
[712]122}
123
124
[662]125//#if 0
[737]126Gamma_RNG::Gamma_RNG ( double a, double b ) {
127        setup ( a, b );
[404]128}
[737]129double  Gamma_RNG::sample() {
[477]130        //A copy of rgamma code from the R package!!
131        //
[32]132
[477]133        /* Constants : */
134        const static double sqrt32 = 5.656854;
135        const static double exp_m1 = 0.36787944117144232159;/* exp(-1) = 1/e */
[276]136
[477]137        /* Coefficients q[k] - for q0 = sum(q[k]*a^(-k))
138         * Coefficients a[k] - for q = q0+(t*t/2)*sum(a[k]*v^k)
139         * Coefficients e[k] - for exp(q)-1 = sum(e[k]*q^k)
140         */
141        const static double q1 = 0.04166669;
142        const static double q2 = 0.02083148;
143        const static double q3 = 0.00801191;
144        const static double q4 = 0.00144121;
145        const static double q5 = -7.388e-5;
146        const static double q6 = 2.4511e-4;
147        const static double q7 = 2.424e-4;
[276]148
[477]149        const static double a1 = 0.3333333;
150        const static double a2 = -0.250003;
151        const static double a3 = 0.2000062;
152        const static double a4 = -0.1662921;
153        const static double a5 = 0.1423657;
154        const static double a6 = -0.1367177;
155        const static double a7 = 0.1233795;
[276]156
[477]157        /* State variables [FIXME for threading!] :*/
158        static double aa = 0.;
159        static double aaa = 0.;
160        static double s, s2, d;    /* no. 1 (step 1) */
161        static double q0, b, si, c;/* no. 2 (step 4) */
[276]162
[477]163        double e, p, q, r, t, u, v, w, x, ret_val;
164        double a = alpha;
165        double scale = 1.0 / beta;
[276]166
[737]167        if ( !R_FINITE ( a ) || !R_FINITE ( scale ) || a < 0.0 || scale <= 0.0 ) {
168                it_error ( "Gamma_RNG wrong parameters" );
[477]169        }
[276]170
[737]171        if ( a < 1. ) { /* GS algorithm for parameters a < 1 */
172                if ( a == 0 )
[477]173                        return 0.;
174                e = 1.0 + exp_m1 * a;
[737]175                for ( ;; ) {  //VS repeat
[477]176                        p = e * unif_rand();
[737]177                        if ( p >= 1.0 ) {
178                                x = -log ( ( e - p ) / a );
179                                if ( exp_rand() >= ( 1.0 - a ) * log ( x ) )
[477]180                                        break;
181                        } else {
[737]182                                x = exp ( log ( p ) / a );
183                                if ( exp_rand() >= x )
[477]184                                        break;
[276]185                        }
186                }
[477]187                return scale * x;
188        }
[276]189
[477]190        /* --- a >= 1 : GD algorithm --- */
[276]191
[477]192        /* Step 1: Recalculations of s2, s, d if a has changed */
[737]193        if ( a != aa ) {
[477]194                aa = a;
195                s2 = a - 0.5;
[737]196                s = sqrt ( s2 );
[477]197                d = sqrt32 - s * 12.0;
198        }
199        /* Step 2: t = standard normal deviate,
200                   x = (s,1/2) -normal deviate. */
[276]201
[477]202        /* immediate acceptance (i) */
203        t = norm_rand();
204        x = s + 0.5 * t;
205        ret_val = x * x;
[737]206        if ( t >= 0.0 )
[477]207                return scale * ret_val;
[276]208
[477]209        /* Step 3: u = 0,1 - uniform sample. squeeze acceptance (s) */
210        u = unif_rand();
[737]211        if ( ( d * u ) <= ( t * t * t ) )
[477]212                return scale * ret_val;
[276]213
[477]214        /* Step 4: recalculations of q0, b, si, c if necessary */
[276]215
[737]216        if ( a != aaa ) {
[477]217                aaa = a;
218                r = 1.0 / a;
[737]219                q0 = ( ( ( ( ( ( q7 * r + q6 ) * r + q5 ) * r + q4 ) * r + q3 ) * r
220                         + q2 ) * r + q1 ) * r;
[276]221
[477]222                /* Approximation depending on size of parameter a */
223                /* The constants in the expressions for b, si and c */
224                /* were established by numerical experiments */
[276]225
[737]226                if ( a <= 3.686 ) {
[477]227                        b = 0.463 + s + 0.178 * s2;
228                        si = 1.235;
229                        c = 0.195 / s - 0.079 + 0.16 * s;
[737]230                } else if ( a <= 13.022 ) {
[477]231                        b = 1.654 + 0.0076 * s2;
232                        si = 1.68 / s + 0.275;
233                        c = 0.062 / s + 0.024;
234                } else {
235                        b = 1.77;
236                        si = 0.75;
237                        c = 0.1515 / s;
[276]238                }
[477]239        }
240        /* Step 5: no quotient test if x not positive */
[276]241
[737]242        if ( x > 0.0 ) {
[477]243                /* Step 6: calculation of v and quotient q */
[737]244                v = t / ( s + s );
245                if ( fabs ( v ) <= 0.25 )
246                        q = q0 + 0.5 * t * t * ( ( ( ( ( ( a7 * v + a6 ) * v + a5 ) * v + a4 ) * v
247                                                     + a3 ) * v + a2 ) * v + a1 ) * v;
[477]248                else
[737]249                        q = q0 - s * t + 0.25 * t * t + ( s2 + s2 ) * log ( 1.0 + v );
[477]250
251
252                /* Step 7: quotient acceptance (q) */
[737]253                if ( log ( 1.0 - u ) <= q )
[477]254                        return scale * ret_val;
255        }
256
[737]257        for ( ;; ) { //VS repeat
[477]258                /* Step 8: e = standard exponential deviate
259                 *      u =  0,1 -uniform deviate
260                 *      t = (b,si)-double exponential (laplace) sample */
261                e = exp_rand();
262                u = unif_rand();
263                u = u + u - 1.0;
[737]264                if ( u < 0.0 )
[477]265                        t = b - si * e;
266                else
267                        t = b + si * e;
268                /* Step  9:  rejection if t < tau(1) = -0.71874483771719 */
[737]269                if ( t >= -0.71874483771719 ) {
[477]270                        /* Step 10:      calculation of v and quotient q */
[737]271                        v = t / ( s + s );
272                        if ( fabs ( v ) <= 0.25 )
[477]273                                q = q0 + 0.5 * t * t *
[737]274                                    ( ( ( ( ( ( a7 * v + a6 ) * v + a5 ) * v + a4 ) * v + a3 ) * v
275                                        + a2 ) * v + a1 ) * v;
[276]276                        else
[737]277                                q = q0 - s * t + 0.25 * t * t + ( s2 + s2 ) * log ( 1.0 + v );
[477]278                        /* Step 11:      hat acceptance (h) */
279                        /* (if q not positive go to step 8) */
[737]280                        if ( q > 0.0 ) {
[477]281                                // TODO: w = expm1(q);
[737]282                                w = exp ( q ) - 1;
[477]283                                /*  ^^^^^ original code had approximation with rel.err < 2e-7 */
284                                /* if t is rejected sample again at step 8 */
[737]285                                if ( ( c * fabs ( u ) ) <= ( w * exp ( e - 0.5 * t * t ) ) )
[477]286                                        break;
287                        }
[276]288                }
[477]289        } /* repeat .. until  `t' is accepted */
290        x = s + 0.5 * t;
291        return scale * x * x;
292}
[276]293
294
[737]295bool qr ( const mat &A, mat &R ) {
[477]296        int info;
297        int m = A.rows();
298        int n = A.cols();
299        int lwork = n;
[737]300        int k = std::min ( m, n );
301        vec tau ( k );
302        vec work ( lwork );
[276]303
[477]304        R = A;
[37]305
[477]306        // perform workspace query for optimum lwork value
307        int lwork_tmp = -1;
[737]308        dgeqrf_ ( &m, &n, R._data(), &m, tau._data(), work._data(), &lwork_tmp,
309                  &info );
310        if ( info == 0 ) {
311                lwork = static_cast<int> ( work ( 0 ) );
312                work.set_size ( lwork, false );
[477]313        }
[737]314        dgeqrf_ ( &m, &n, R._data(), &m, tau._data(), work._data(), &lwork, &info );
[37]315
[477]316        // construct R
[737]317        for ( int i = 0; i < m; i++ )
318                for ( int j = 0; j < std::min ( i, n ); j++ )
319                        R ( i, j ) = 0;
[37]320
[737]321        return ( info == 0 );
[477]322}
[37]323
[662]324//#endif
[737]325std::string num2str ( double d ) {
[477]326        char tmp[20];//that should do
[737]327        sprintf ( tmp, "%f", d );
328        return std::string ( tmp );
[477]329};
[737]330std::string num2str ( int i ) {
[477]331        char tmp[10];//that should do
[737]332        sprintf ( tmp, "%d", i );
333        return std::string ( tmp );
[477]334};
[328]335
336// digamma
337// copied from C. Bonds' source
338#include <math.h>
339#define el 0.5772156649015329
340
[737]341double psi ( double x ) {
[477]342        double s, ps, xa, x2;
343        int n, k;
344        static double a[] = {
345                -0.8333333333333e-01,
346                0.83333333333333333e-02,
347                -0.39682539682539683e-02,
348                0.41666666666666667e-02,
349                -0.75757575757575758e-02,
350                0.21092796092796093e-01,
351                -0.83333333333333333e-01,
352                0.4432598039215686
353        };
[328]354
[737]355        xa = fabs ( x );
[477]356        s = 0.0;
[737]357        if ( ( x == ( int ) x ) && ( x <= 0.0 ) ) {
[477]358                ps = 1e308;
359                return ps;
360        }
[737]361        if ( xa == ( int ) xa ) {
[766]362                n = (int) xa;
[737]363                for ( k = 1; k < n; k++ ) {
[477]364                        s += 1.0 / k;
365                }
366                ps =  s - el;
[737]367        } else if ( ( xa + 0.5 ) == ( ( int ) ( xa + 0.5 ) ) ) {
[766]368                n = (int) (xa - 0.5);
[737]369                for ( k = 1; k <= n; k++ ) {
370                        s += 1.0 / ( 2.0 * k - 1.0 );
[477]371                }
372                ps = 2.0 * s - el - 1.386294361119891;
373        } else {
[737]374                if ( xa < 10.0 ) {
375                        n = 10 - ( int ) xa;
376                        for ( k = 0; k < n; k++ ) {
377                                s += 1.0 / ( xa + k );
[477]378                        }
379                        xa += n;
380                }
[737]381                x2 = 1.0 / ( xa * xa );
382                ps = log ( xa ) - 0.5 / xa + x2 * ( ( ( ( ( ( ( a[7] * x2 + a[6] ) * x2 + a[5] ) * x2 +
383                                                        a[4] ) * x2 + a[3] ) * x2 + a[2] ) * x2 + a[1] ) * x2 + a[0] );
[477]384                ps -= s;
385        }
[737]386        if ( x < 0.0 )
387                ps = ps - M_PI * std::cos ( M_PI * x ) / std::sin ( M_PI * x ) - 1.0 / x;
[477]388        return ps;
[32]389}
[328]390
[737]391void triu ( mat &A ) {
392        for ( int i = 1; i < A.rows(); i++ ) { // row cycle
[787]393                for ( int j = 0; (j < i) && (j<A.cols()); j++ ) {
[737]394                        A ( i, j ) = 0;
395                }
[508]396        }
[328]397}
[579]398
[661]399//! Storage of randun() internals
[737]400class RandunStorage {
401        const int A;
402        const int M;
403        static double seed;
404        static int counter;
405public:
406        RandunStorage() : A ( 16807 ), M ( 2147483647 ) {};
407        //!set seed of the randun() generator
408        void set_seed ( double seed0 ) {
409                seed = seed0;
410        }
411        //! generate randun() sample
412        double get() {
[869]413                long long tmp = (long long) (A * seed);
[737]414                tmp = tmp % M;
[942]415                seed = (double) tmp;
[737]416                counter++;
417                return seed / M;
418        }
[579]419};
420static RandunStorage randun_global_storage;
[622]421double RandunStorage::seed = 1111111;
422int RandunStorage::counter = 0;
[737]423double randun() {
424        return randun_global_storage.get();
425};
426vec randun ( int n ) {
427        vec res ( n );
428        for ( int i = 0; i < n; i++ ) {
429                res ( i ) = randun();
430        };
431        return res;
432};
433mat randun ( int n, int m ) {
434        mat res ( n, m );
435        for ( int i = 0; i < n*m; i++ ) {
436                res ( i ) = randun();
437        };
438        return res;
439};
[622]440
[737]441ivec unique ( const ivec &in ) {
442        ivec uniq ( 0 );
[584]443        int j = 0;
444        bool found = false;
[737]445        for ( int i = 0; i < in.length(); i++ ) {
[584]446                found = false;
447                j = 0;
[737]448                while ( ( !found ) && ( j < uniq.length() ) ) {
449                        if ( in ( i ) == uniq ( j ) ) found = true;
[584]450                        j++;
451                }
[737]452                if ( !found ) uniq = concat ( uniq, in ( i ) );
[584]453        }
454        return uniq;
455}
[579]456
[737]457ivec unique_complement ( const ivec &in, const ivec &base ) {
[586]458        // almost a copy of unique
[737]459        ivec uniq ( 0 );
[586]460        int j = 0;
461        bool found = false;
[737]462        for ( int i = 0; i < in.length(); i++ ) {
[586]463                found = false;
464                j = 0;
[737]465                while ( ( !found ) && ( j < uniq.length() ) ) {
466                        if ( in ( i ) == uniq ( j ) ) found = true;
[586]467                        j++;
468                }
[737]469                j = 0;
470                while ( ( !found ) && ( j < base.length() ) ) {
471                        if ( in ( i ) == base ( j ) ) found = true;
[586]472                        j++;
473                }
[737]474                if ( !found ) uniq = concat ( uniq, in ( i ) );
[586]475        }
476        return uniq;
[508]477}
[586]478
479}
Note: See TracBrowser for help on using the browser.