root/library/bdm/itpp_ext.cpp @ 723

Revision 723, 11.0 kB (checked in by smidl, 15 years ago)

Big commit of LQG stuff

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