[0b990d] | 1 |
|
---|
| 2 | /*
|
---|
| 3 | * These routines are based on the work of Edward T. Seidl at the
|
---|
| 4 | * National Institutes of Health.
|
---|
| 5 | */
|
---|
| 6 |
|
---|
| 7 | #include <stdio.h>
|
---|
| 8 | #include <stdlib.h>
|
---|
| 9 | #include <math.h>
|
---|
| 10 | #include <string.h>
|
---|
| 11 | #include <math/scmat/cmatrix.h>
|
---|
| 12 |
|
---|
| 13 | static void ludcmp(double**, int, int*, double*);
|
---|
| 14 | static void lubksb(double**, int, int*, double*);
|
---|
| 15 | static void symm_lu_decomp(double**, int, double*);
|
---|
| 16 | static void symm_lu_back_sub(double**, int, double*);
|
---|
| 17 |
|
---|
| 18 | static void tred2(int dim,double**,double*,double*,int);
|
---|
| 19 | static void tqli(int dim,double*,double**,double*,int,double);
|
---|
| 20 | static void eigsort(int dim,double*,double**);
|
---|
| 21 |
|
---|
| 22 | double**
|
---|
| 23 | cmat_new_square_matrix(int n)
|
---|
| 24 | {
|
---|
| 25 | double *mat;
|
---|
| 26 | double **r;
|
---|
| 27 | if (n == 0) return 0;
|
---|
| 28 | mat = (double*) malloc(sizeof(double)*n*n);
|
---|
| 29 | if (!mat) return 0;
|
---|
| 30 | r = (double**) malloc(sizeof(double*)*n);
|
---|
| 31 | if (!r) {
|
---|
| 32 | free(mat);
|
---|
| 33 | return 0;
|
---|
| 34 | }
|
---|
| 35 | cmat_matrix_pointers(r,mat,n,n);
|
---|
| 36 | return r;
|
---|
| 37 | }
|
---|
| 38 |
|
---|
| 39 | double**
|
---|
| 40 | cmat_new_rect_matrix(int n,int m)
|
---|
| 41 | {
|
---|
| 42 | double *mat;
|
---|
| 43 | double **r;
|
---|
| 44 | if (n == 0 || m == 0) return 0;
|
---|
| 45 | mat = (double*) malloc(sizeof(double)*n*m);
|
---|
| 46 | if (!mat) return 0;
|
---|
| 47 | r = (double**) malloc(sizeof(double*)*n);
|
---|
| 48 | if (!r) {
|
---|
| 49 | free(mat);
|
---|
| 50 | return 0;
|
---|
| 51 | }
|
---|
| 52 | cmat_matrix_pointers(r,mat,n,m);
|
---|
| 53 | return r;
|
---|
| 54 | }
|
---|
| 55 |
|
---|
| 56 | /* this deletes both square and triangular matrices */
|
---|
| 57 | void
|
---|
| 58 | cmat_delete_matrix(double**m)
|
---|
| 59 | {
|
---|
| 60 | if (m) {
|
---|
| 61 | free(m[0]);
|
---|
| 62 | free(m);
|
---|
| 63 | }
|
---|
| 64 | }
|
---|
| 65 |
|
---|
| 66 | void
|
---|
| 67 | cmat_transpose_square_matrix(double**matrix, int n)
|
---|
| 68 | {
|
---|
| 69 | int i,j;
|
---|
| 70 | for (i=0; i<n; i++) {
|
---|
| 71 | for (j=0; j<i; j++) {
|
---|
| 72 | double tmp = matrix[i][j];
|
---|
| 73 | matrix[i][j] = matrix[j][i];
|
---|
| 74 | matrix[j][i] = tmp;
|
---|
| 75 | }
|
---|
| 76 | }
|
---|
| 77 | }
|
---|
| 78 |
|
---|
| 79 | void
|
---|
| 80 | cmat_matrix_pointers(double**ptrs,double*matrix,int nrow, int ncol)
|
---|
| 81 | {
|
---|
| 82 | int i;
|
---|
| 83 | for (i=0; i<nrow; i++) ptrs[i] = &matrix[i*ncol];
|
---|
| 84 | }
|
---|
| 85 |
|
---|
| 86 | /*
|
---|
| 87 | * a contains pointers to the an area of contiguous storage.
|
---|
| 88 | * Its dimensions are nr by nc. On exit it will be transposed,
|
---|
| 89 | * however the a vector of double* is itself unchanged. Another
|
---|
| 90 | * vector is needed to access the storage or a must be updated
|
---|
| 91 | * after this routine is called.
|
---|
| 92 | */
|
---|
| 93 | void
|
---|
| 94 | cmat_transpose_matrix(double**a, int nr, int nc)
|
---|
| 95 | {
|
---|
| 96 | int i,j;
|
---|
| 97 | double* tmpp;
|
---|
| 98 | double* tmp;
|
---|
| 99 |
|
---|
| 100 | if (nr == 0 || nc == 0) return;
|
---|
| 101 |
|
---|
| 102 | if (nr == nc) {
|
---|
| 103 | cmat_transpose_square_matrix(a,nr);
|
---|
| 104 | return;
|
---|
| 105 | };
|
---|
| 106 |
|
---|
| 107 | tmp = (double*) malloc(sizeof(double)*nr*nc);
|
---|
| 108 | if (!tmp && nr && nc) {
|
---|
| 109 | fprintf(stderr,"cmat_transpose_matrix: malloc failed\n");
|
---|
| 110 | abort();
|
---|
| 111 | }
|
---|
| 112 |
|
---|
| 113 | tmpp = tmp;
|
---|
| 114 | for (i=0; i<nc; i++) {
|
---|
| 115 | for (j=0; j<nr; j++) {
|
---|
| 116 | *tmpp = a[j][i];
|
---|
| 117 | tmpp++;
|
---|
| 118 | }
|
---|
| 119 | }
|
---|
| 120 |
|
---|
| 121 | memcpy(a[0],tmp,sizeof(double)*nr*nc);
|
---|
| 122 |
|
---|
| 123 | if (tmp) free(tmp);
|
---|
| 124 | }
|
---|
| 125 |
|
---|
| 126 | /* a is symmetric if sym is true */
|
---|
| 127 | double
|
---|
| 128 | cmat_determ(double** a, int sym, int dim)
|
---|
| 129 | {
|
---|
| 130 | int i;
|
---|
| 131 | double det=0;
|
---|
| 132 |
|
---|
| 133 | if (sym) {
|
---|
| 134 | symm_lu_decomp(a,dim,&det);
|
---|
| 135 | } else {
|
---|
| 136 | int *indx= (int*) malloc(sizeof(int)*dim);
|
---|
| 137 | ludcmp(a,dim,indx,&det);
|
---|
| 138 | free(indx);
|
---|
| 139 | }
|
---|
| 140 |
|
---|
| 141 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 142 |
|
---|
| 143 | for (i=0; i < dim; i++) det *= a[i][i];
|
---|
| 144 |
|
---|
| 145 | return det;
|
---|
| 146 | }
|
---|
| 147 |
|
---|
| 148 | /* a is symmetric if sym is true */
|
---|
| 149 | double
|
---|
| 150 | cmat_solve_lin(double** a, int sym, double* b, int dim)
|
---|
| 151 | {
|
---|
| 152 | int i;
|
---|
| 153 | double det=0;
|
---|
| 154 |
|
---|
| 155 | if (sym) {
|
---|
| 156 | symm_lu_decomp(a,dim,&det);
|
---|
| 157 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 158 | symm_lu_back_sub(a,dim,b);
|
---|
| 159 | } else {
|
---|
| 160 | int *indx= (int*) malloc(sizeof(int)*dim);
|
---|
| 161 | ludcmp(a,dim,indx,&det);
|
---|
| 162 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 163 | lubksb(a,dim,indx,b);
|
---|
| 164 | free(indx);
|
---|
| 165 | }
|
---|
| 166 |
|
---|
| 167 | for(i=0; i < dim; i++) det *= a[i][i];
|
---|
| 168 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 169 |
|
---|
| 170 | return det;
|
---|
| 171 | }
|
---|
| 172 |
|
---|
| 173 | double
|
---|
| 174 | cmat_invert(double**a, int sym, int dim)
|
---|
| 175 | {
|
---|
| 176 | int i,j;
|
---|
| 177 | double det=0;
|
---|
| 178 | double **y;
|
---|
| 179 | double *b;
|
---|
| 180 |
|
---|
| 181 | b = (double*) malloc(sizeof(double)*dim);
|
---|
| 182 | y = cmat_new_square_matrix(dim);
|
---|
| 183 |
|
---|
| 184 | if (sym) {
|
---|
| 185 | symm_lu_decomp(a,dim,&det);
|
---|
| 186 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 187 |
|
---|
| 188 | for (i=0; i < dim; i++) det *= a[i][i];
|
---|
| 189 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 190 |
|
---|
| 191 | for (i=0; i < dim; i++) {
|
---|
| 192 | for (j=0; j < dim; j++) b[j]=0;
|
---|
| 193 | b[i]=1;
|
---|
| 194 | symm_lu_back_sub(a,dim,b);
|
---|
| 195 | for (j=0; j < dim; j++) y[j][i]=b[j];
|
---|
| 196 | }
|
---|
| 197 |
|
---|
| 198 | for (i=0; i < dim; i++)
|
---|
| 199 | for (j=0; j <= i; j++)
|
---|
| 200 | a[i][j] = y[i][j];
|
---|
| 201 |
|
---|
| 202 | } else {
|
---|
| 203 | int *indx= (int*) malloc(sizeof(int)*dim);
|
---|
| 204 |
|
---|
| 205 | ludcmp(a,dim,indx,&det);
|
---|
| 206 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 207 |
|
---|
| 208 | for (i=0; i < dim; i++) det *= a[i][i];
|
---|
| 209 | if (fabs(det) < 1.0e-16) return 0;
|
---|
| 210 |
|
---|
| 211 | for (i=0; i < dim; i++) {
|
---|
| 212 | memset(b,0,sizeof(double)*dim);
|
---|
| 213 | b[i]=1;
|
---|
| 214 | lubksb(a,dim,indx,b);
|
---|
| 215 | for (j=0; j < dim; j++) y[j][i]=b[j];
|
---|
| 216 | }
|
---|
| 217 |
|
---|
| 218 | for (i=0; i < dim; i++)
|
---|
| 219 | for (j=0; j < dim; j++)
|
---|
| 220 | a[i][j] = y[i][j];
|
---|
| 221 | free(indx);
|
---|
| 222 | }
|
---|
| 223 |
|
---|
| 224 | free(b);
|
---|
| 225 | cmat_delete_matrix(y);
|
---|
| 226 |
|
---|
| 227 | return det;
|
---|
| 228 | }
|
---|
| 229 |
|
---|
| 230 | static void
|
---|
| 231 | ludcmp(double** a, int n, int *indx, double *d)
|
---|
| 232 | {
|
---|
| 233 | int i,imax=0,j,k;
|
---|
| 234 | double big,dum,sum,temp;
|
---|
| 235 |
|
---|
| 236 | double* vv = (double*) malloc(sizeof(double)*n);
|
---|
| 237 |
|
---|
| 238 | *d = 1.0;
|
---|
| 239 |
|
---|
| 240 | for (i=0; i < n ; i++) {
|
---|
| 241 | big=0.0;
|
---|
| 242 | for (j=0; j < n; j++) if ((temp=fabs(a[i][j])) > big) big=temp;
|
---|
| 243 | #if 1
|
---|
| 244 | if (big == 0.0) {
|
---|
| 245 | *d = 0.0;
|
---|
| 246 | free(vv);
|
---|
| 247 | return;
|
---|
| 248 | }
|
---|
| 249 | #else
|
---|
| 250 | if(big==0.0) big=1.0e-16;
|
---|
| 251 | #endif
|
---|
| 252 | vv[i] = 1.0/big;
|
---|
| 253 | }
|
---|
| 254 |
|
---|
| 255 | for (j=0; j < n ; j++) {
|
---|
| 256 | for (i=0; i < j ; i++) {
|
---|
| 257 | sum = a[i][j];
|
---|
| 258 | for (k=0; k < i ; k++) sum -= a[i][k]*a[k][j];
|
---|
| 259 | a[i][j] = sum;
|
---|
| 260 | }
|
---|
| 261 |
|
---|
| 262 | big = 0.0;
|
---|
| 263 | for (i=j ; i < n ; i++) {
|
---|
| 264 | sum=a[i][j];
|
---|
| 265 | for (k=0; k < j ; k++) sum -= a[i][k]*a[k][j];
|
---|
| 266 | a[i][j] = sum;
|
---|
| 267 | if ((dum=vv[i]*fabs(sum)) >= big) {
|
---|
| 268 | big = dum;
|
---|
| 269 | imax = i;
|
---|
| 270 | }
|
---|
| 271 | }
|
---|
| 272 |
|
---|
| 273 | if (j != imax) {
|
---|
| 274 | for (k=0; k < n; k++) {
|
---|
| 275 | dum=a[imax][k];
|
---|
| 276 | a[imax][k]=a[j][k];
|
---|
| 277 | a[j][k]=dum;
|
---|
| 278 | }
|
---|
| 279 | *d = -(*d);
|
---|
| 280 | vv[imax]=vv[j];
|
---|
| 281 | }
|
---|
| 282 |
|
---|
| 283 | indx[j]=imax;
|
---|
| 284 | if (a[j][j] == 0.0) a[j][j] = 1.0e-20;
|
---|
| 285 | if (j != n-1) {
|
---|
| 286 | dum = 1.0/a[j][j];
|
---|
| 287 | for (i=j+1; i < n ; i++) a[i][j] *= dum;
|
---|
| 288 | }
|
---|
| 289 | }
|
---|
| 290 | free(vv);
|
---|
| 291 | }
|
---|
| 292 |
|
---|
| 293 | static void
|
---|
| 294 | lubksb(double** a, int n, int *indx, double* b)
|
---|
| 295 | {
|
---|
| 296 | int i,ii=0,ip,j;
|
---|
| 297 | int t=0;
|
---|
| 298 | double sum;
|
---|
| 299 |
|
---|
| 300 | for (i=0; i < n ; i++) {
|
---|
| 301 | ip = indx[i];
|
---|
| 302 | sum = b[ip];
|
---|
| 303 | b[ip]=b[i];
|
---|
| 304 |
|
---|
| 305 | if(t) {
|
---|
| 306 | for (j=ii; j <= i-1 ; j++)
|
---|
| 307 | sum -= a[i][j]*b[j];
|
---|
| 308 | }
|
---|
| 309 | else if(sum) {
|
---|
| 310 | ii=i;
|
---|
| 311 | t++;
|
---|
| 312 | }
|
---|
| 313 |
|
---|
| 314 | b[i]=sum;
|
---|
| 315 | }
|
---|
| 316 |
|
---|
| 317 | for (i=n-1; i >= 0 ; i--) {
|
---|
| 318 | sum = b[i];
|
---|
| 319 | for (j=i+1; j < n ; j++) sum -= a[i][j]*b[j];
|
---|
| 320 | b[i] = sum/a[i][i];
|
---|
| 321 | }
|
---|
| 322 | }
|
---|
| 323 |
|
---|
| 324 | /*
|
---|
| 325 | * this is LU decomposition where A is a symmetric matrix
|
---|
| 326 | * when A is symmetric, then
|
---|
| 327 | * beta(i,j) = A(i,j) - sum_k(i-1) beta(k,i)*beta(k,j)/beta(k,k)
|
---|
| 328 | * alpha(i,j) = beta(j,i)/beta(j,j)
|
---|
| 329 | *
|
---|
| 330 | * since we're storing beta in a, the indices of beta will be switched
|
---|
| 331 | * since alpha is expressed in terms of beta, we don't store it
|
---|
| 332 | *
|
---|
| 333 | * so we have
|
---|
| 334 | * beta(i,j) = A(i,j) - sum_k(i-1) beta(i,k)*beta(j,k)/beta(k,k)
|
---|
| 335 | * alpha(i,j) = beta(i,j)/beta(j,j)
|
---|
| 336 | */
|
---|
| 337 |
|
---|
| 338 | static void
|
---|
| 339 | symm_lu_decomp(double** a, int n, double *d)
|
---|
| 340 | {
|
---|
| 341 | int i,j,k;
|
---|
| 342 | double tmp;
|
---|
| 343 |
|
---|
| 344 | double* v = (double*) malloc(sizeof(double)*n);
|
---|
| 345 | memset(v,0,sizeof(double)*n);
|
---|
| 346 |
|
---|
| 347 | /* check for singular matrix */
|
---|
| 348 | for (i=0; i < n; i++) {
|
---|
| 349 | for (j=0; j < i; j++) {
|
---|
| 350 | v[i] = ((tmp=fabs(a[i][j])) > v[i]) ? tmp : v[i];
|
---|
| 351 | v[j] = (tmp > v[j]) ? tmp : v[j];
|
---|
| 352 | }
|
---|
| 353 | v[i] = ((tmp=fabs(a[i][i])) > v[i]) ? tmp : v[i];
|
---|
| 354 | }
|
---|
| 355 |
|
---|
| 356 | for (i=0; i < n; i++) {
|
---|
| 357 | if (fabs(v[i]) < 1.0e-16) {
|
---|
| 358 | fprintf(stderr,"\n warning: singular matrix in symm_lu_decomp\n");
|
---|
| 359 | *d = 0.0;
|
---|
| 360 | return;
|
---|
| 361 | }
|
---|
| 362 | }
|
---|
| 363 |
|
---|
| 364 | free(v);
|
---|
| 365 |
|
---|
| 366 | *d = 1.0;
|
---|
| 367 |
|
---|
| 368 | for (i=0; i < n ; i++) {
|
---|
| 369 | /* check to make sure we're not going to blow up */
|
---|
| 370 | if (i < n-1) {
|
---|
| 371 | tmp = 0;
|
---|
| 372 | for (k=0; k < i-1; k++)
|
---|
| 373 | tmp += a[i][k]*a[i][k]/a[k][k];
|
---|
| 374 | if (fabs(a[i][i]-tmp) < 1.0e-16) {
|
---|
| 375 | fprintf(stderr,"\n warning: singular matrix in symm_lu_decomp 2\n");
|
---|
| 376 | *d = 0;
|
---|
| 377 | return;
|
---|
| 378 | }
|
---|
| 379 | }
|
---|
| 380 | for (j=i; j < n; j++) {
|
---|
| 381 | tmp = 0;
|
---|
| 382 | for (k=0; k <= i-1; k++)
|
---|
| 383 | tmp -= a[i][k]*a[j][k]/a[k][k];
|
---|
| 384 | a[j][i] += tmp;
|
---|
| 385 | }
|
---|
| 386 | }
|
---|
| 387 | }
|
---|
| 388 |
|
---|
| 389 | static void
|
---|
| 390 | symm_lu_back_sub(double** a, int n, double* b)
|
---|
| 391 | {
|
---|
| 392 | int i,j;
|
---|
| 393 | double sum;
|
---|
| 394 |
|
---|
| 395 | /* form y(i) = bi - sum_j(i-1) alpha(i,j)*y(j)
|
---|
| 396 | * alpha(i,j) = beta(j,i)/beta(j,j), but beta is stored lower instead of
|
---|
| 397 | * upper triangle, so alpha(i,j) = beta(i,j)/beta(j,j)
|
---|
| 398 | */
|
---|
| 399 | for (i=0; i < n ; i++) {
|
---|
| 400 | sum = 0;
|
---|
| 401 | for (j=0; j < i; j++)
|
---|
| 402 | sum += (a[i][j]/a[j][j]) * b[j];
|
---|
| 403 | b[i] -= sum;
|
---|
| 404 | }
|
---|
| 405 |
|
---|
| 406 | /* now form x(i) = 1/beta(i,i)*[y(i) - sum_j=i+1(N) beta(i,j)*x(j)]
|
---|
| 407 | * is really ...[...beta(j,i)*x(j)]
|
---|
| 408 | */
|
---|
| 409 | for (i=n-1; i >= 0 ; i--) {
|
---|
| 410 | sum = b[i];
|
---|
| 411 | for (j=i+1; j < n ; j++) sum -= a[j][i]*b[j];
|
---|
| 412 | b[i] = sum/a[i][i];
|
---|
| 413 | }
|
---|
| 414 | }
|
---|
| 415 |
|
---|
| 416 | /*
|
---|
| 417 | * This does c(t) (+)= a(t) * b(t), where the (t) means the transpose
|
---|
| 418 | * of the matrix can be optionally used and the (+) means that accumulation
|
---|
| 419 | * is optional. The dimensions of the matrices is as follows:
|
---|
| 420 | * a(nr,nl) (if ta then a(nl,nr))
|
---|
| 421 | * b(nl,nc) (if tb then b(nc,nl))
|
---|
| 422 | * c(nr,nc) (if tc then c(nc,nr))
|
---|
| 423 | */
|
---|
| 424 | void
|
---|
| 425 | cmat_mxm(double** a, int ta, double** b, int tb, double** c, int tc,
|
---|
| 426 | int nr, int nl, int nc, int add)
|
---|
| 427 | {
|
---|
| 428 | int odd_nr,odd_nc;
|
---|
| 429 | int i,j,k;
|
---|
| 430 | double t00,t01,t10,t11;
|
---|
| 431 | double *att,*bt;
|
---|
| 432 | double *at1,*bt1;
|
---|
| 433 | double** old_a = 0;
|
---|
| 434 | double** old_b = 0;
|
---|
| 435 |
|
---|
| 436 | odd_nr = (nr)%2;
|
---|
| 437 | odd_nc = (nc)%2;
|
---|
| 438 |
|
---|
| 439 | if(ta) {
|
---|
| 440 | cmat_transpose_matrix(a,nl,nr);
|
---|
| 441 | if (nr > nl) {
|
---|
| 442 | old_a = a;
|
---|
| 443 | a = (double**) malloc(nr*sizeof(double*));
|
---|
| 444 | if (!a) {
|
---|
| 445 | fprintf(stderr,"cmat_mxm: malloc a failed\n");
|
---|
| 446 | abort();
|
---|
| 447 | }
|
---|
| 448 | a[0] = old_a[0];
|
---|
| 449 | }
|
---|
| 450 | cmat_matrix_pointers(a,a[0],nr,nl);
|
---|
| 451 | }
|
---|
| 452 | if(!tb) {
|
---|
| 453 | cmat_transpose_matrix(b,nl,nc);
|
---|
| 454 | if (nc > nl) {
|
---|
| 455 | old_b = b;
|
---|
| 456 | b = (double**) malloc(nc*sizeof(double*));
|
---|
| 457 | if (!b) {
|
---|
| 458 | fprintf(stderr,"cmat_mxm: malloc b failed\n");
|
---|
| 459 | abort();
|
---|
| 460 | }
|
---|
| 461 | b[0] = old_b[0];
|
---|
| 462 | }
|
---|
| 463 | cmat_matrix_pointers(b,b[0],nc,nl);
|
---|
| 464 | }
|
---|
| 465 |
|
---|
| 466 | for(j=0; j < nc-1 ; j+=2) {
|
---|
| 467 | for(i=0; i < nr-1 ; i+=2) {
|
---|
| 468 | att=a[i]; bt=b[j];
|
---|
| 469 | at1=a[i+1]; bt1=b[j+1];
|
---|
| 470 | if(add) {
|
---|
| 471 | if(tc) {
|
---|
| 472 | t00 = c[j][i];
|
---|
| 473 | t01 = c[j+1][i];
|
---|
| 474 | t10 = c[j][i+1];
|
---|
| 475 | t11 = c[j+1][i+1];
|
---|
| 476 | }
|
---|
| 477 | else {
|
---|
| 478 | t00 = c[i][j];
|
---|
| 479 | t01 = c[i][j+1];
|
---|
| 480 | t10 = c[i+1][j];
|
---|
| 481 | t11 = c[i+1][j+1];
|
---|
| 482 | }
|
---|
| 483 | }
|
---|
| 484 | else
|
---|
| 485 | t00=t01=t10=t11=0.0;
|
---|
| 486 | for(k=nl; k ; k--,att++,bt++,at1++,bt1++) {
|
---|
| 487 | t00 += *att * *bt;
|
---|
| 488 | t01 += *att * *bt1;
|
---|
| 489 | t10 += *at1 * *bt;
|
---|
| 490 | t11 += *at1 * *bt1;
|
---|
| 491 | }
|
---|
| 492 | if(tc) {
|
---|
| 493 | c[j][i]=t00;
|
---|
| 494 | c[j+1][i]=t01;
|
---|
| 495 | c[j][i+1]=t10;
|
---|
| 496 | c[j+1][i+1]=t11;
|
---|
| 497 | }
|
---|
| 498 | else {
|
---|
| 499 | c[i][j]=t00;
|
---|
| 500 | c[i][j+1]=t01;
|
---|
| 501 | c[i+1][j]=t10;
|
---|
| 502 | c[i+1][j+1]=t11;
|
---|
| 503 | }
|
---|
| 504 | }
|
---|
| 505 | if(odd_nr) {
|
---|
| 506 | att=a[i]; bt=b[j];
|
---|
| 507 | bt1=b[j+1];
|
---|
| 508 | if(add) {
|
---|
| 509 | if(tc) {
|
---|
| 510 | t00 = c[j][i];
|
---|
| 511 | t01 = c[j+1][i];
|
---|
| 512 | }
|
---|
| 513 | else {
|
---|
| 514 | t00 = c[i][j];
|
---|
| 515 | t01 = c[i][j+1];
|
---|
| 516 | }
|
---|
| 517 | }
|
---|
| 518 | else t00=t01=0.0;
|
---|
| 519 | for(k= nl; k ; k--,att++,bt++,bt1++) {
|
---|
| 520 | t00 += *att * *bt;
|
---|
| 521 | t01 += *att * *bt1;
|
---|
| 522 | }
|
---|
| 523 | if(tc) {
|
---|
| 524 | c[j][i]=t00;
|
---|
| 525 | c[j+1][i]=t01;
|
---|
| 526 | }
|
---|
| 527 | else {
|
---|
| 528 | c[i][j]=t00;
|
---|
| 529 | c[i][j+1]=t01;
|
---|
| 530 | }
|
---|
| 531 | }
|
---|
| 532 | }
|
---|
| 533 | if(odd_nc) {
|
---|
| 534 | for(i=0; i < nr-1 ; i+=2) {
|
---|
| 535 | att=a[i]; bt=b[j];
|
---|
| 536 | at1=a[i+1];
|
---|
| 537 | if(add) {
|
---|
| 538 | if(tc) {
|
---|
| 539 | t00 = c[j][i];
|
---|
| 540 | t10 = c[j][i+1];
|
---|
| 541 | }
|
---|
| 542 | else {
|
---|
| 543 | t00 = c[i][j];
|
---|
| 544 | t10 = c[i+1][j];
|
---|
| 545 | }
|
---|
| 546 | }
|
---|
| 547 | else t00=t10=0.0;
|
---|
| 548 | for(k= nl; k ; k--,att++,bt++,at1++) {
|
---|
| 549 | t00 += *att * *bt;
|
---|
| 550 | t10 += *at1 * *bt;
|
---|
| 551 | }
|
---|
| 552 | if(tc) {
|
---|
| 553 | c[j][i]=t00;
|
---|
| 554 | c[j][i+1]=t10;
|
---|
| 555 | }
|
---|
| 556 | else {
|
---|
| 557 | c[i][j]=t00;
|
---|
| 558 | c[i+1][j]=t10;
|
---|
| 559 | }
|
---|
| 560 | }
|
---|
| 561 | if(odd_nr) {
|
---|
| 562 | att=a[i]; bt=b[j];
|
---|
| 563 | if(add) t00 = (tc) ? c[j][i] : c[i][j];
|
---|
| 564 | else t00=0.0;
|
---|
| 565 | for(k=nl; k ; k--,att++,bt++) t00 += *att * *bt;
|
---|
| 566 | if(tc) c[j][i]=t00;
|
---|
| 567 | else c[i][j]=t00;
|
---|
| 568 | }
|
---|
| 569 | }
|
---|
| 570 |
|
---|
| 571 | if(ta) {
|
---|
| 572 | cmat_transpose_matrix(a,nr,nl);
|
---|
| 573 | if (old_a) {
|
---|
| 574 | free(a);
|
---|
| 575 | a = old_a;
|
---|
| 576 | }
|
---|
| 577 | cmat_matrix_pointers(a,a[0],nr,nl);
|
---|
| 578 | }
|
---|
| 579 | if(!tb) {
|
---|
| 580 | cmat_transpose_matrix(b,nc,nl);
|
---|
| 581 | if (old_b) {
|
---|
| 582 | free(b);
|
---|
| 583 | b = old_b;
|
---|
| 584 | }
|
---|
| 585 | cmat_matrix_pointers(b,b[0],nl,nc);
|
---|
| 586 | }
|
---|
| 587 | }
|
---|
| 588 |
|
---|
| 589 | /*
|
---|
| 590 | * a is symmetric (na,na) in a triangular storage format
|
---|
| 591 | * b is rectangular (na,nb)
|
---|
| 592 | * a (+)= b * transpose(b) (+= if add)
|
---|
| 593 | */
|
---|
| 594 | void
|
---|
| 595 | cmat_symmetric_mxm(double**a,int na, /* a is (na,na) */
|
---|
| 596 | double**b,int nb, /* b is (na,nb) */
|
---|
| 597 | int add)
|
---|
| 598 | {
|
---|
| 599 | int i,j,k;
|
---|
| 600 | for (i=0; i<na; i++) {
|
---|
| 601 | double*ai=a[i];
|
---|
| 602 | for (j=0; j<=i; j++) {
|
---|
| 603 | double*bi=b[i];
|
---|
| 604 | double*bj=b[j];
|
---|
| 605 | double tmp;
|
---|
| 606 | if (add) tmp = 0.0;
|
---|
| 607 | else tmp = ai[j];
|
---|
| 608 | for (k=nb; k; k--,bi++,bj++) {
|
---|
| 609 | tmp += *bi * *bj;
|
---|
| 610 | }
|
---|
| 611 | ai[j] = tmp;
|
---|
| 612 | }
|
---|
| 613 | }
|
---|
| 614 | }
|
---|
| 615 |
|
---|
| 616 | /*
|
---|
| 617 | * a is symmetric (na,na) in a triangular storage format
|
---|
| 618 | * b is symmetric (nb,nb) in a triangular storage format
|
---|
| 619 | * a (+)= c * b * transpose(c) (+= if add)
|
---|
| 620 | */
|
---|
| 621 | void
|
---|
| 622 | cmat_transform_symmetric_matrix(double**a,int na, /* a is (na,na) */
|
---|
| 623 | double**b,int nb, /* b is (nb,nb) */
|
---|
| 624 | double**c, /* c is (na,nb) */
|
---|
| 625 | int add)
|
---|
| 626 | {
|
---|
| 627 | int i,j,k;
|
---|
| 628 | double**t;
|
---|
| 629 | double* brow;
|
---|
| 630 |
|
---|
| 631 | /* create a temporary matrix, t */
|
---|
| 632 | t = cmat_new_rect_matrix(na,nb);
|
---|
| 633 |
|
---|
| 634 | /* t = transpose(b * transpose(c)) */
|
---|
| 635 | brow = (double*) malloc(sizeof(double)*nb);
|
---|
| 636 | if (!brow) {
|
---|
| 637 | fprintf(stderr,"cmat_transform_symmetric_matrix: malloc brow failed\n");
|
---|
| 638 | abort();
|
---|
| 639 | }
|
---|
| 640 | for (i=0; i<nb; i++) {
|
---|
| 641 | for (k=0; k<=i; k++) brow[k] = b[i][k];
|
---|
| 642 | for ( ; k<nb; k++) brow[k] = b[k][i];
|
---|
| 643 | for (j=0; j<na; j++) {
|
---|
| 644 | double*bi = brow;
|
---|
| 645 | double*cj = c[j];
|
---|
| 646 | double tmp = 0.0;
|
---|
| 647 | for (k=nb; k; k--,bi++,cj++) tmp += *bi * *cj;
|
---|
| 648 | t[j][i] = tmp;
|
---|
| 649 | }
|
---|
| 650 | }
|
---|
| 651 | free(brow);
|
---|
| 652 |
|
---|
| 653 | /* a = c * transpose(t) */
|
---|
| 654 | for (i=0; i<na; i++) {
|
---|
| 655 | for (j=0; j<=i; j++) {
|
---|
| 656 | double*ci = c[i];
|
---|
| 657 | double*tj = t[j];
|
---|
| 658 | double tmp;
|
---|
| 659 | if (add) tmp = a[i][j];
|
---|
| 660 | else tmp = 0.0;
|
---|
| 661 | for (k=nb; k; k--,ci++,tj++) tmp += *ci * *tj;
|
---|
| 662 | a[i][j] = tmp;
|
---|
| 663 | }
|
---|
| 664 | }
|
---|
| 665 |
|
---|
| 666 | /* delete the temporary */
|
---|
| 667 | cmat_delete_matrix(t);
|
---|
| 668 | }
|
---|
| 669 |
|
---|
| 670 | /*
|
---|
| 671 | * a is symmetric (na,na) in a triangular storage format
|
---|
| 672 | * b is diagonal (nb,nb) in a vector storage format
|
---|
| 673 | * a (+)= c * b * transpose(c) (+= if add)
|
---|
| 674 | */
|
---|
| 675 | void
|
---|
| 676 | cmat_transform_diagonal_matrix(double**a,int na, /* a is (na,na) */
|
---|
| 677 | double*b,int nb, /* b is (nb,nb) */
|
---|
| 678 | double**c, /* c is (na,nb) */
|
---|
| 679 | int add)
|
---|
| 680 | {
|
---|
| 681 | int i,j,k;
|
---|
| 682 | double t;
|
---|
| 683 |
|
---|
| 684 | for (i=0; i < na; i++) {
|
---|
| 685 | for (j=0; j <= i; j++) {
|
---|
| 686 | t=0;
|
---|
| 687 | for (k=0; k < nb; k++)
|
---|
| 688 | t += c[i][k] * c[j][k] * b[k];
|
---|
| 689 | if (add)
|
---|
| 690 | a[i][j] += t;
|
---|
| 691 | else
|
---|
| 692 | a[i][j] = t;
|
---|
| 693 | }
|
---|
| 694 | }
|
---|
| 695 | }
|
---|
| 696 |
|
---|
| 697 | /*
|
---|
| 698 | * Argument a contains pointers to the rows of a symmetrix matrix. The
|
---|
| 699 | * in each row is the row number + 1. These rows are stored in
|
---|
| 700 | * contiguous memory starting with 0. Evecs also contains pointers to
|
---|
| 701 | * contiguous memory. N is the dimension.
|
---|
| 702 | */
|
---|
| 703 | void
|
---|
| 704 | cmat_diag(double**a, double*evals, double**evecs, int n,
|
---|
| 705 | int matz, double tol)
|
---|
| 706 | {
|
---|
| 707 | int i,j;
|
---|
| 708 | int diagonal=1;
|
---|
| 709 | double*fv1;
|
---|
| 710 |
|
---|
| 711 | /* I'm having problems with diagonalizing matrices which are already
|
---|
| 712 | * diagonal. So let's first check to see if _a_ is diagonal, and if it
|
---|
| 713 | * is, then just return the diagonal elements in evals and a unit matrix
|
---|
| 714 | * in evecs
|
---|
| 715 | */
|
---|
| 716 |
|
---|
| 717 | for (i=1; i < n; i++) {
|
---|
| 718 | for (j=0; j < i; j++) {
|
---|
| 719 | if (fabs(a[i][j]) > tol) diagonal=0;
|
---|
| 720 | }
|
---|
| 721 | }
|
---|
| 722 |
|
---|
| 723 | if (diagonal) {
|
---|
| 724 | for(i=0; i < n; i++) {
|
---|
| 725 | evals[i] = a[i][i];
|
---|
| 726 | evecs[i][i] = 1.0;
|
---|
| 727 |
|
---|
| 728 | for(j=0; j < i; j++) {
|
---|
| 729 | evecs[i][j] = evecs[j][i] = 0.0;
|
---|
| 730 | }
|
---|
| 731 | }
|
---|
| 732 | eigsort(n,evals,evecs);
|
---|
| 733 | return;
|
---|
| 734 | }
|
---|
| 735 |
|
---|
| 736 | fv1 = (double*) malloc(sizeof(double)*n);
|
---|
| 737 | if (!fv1) {
|
---|
| 738 | fprintf(stderr,"cmat_diag: malloc fv1 failed\n");
|
---|
| 739 | abort();
|
---|
| 740 | }
|
---|
| 741 |
|
---|
| 742 | for(i=0; i < n; i++) {
|
---|
| 743 | for(j=0; j <= i; j++) {
|
---|
| 744 | evecs[i][j] = evecs[j][i] = a[i][j];
|
---|
| 745 | }
|
---|
| 746 | }
|
---|
| 747 |
|
---|
| 748 | tred2(n,evecs,evals,fv1,1);
|
---|
| 749 |
|
---|
| 750 | cmat_transpose_square_matrix(evecs,n);
|
---|
| 751 | tqli(n,evals,evecs,fv1,1,tol);
|
---|
| 752 | cmat_transpose_square_matrix(evecs,n);
|
---|
| 753 |
|
---|
| 754 | eigsort(n,evals,evecs);
|
---|
| 755 |
|
---|
| 756 | free(fv1);
|
---|
| 757 | }
|
---|
| 758 |
|
---|
| 759 | #define dsign(a,b) (((b) >= 0.0) ? fabs(a) : -fabs(a))
|
---|
| 760 |
|
---|
| 761 | static void
|
---|
| 762 | tred2(int n,double** a,double* d,double* e,int matz)
|
---|
| 763 | {
|
---|
| 764 | int i,j,k,l;
|
---|
| 765 | double f,g,h,hh,scale,scale_inv,h_inv;
|
---|
| 766 | if (n == 1) return;
|
---|
| 767 |
|
---|
| 768 | for(i=n-1; i > 0; i--) {
|
---|
| 769 | l = i-1;
|
---|
| 770 | h = 0.0;
|
---|
| 771 | scale = 0.0;
|
---|
| 772 | if(l) {
|
---|
| 773 | for(k=0; k <= l; k++) scale += fabs(a[i][k]);
|
---|
| 774 | if (scale == 0.0) e[i] = a[i][l];
|
---|
| 775 | else {
|
---|
| 776 | scale_inv=1.0/scale;
|
---|
| 777 | for (k=0; k <= l; k++) {
|
---|
| 778 | a[i][k] *= scale_inv;
|
---|
| 779 | h += a[i][k]*a[i][k];
|
---|
| 780 | }
|
---|
| 781 | f=a[i][l];
|
---|
| 782 | g= -(dsign(sqrt(h),f));
|
---|
| 783 | e[i] = scale*g;
|
---|
| 784 | h -= f*g;
|
---|
| 785 | a[i][l] = f-g;
|
---|
| 786 | f = 0.0;
|
---|
| 787 | h_inv=1.0/h;
|
---|
| 788 | for (j=0; j <= l; j++) {
|
---|
| 789 | if (matz) a[j][i] = a[i][j]*h_inv;
|
---|
| 790 | g = 0.0;
|
---|
| 791 | for (k=0; k <= j; k++) g += a[j][k]*a[i][k];
|
---|
| 792 | if (l > j) for (k=j+1; k <= l; k++) g += a[k][j]*a[i][k];
|
---|
| 793 | e[j] = g*h_inv;
|
---|
| 794 | f += e[j]*a[i][j];
|
---|
| 795 | }
|
---|
| 796 | hh = f/(h+h);
|
---|
| 797 | for (j=0; j <= l; j++) {
|
---|
| 798 | f = a[i][j];
|
---|
| 799 | g = e[j] - hh*f;
|
---|
| 800 | e[j] = g;
|
---|
| 801 | for (k=0; k <= j; k++) a[j][k] -= (f*e[k] + g*a[i][k]);
|
---|
| 802 | }
|
---|
| 803 | }
|
---|
| 804 | }
|
---|
| 805 | else {
|
---|
| 806 | e[i] = a[i][l];
|
---|
| 807 | }
|
---|
| 808 | d[i] = h;
|
---|
| 809 | }
|
---|
| 810 | if(matz) d[0] = 0.0;
|
---|
| 811 | e[0] = 0.0;
|
---|
| 812 |
|
---|
| 813 | for(i=0; i < n; i++) {
|
---|
| 814 | l = i-1;
|
---|
| 815 | if (matz) {
|
---|
| 816 | if(d[i]) {
|
---|
| 817 | for(j=0; j <= l; j++) {
|
---|
| 818 | g = 0.0;
|
---|
| 819 | for(k=0; k <= l; k++) g += a[i][k]*a[k][j];
|
---|
| 820 | for(k=0; k <= l; k++) a[k][j] -= g*a[k][i];
|
---|
| 821 | }
|
---|
| 822 | }
|
---|
| 823 | }
|
---|
| 824 | d[i] = a[i][i];
|
---|
| 825 | if(matz) {
|
---|
| 826 | a[i][i] = 1.0;
|
---|
| 827 | if(l >= 0) for (j=0; j<= l; j++) a[i][j] = a[j][i] = 0.0;
|
---|
| 828 | }
|
---|
| 829 | }
|
---|
| 830 | }
|
---|
| 831 |
|
---|
| 832 | static void
|
---|
| 833 | tqli(int n, double* d, double** z, double* e, int matz, double toler)
|
---|
| 834 | {
|
---|
| 835 | register int k;
|
---|
| 836 | int i,l,m,iter;
|
---|
| 837 | double g,r,s,c,p,f,b;
|
---|
| 838 | double azi;
|
---|
| 839 |
|
---|
| 840 | f=0.0;
|
---|
| 841 | if (n == 1) {
|
---|
| 842 | d[0]=z[0][0];
|
---|
| 843 | z[0][0] = 1.0;
|
---|
| 844 | return;
|
---|
| 845 | }
|
---|
| 846 |
|
---|
| 847 | for (i=1; i < n ; i++) e[i-1] = e[i];
|
---|
| 848 | e[n-1] = 0.0;
|
---|
| 849 | for (l=0; l < n; l++) {
|
---|
| 850 | iter = 0;
|
---|
| 851 | L1:
|
---|
| 852 | for (m=l; m < n-1;m++) if (fabs(e[m]) < toler) goto L2;
|
---|
| 853 | m=n-1;
|
---|
| 854 | L2:
|
---|
| 855 | if (m != l) {
|
---|
| 856 | if (iter++ == 30) {
|
---|
| 857 | fprintf (stderr,"tqli not converging %d %g\n",l,e[l]);
|
---|
| 858 | continue;
|
---|
| 859 | }
|
---|
| 860 |
|
---|
| 861 | g = (d[l+1]-d[l])/(2.0*e[l]);
|
---|
| 862 | r = sqrt(g*g + 1.0);
|
---|
| 863 | g = d[m] - d[l] + e[l]/((g + dsign(r,g)));
|
---|
| 864 | s=1.0;
|
---|
| 865 | c=1.0;
|
---|
| 866 | p=0.0;
|
---|
| 867 | for (i=m-1; i >= l; i--) {
|
---|
| 868 | f = s*e[i];
|
---|
| 869 | b = c*e[i];
|
---|
| 870 | if (fabs(f) >= fabs(g)) {
|
---|
| 871 | c = g/f;
|
---|
| 872 | r = sqrt(c*c + 1.0);
|
---|
| 873 | e[i+1] = f*r;
|
---|
| 874 | s=1.0/r;
|
---|
| 875 | c *= s;
|
---|
| 876 | }
|
---|
| 877 | else {
|
---|
| 878 | s = f/g;
|
---|
| 879 | r = sqrt(s*s + 1.0);
|
---|
| 880 | e[i+1] = g*r;
|
---|
| 881 | c = 1.0/r;
|
---|
| 882 | s *= c;
|
---|
| 883 | }
|
---|
| 884 | g = d[i+1] - p;
|
---|
| 885 | r = (d[i]-g)*s + 2.0*c*b;
|
---|
| 886 | p = s*r;
|
---|
| 887 | d[i+1] = g+p;
|
---|
| 888 | g = c*r-b;
|
---|
| 889 |
|
---|
| 890 | if (matz) {
|
---|
| 891 | double *zi = z[i];
|
---|
| 892 | double *zi1 = z[i+1];
|
---|
| 893 | for (k=n; k ; k--,zi++,zi1++) {
|
---|
| 894 | azi = *zi;
|
---|
| 895 | f = *zi1;
|
---|
| 896 | *zi1 = azi*s + c*f;
|
---|
| 897 | *zi = azi*c - s*f;
|
---|
| 898 | }
|
---|
| 899 | }
|
---|
| 900 | }
|
---|
| 901 |
|
---|
| 902 | d[l] -= p;
|
---|
| 903 | e[l] = g;
|
---|
| 904 | e[m] = 0.0;
|
---|
| 905 | goto L1;
|
---|
| 906 | }
|
---|
| 907 | }
|
---|
| 908 | }
|
---|
| 909 |
|
---|
| 910 | static void
|
---|
| 911 | eigsort(int n, double* d, double** v)
|
---|
| 912 | {
|
---|
| 913 | int i,j,k;
|
---|
| 914 | double p;
|
---|
| 915 |
|
---|
| 916 | for(i=0; i < n-1 ; i++) {
|
---|
| 917 | k=i;
|
---|
| 918 | p=d[i];
|
---|
| 919 | for(j=i+1; j < n; j++) {
|
---|
| 920 | if(d[j] < p) {
|
---|
| 921 | k=j;
|
---|
| 922 | p=d[j];
|
---|
| 923 | }
|
---|
| 924 | }
|
---|
| 925 | if(k != i) {
|
---|
| 926 | d[k]=d[i];
|
---|
| 927 | d[i]=p;
|
---|
| 928 | for(j=0; j < n; j++) {
|
---|
| 929 | p=v[j][i];
|
---|
| 930 | v[j][i]=v[j][k];
|
---|
| 931 | v[j][k]=p;
|
---|
| 932 | }
|
---|
| 933 | }
|
---|
| 934 | }
|
---|
| 935 | }
|
---|
| 936 |
|
---|
| 937 | void
|
---|
| 938 | cmat_schmidt(double **C, double *S, int nrow, int nc)
|
---|
| 939 | {
|
---|
| 940 | int i,j,ij;
|
---|
| 941 | int m;
|
---|
| 942 | double vtmp;
|
---|
| 943 | double *v = (double*) malloc(sizeof(double)*nrow);
|
---|
| 944 |
|
---|
| 945 | if (!v) {
|
---|
| 946 | fprintf(stderr,"cmat_schmidt: could not malloc v(%d)\n",nrow);
|
---|
| 947 | abort();
|
---|
| 948 | }
|
---|
| 949 |
|
---|
| 950 | for (m=0; m < nc; m++) {
|
---|
| 951 | v[0] = C[0][m] * S[0];
|
---|
| 952 |
|
---|
| 953 | for (i=ij=1; i < nrow; i++) {
|
---|
| 954 | for (j=0,vtmp=0.0; j < i; j++,ij++) {
|
---|
| 955 | vtmp += C[j][m]*S[ij];
|
---|
| 956 | v[j] += C[i][m]*S[ij];
|
---|
| 957 | }
|
---|
| 958 | v[i] = vtmp + C[i][m]*S[ij];
|
---|
| 959 | ij++;
|
---|
| 960 | }
|
---|
| 961 |
|
---|
| 962 | for (i=0,vtmp=0.0; i < nrow; i++)
|
---|
| 963 | vtmp += v[i]*C[i][m];
|
---|
| 964 |
|
---|
| 965 | if (!vtmp) {
|
---|
| 966 | fprintf(stderr,"cmat_schmidt: bogus\n");
|
---|
| 967 | abort();
|
---|
| 968 | }
|
---|
| 969 |
|
---|
| 970 | if (vtmp < 1.0e-15)
|
---|
| 971 | vtmp = 1.0e-15;
|
---|
| 972 |
|
---|
| 973 | vtmp = 1.0/sqrt(vtmp);
|
---|
| 974 |
|
---|
| 975 | for (i=0; i < nrow; i++) {
|
---|
| 976 | v[i] *= vtmp;
|
---|
| 977 | C[i][m] *= vtmp;
|
---|
| 978 | }
|
---|
| 979 |
|
---|
| 980 | if (m < nc-1) {
|
---|
| 981 | for (i=m+1,vtmp=0.0; i < nc; i++) {
|
---|
| 982 | for (j=0,vtmp=0.0; j < nrow; j++)
|
---|
| 983 | vtmp += v[j] * C[j][i];
|
---|
| 984 | for (j=0; j < nrow; j++)
|
---|
| 985 | C[j][i] -= vtmp * C[j][m];
|
---|
| 986 | }
|
---|
| 987 | }
|
---|
| 988 | }
|
---|
| 989 | }
|
---|
| 990 |
|
---|
| 991 | /* Returns the number of linearly independent vectors
|
---|
| 992 | orthogonal wrt S. */
|
---|
| 993 | int
|
---|
| 994 | cmat_schmidt_tol(double **C, double *S, int nrow, int ncol,
|
---|
| 995 | double tolerance, double *res)
|
---|
| 996 | {
|
---|
| 997 | int i,j,ij;
|
---|
| 998 | int m;
|
---|
| 999 | double vtmp;
|
---|
| 1000 | int northog = 0;
|
---|
| 1001 | double *v = (double*) malloc(sizeof(double)*nrow);
|
---|
| 1002 |
|
---|
| 1003 | if (res) *res = 1.0;
|
---|
| 1004 |
|
---|
| 1005 | if (!v) {
|
---|
| 1006 | fprintf(stderr,"cmat_schmidt_tol: could not malloc v(%d)\n",nrow);
|
---|
| 1007 | abort();
|
---|
| 1008 | }
|
---|
| 1009 |
|
---|
| 1010 | /* Orthonormalize the columns of C wrt S. */
|
---|
| 1011 | for (m=0; m < ncol; m++) {
|
---|
| 1012 | v[0] = C[0][m] * S[0];
|
---|
| 1013 |
|
---|
| 1014 | for (i=ij=1; i < nrow; i++) {
|
---|
| 1015 | for (j=0,vtmp=0.0; j < i; j++,ij++) {
|
---|
| 1016 | vtmp += C[j][m]*S[ij];
|
---|
| 1017 | v[j] += C[i][m]*S[ij];
|
---|
| 1018 | }
|
---|
| 1019 | v[i] = vtmp + C[i][m]*S[ij];
|
---|
| 1020 | ij++;
|
---|
| 1021 | }
|
---|
| 1022 |
|
---|
| 1023 | for (i=0,vtmp=0.0; i < nrow; i++)
|
---|
| 1024 | vtmp += v[i]*C[i][m];
|
---|
| 1025 |
|
---|
| 1026 | if (vtmp < tolerance) continue;
|
---|
| 1027 |
|
---|
| 1028 | if (res && (m == 0 || vtmp < *res)) *res = vtmp;
|
---|
| 1029 |
|
---|
| 1030 | vtmp = 1.0/sqrt(vtmp);
|
---|
| 1031 |
|
---|
| 1032 | for (i=0; i < nrow; i++) {
|
---|
| 1033 | v[i] *= vtmp;
|
---|
| 1034 | C[i][northog] = C[i][m] * vtmp;
|
---|
| 1035 | }
|
---|
| 1036 |
|
---|
| 1037 | for (i=m+1,vtmp=0.0; i < ncol; i++) {
|
---|
| 1038 | for (j=0,vtmp=0.0; j < nrow; j++)
|
---|
| 1039 | vtmp += v[j] * C[j][i];
|
---|
| 1040 | for (j=0; j < nrow; j++)
|
---|
| 1041 | C[j][i] -= vtmp * C[j][northog];
|
---|
| 1042 | }
|
---|
| 1043 | northog++;
|
---|
| 1044 | }
|
---|
| 1045 | return northog;
|
---|
| 1046 | }
|
---|