root/bdm/itpp_ext.cpp @ 145

Revision 145, 7.1 kB (checked in by smidl, 16 years ago)

Oprava dokumentace

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