root/library/bdm/itpp_ext.cpp @ 869

Revision 869, 11.5 kB (checked in by sarka, 14 years ago)

assertion failure on windows eliminated

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