root/bdm/itpp_ext.cpp @ 86

Revision 86, 6.0 kB (checked in by smidl, 16 years ago)

linspace ITpp extenion

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/itbase.h>
14#include "itpp_ext.h"
15
16// from  algebra/lapack.h
17
18extern "C" {  /* QR factorization of a general matrix A  */
19  void dgeqrf_(int *m, int *n, double *a, int *lda, double *tau, double *work,
20               int *lwork, int *info);
21};
22
23namespace itpp {
24        Array<int> to_Arr(const ivec &indices)
25    {
26    Array<int> a(indices.size());
27    for (int i = 0; i < a.size(); i++) {
28      a(i) = indices(i);
29    }
30    return a;
31  }
32
33ivec linspace ( int from, int to ) {
34        int n=to-from+1; 
35        int i;
36        it_assert_debug ( n>0,"wrong linspace" ); 
37        ivec iv ( n ); for (i=0;i<n;i++) iv(i)=from+i;
38        return iv;
39};
40
41//Gamma
42
43  Gamma_RNG::Gamma_RNG(double a, double b)
44  {
45    setup(a,b);
46  }
47
48#define log std::log
49#define exp std::exp
50#define sqrt std::sqrt
51#define R_FINITE std::isfinite
52
53double  Gamma_RNG::sample()
54  {
55  //A copy of rgamma code from the R package!!
56  //
57 
58/* Constants : */
59    const static double sqrt32 = 5.656854;
60    const static double exp_m1 = 0.36787944117144232159;/* exp(-1) = 1/e */
61
62    /* Coefficients q[k] - for q0 = sum(q[k]*a^(-k))
63     * Coefficients a[k] - for q = q0+(t*t/2)*sum(a[k]*v^k)
64     * Coefficients e[k] - for exp(q)-1 = sum(e[k]*q^k)
65     */
66    const static double q1 = 0.04166669;
67    const static double q2 = 0.02083148;
68    const static double q3 = 0.00801191;
69    const static double q4 = 0.00144121;
70    const static double q5 = -7.388e-5;
71    const static double q6 = 2.4511e-4;
72    const static double q7 = 2.424e-4;
73
74    const static double a1 = 0.3333333;
75    const static double a2 = -0.250003;
76    const static double a3 = 0.2000062;
77    const static double a4 = -0.1662921;
78    const static double a5 = 0.1423657;
79    const static double a6 = -0.1367177;
80    const static double a7 = 0.1233795;
81
82    /* State variables [FIXME for threading!] :*/
83    static double aa = 0.;
84    static double aaa = 0.;
85    static double s, s2, d;    /* no. 1 (step 1) */
86    static double q0, b, si, c;/* no. 2 (step 4) */
87
88    double e, p, q, r, t, u, v, w, x, ret_val;
89    double a=alpha;
90    double scale=1.0/beta;
91
92    if (!R_FINITE(a) || !R_FINITE(scale) || a < 0.0 || scale <= 0.0)
93        it_error("Gamma_RNG wrong parameters");
94
95    if (a < 1.) { /* GS algorithm for parameters a < 1 */
96        if(a == 0)
97            return 0.;
98        e = 1.0 + exp_m1 * a;
99        for(;;) {  //VS repeat
100            p = e * unif_rand();
101            if (p >= 1.0) {
102                x = -log((e - p) / a);
103                if (exp_rand() >= (1.0 - a) * log(x))
104                    break;
105            } else {
106                x = exp(log(p) / a);
107                if (exp_rand() >= x)
108                    break;
109            }
110        }
111        return scale * x;
112    }
113
114    /* --- a >= 1 : GD algorithm --- */
115
116    /* Step 1: Recalculations of s2, s, d if a has changed */
117    if (a != aa) {
118        aa = a;
119        s2 = a - 0.5;
120        s = sqrt(s2);
121        d = sqrt32 - s * 12.0;
122    }
123    /* Step 2: t = standard normal deviate,
124               x = (s,1/2) -normal deviate. */
125
126    /* immediate acceptance (i) */
127    t = norm_rand();
128    x = s + 0.5 * t;
129    ret_val = x * x;
130    if (t >= 0.0)
131        return scale * ret_val;
132
133    /* Step 3: u = 0,1 - uniform sample. squeeze acceptance (s) */
134    u = unif_rand();
135    if ((d * u) <= (t * t * t))
136        return scale * ret_val;
137
138    /* Step 4: recalculations of q0, b, si, c if necessary */
139
140    if (a != aaa) {
141        aaa = a;
142        r = 1.0 / a;
143        q0 = ((((((q7 * r + q6) * r + q5) * r + q4) * r + q3) * r
144               + q2) * r + q1) * r;
145
146        /* Approximation depending on size of parameter a */
147        /* The constants in the expressions for b, si and c */
148        /* were established by numerical experiments */
149
150        if (a <= 3.686) {
151            b = 0.463 + s + 0.178 * s2;
152            si = 1.235;
153            c = 0.195 / s - 0.079 + 0.16 * s;
154        } else if (a <= 13.022) {
155            b = 1.654 + 0.0076 * s2;
156            si = 1.68 / s + 0.275;
157            c = 0.062 / s + 0.024;
158        } else {
159            b = 1.77;
160            si = 0.75;
161            c = 0.1515 / s;
162        }
163    }
164    /* Step 5: no quotient test if x not positive */
165
166    if (x > 0.0) {
167        /* Step 6: calculation of v and quotient q */
168        v = t / (s + s);
169        if (fabs(v) <= 0.25)
170            q = q0 + 0.5 * t * t * ((((((a7 * v + a6) * v + a5) * v + a4) * v
171                                      + a3) * v + a2) * v + a1) * v;
172        else
173            q = q0 - s * t + 0.25 * t * t + (s2 + s2) * log(1.0 + v);
174
175
176        /* Step 7: quotient acceptance (q) */
177        if (log(1.0 - u) <= q)
178            return scale * ret_val;
179    }
180
181    for(;;) { //VS repeat
182        /* Step 8: e = standard exponential deviate
183         *      u =  0,1 -uniform deviate
184         *      t = (b,si)-double exponential (laplace) sample */
185        e = exp_rand();
186        u = unif_rand();
187        u = u + u - 1.0;
188        if (u < 0.0)
189            t = b - si * e;
190        else
191            t = b + si * e;
192        /* Step  9:  rejection if t < tau(1) = -0.71874483771719 */
193        if (t >= -0.71874483771719) {
194            /* Step 10:  calculation of v and quotient q */
195            v = t / (s + s);
196            if (fabs(v) <= 0.25)
197                q = q0 + 0.5 * t * t *
198                    ((((((a7 * v + a6) * v + a5) * v + a4) * v + a3) * v
199                      + a2) * v + a1) * v;
200            else
201                q = q0 - s * t + 0.25 * t * t + (s2 + s2) * log(1.0 + v);
202            /* Step 11:  hat acceptance (h) */
203            /* (if q not positive go to step 8) */
204            if (q > 0.0) {
205                // TODO: w = expm1(q);
206                w = exp(q)-1;
207                /*  ^^^^^ original code had approximation with rel.err < 2e-7 */
208                /* if t is rejected sample again at step 8 */
209                if ((c * fabs(u)) <= (w * exp(e - 0.5 * t * t)))
210                    break;
211            }
212        }
213    } /* repeat .. until  `t' is accepted */
214    x = s + 0.5 * t;
215    return scale * x * x;
216}
217
218
219  bool qr(const mat &A, mat &R)
220  {
221    int info;
222    int m = A.rows();
223    int n = A.cols();
224    int lwork = n;
225    int k = std::min(m, n);
226    vec tau(k);
227    vec work(lwork);
228
229    R = A;
230
231    // perform workspace query for optimum lwork value
232    int lwork_tmp = -1;
233    dgeqrf_(&m, &n, R._data(), &m, tau._data(), work._data(), &lwork_tmp,
234            &info);
235    if (info == 0) {
236      lwork = static_cast<int>(work(0));
237      work.set_size(lwork, false);
238    }
239    dgeqrf_(&m, &n, R._data(), &m, tau._data(), work._data(), &lwork, &info);
240
241    // construct R
242    for (int i=0; i<m; i++)
243      for (int j=0; j<std::min(i,n); j++)
244        R(i,j) = 0;
245
246    return (info==0);
247  }
248
249}
Note: See TracBrowser for help on using the browser.