root/bdm/itpp_ext.cpp @ 37

Revision 37, 5.8 kB (checked in by smidl, 16 years ago)

Matrix in Cholesky decomposition, Square-root Kalman and many bug fixes

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