root/library/bdm/itpp_ext.cpp @ 1274

Revision 1200, 13.4 kB (checked in by smidl, 14 years ago)

vec ivec conversion

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