root/library/bdm/itpp_ext.cpp @ 679

Revision 679, 10.6 kB (checked in by smidl, 15 years ago)

Major changes in BM -- OK is only test suite and tests/tutorial -- the rest is broken!!!

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