1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.11 1995/03/10 04:44:45 bsmith Exp bsmith $"; 3 #endif 4 5 /* 6 Standard Fortran style matrices 7 */ 8 #include "ptscimpl.h" 9 #include "plapack.h" 10 #include "matimpl.h" 11 #include "math.h" 12 #include "vec/vecimpl.h" 13 14 typedef struct { 15 Scalar *v; 16 int roworiented; 17 int m,n,pad; 18 int *pivots; /* pivots in LU factorization */ 19 } MatiSD; 20 21 22 static int MatiSDnz(Mat matin,int *nz) 23 { 24 MatiSD *mat = (MatiSD *) matin->data; 25 int i,N = mat->m*mat->n,count = 0; 26 Scalar *v = mat->v; 27 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 28 *nz = count; return 0; 29 } 30 static int MatiSDmemory(Mat matin,int *mem) 31 { 32 MatiSD *mat = (MatiSD *) matin->data; 33 *mem = mat->m*mat->n*sizeof(Scalar); return 0; 34 } 35 36 /* ---------------------------------------------------------------*/ 37 /* COMMENT: I have chosen to hide column permutation in the pivots, 38 rather than put it in the Mat->col slot.*/ 39 static int MatiSDlufactor(Mat matin,IS row,IS col) 40 { 41 MatiSD *mat = (MatiSD *) matin->data; 42 int info; 43 if (!mat->pivots) { 44 mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 45 CHKPTR(mat->pivots); 46 } 47 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 48 if (info) SETERR(1,"Bad LU factorization"); 49 matin->factor = FACTOR_LU; 50 return 0; 51 } 52 static int MatiSDlufactorsymbolic(Mat matin,IS row,IS col,Mat *fact) 53 { 54 int ierr; 55 if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 56 return 0; 57 } 58 static int MatiSDlufactornumeric(Mat matin,Mat *fact) 59 { 60 return MatLUFactor(*fact,0,0); 61 } 62 static int MatiSDchfactorsymbolic(Mat matin,IS row,Mat *fact) 63 { 64 int ierr; 65 if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 66 return 0; 67 } 68 static int MatiSDchfactornumeric(Mat matin,Mat *fact) 69 { 70 return MatCholeskyFactor(*fact,0); 71 } 72 static int MatiSDchfactor(Mat matin,IS perm) 73 { 74 MatiSD *mat = (MatiSD *) matin->data; 75 int info; 76 if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 77 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 78 if (info) SETERR(1,"Bad Cholesky factorization"); 79 matin->factor = FACTOR_CHOLESKY; 80 return 0; 81 } 82 83 static int MatiSDsolve(Mat matin,Vec xx,Vec yy) 84 { 85 MatiSD *mat = (MatiSD *) matin->data; 86 int one = 1, info; 87 Scalar *x, *y; 88 VecGetArray(xx,&x); VecGetArray(yy,&y); 89 MEMCPY(y,x,mat->m*sizeof(Scalar)); 90 if (matin->factor == FACTOR_LU) { 91 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 92 y, &mat->m, &info ); 93 } 94 else if (matin->factor == FACTOR_CHOLESKY){ 95 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 96 y, &mat->m, &info ); 97 } 98 else SETERR(1,"Matrix must be factored to solve"); 99 if (info) SETERR(1,"Bad solve"); 100 return 0; 101 } 102 static int MatiSDsolvetrans(Mat matin,Vec xx,Vec yy) 103 { 104 MatiSD *mat = (MatiSD *) matin->data; 105 int one = 1, info; 106 Scalar *x, *y; 107 VecGetArray(xx,&x); VecGetArray(yy,&y); 108 MEMCPY(y,x,mat->m*sizeof(Scalar)); 109 /* assume if pivots exist then LU else Cholesky */ 110 if (mat->pivots) { 111 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 112 y, &mat->m, &info ); 113 } 114 else { 115 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 116 y, &mat->m, &info ); 117 } 118 if (info) SETERR(1,"Bad solve"); 119 return 0; 120 } 121 static int MatiSDsolveadd(Mat matin,Vec xx,Vec zz,Vec yy) 122 { 123 MatiSD *mat = (MatiSD *) matin->data; 124 int one = 1, info,ierr; 125 Scalar *x, *y, sone = 1.0; 126 Vec tmp = 0; 127 VecGetArray(xx,&x); VecGetArray(yy,&y); 128 if (yy == zz) { 129 ierr = VecCreate(yy,&tmp); CHKERR(ierr); 130 ierr = VecCopy(yy,tmp); CHKERR(ierr); 131 } 132 MEMCPY(y,x,mat->m*sizeof(Scalar)); 133 /* assume if pivots exist then LU else Cholesky */ 134 if (mat->pivots) { 135 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 136 y, &mat->m, &info ); 137 } 138 else { 139 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 140 y, &mat->m, &info ); 141 } 142 if (info) SETERR(1,"Bad solve"); 143 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 144 else VecAXPY(&sone,zz,yy); 145 return 0; 146 } 147 static int MatiSDsolvetransadd(Mat matin,Vec xx,Vec zz, Vec yy) 148 { 149 MatiSD *mat = (MatiSD *) matin->data; 150 int one = 1, info,ierr; 151 Scalar *x, *y, sone = 1.0; 152 Vec tmp; 153 VecGetArray(xx,&x); VecGetArray(yy,&y); 154 if (yy == zz) { 155 ierr = VecCreate(yy,&tmp); CHKERR(ierr); 156 ierr = VecCopy(yy,tmp); CHKERR(ierr); 157 } 158 MEMCPY(y,x,mat->m*sizeof(Scalar)); 159 /* assume if pivots exist then LU else Cholesky */ 160 if (mat->pivots) { 161 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 162 y, &mat->m, &info ); 163 } 164 else { 165 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 166 y, &mat->m, &info ); 167 } 168 if (info) SETERR(1,"Bad solve"); 169 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 170 else VecAXPY(&sone,zz,yy); 171 return 0; 172 } 173 /* ------------------------------------------------------------------*/ 174 static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,double shift, 175 int its,Vec xx) 176 { 177 MatiSD *mat = (MatiSD *) matin->data; 178 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 179 int o = 1,ierr, m = mat->m, i; 180 181 if (flag & SOR_ZERO_INITIAL_GUESS) { 182 /* this is a hack fix, should have another version without 183 the second BLdot */ 184 if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0); 185 } 186 VecGetArray(xx,&x); VecGetArray(bb,&b); 187 while (its--) { 188 if (flag & SOR_FORWARD_SWEEP){ 189 for ( i=0; i<m; i++ ) { 190 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 191 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 192 } 193 } 194 if (flag & SOR_BACKWARD_SWEEP) { 195 for ( i=m-1; i>=0; i-- ) { 196 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 197 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 198 } 199 } 200 } 201 return 0; 202 } 203 204 /* -----------------------------------------------------------------*/ 205 static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy) 206 { 207 MatiSD *mat = (MatiSD *) matin->data; 208 Scalar *v = mat->v, *x, *y; 209 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 210 VecGetArray(xx,&x), VecGetArray(yy,&y); 211 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 212 x, &_One, &_DZero, y, &_One ); 213 return 0; 214 } 215 static int MatiSDmult(Mat matin,Vec xx,Vec yy) 216 { 217 MatiSD *mat = (MatiSD *) matin->data; 218 Scalar *v = mat->v, *x, *y; 219 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 220 VecGetArray(xx,&x); VecGetArray(yy,&y); 221 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 222 x, &_One, &_DZero, y, &_One ); 223 return 0; 224 } 225 static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy) 226 { 227 MatiSD *mat = (MatiSD *) matin->data; 228 Scalar *v = mat->v, *x, *y, *z; 229 int _One=1; Scalar _DOne=1.0; 230 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 231 if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 232 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 233 x, &_One, &_DOne, y, &_One ); 234 return 0; 235 } 236 static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy) 237 { 238 MatiSD *mat = (MatiSD *) matin->data; 239 Scalar *v = mat->v, *x, *y, *z; 240 int _One=1; 241 Scalar _DOne=1.0; 242 VecGetArray(xx,&x); VecGetArray(yy,&y); 243 VecGetArray(zz,&z); 244 if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 245 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 246 x, &_One, &_DOne, y, &_One ); 247 return 0; 248 } 249 250 /* -----------------------------------------------------------------*/ 251 static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols, 252 Scalar **vals) 253 { 254 MatiSD *mat = (MatiSD *) matin->data; 255 Scalar *v; 256 int i; 257 *ncols = mat->n; 258 if (cols) { 259 *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 260 for ( i=0; i<mat->n; i++ ) *cols[i] = i; 261 } 262 if (vals) { 263 *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 264 v = mat->v + row; 265 for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 266 } 267 return 0; 268 } 269 static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols, 270 Scalar **vals) 271 { 272 if (cols) { FREE(*cols); } 273 if (vals) { FREE(*vals); } 274 return 0; 275 } 276 /* ----------------------------------------------------------------*/ 277 static int MatiSDinsert(Mat matin,int m,int *indexm,int n, 278 int *indexn,Scalar *v,InsertMode addv) 279 { 280 MatiSD *mat = (MatiSD *) matin->data; 281 int i,j; 282 283 if (!mat->roworiented) { 284 if (addv == InsertValues) { 285 for ( j=0; j<n; j++ ) { 286 if (indexn[j] < 0) {*v += m; continue;} 287 for ( i=0; i<m; i++ ) { 288 if (indexm[i] < 0) {*v++; continue;} 289 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 290 } 291 } 292 } 293 else { 294 for ( j=0; j<n; j++ ) { 295 if (indexn[j] < 0) {*v += m; continue;} 296 for ( i=0; i<m; i++ ) { 297 if (indexm[i] < 0) {*v++; continue;} 298 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 299 } 300 } 301 } 302 } 303 else { 304 if (addv == InsertValues) { 305 for ( i=0; i<m; i++ ) { 306 if (indexm[i] < 0) {*v += n; continue;} 307 for ( j=0; j<n; j++ ) { 308 if (indexn[j] < 0) {*v++; continue;} 309 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 310 } 311 } 312 } 313 else { 314 for ( i=0; i<m; i++ ) { 315 if (indexm[i] < 0) {*v += n; continue;} 316 for ( j=0; j<n; j++ ) { 317 if (indexn[j] < 0) {*v++; continue;} 318 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 319 } 320 } 321 } 322 } 323 return 0; 324 } 325 326 /* -----------------------------------------------------------------*/ 327 static int MatiSDcopy(Mat matin,Mat *newmat) 328 { 329 MatiSD *mat = (MatiSD *) matin->data; 330 int ierr; 331 Mat newi; 332 MatiSD *l; 333 if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0); 334 l = (MatiSD *) newi->data; 335 MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 336 *newmat = newi; 337 return 0; 338 } 339 #include "viewer.h" 340 341 int MatiSDview(PetscObject obj,Viewer ptr) 342 { 343 Mat matin = (Mat) obj; 344 MatiSD *mat = (MatiSD *) matin->data; 345 Scalar *v; 346 int i,j; 347 PetscObject ojb = (PetscObject) ptr; 348 349 if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) { 350 return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v); 351 } 352 else { 353 for ( i=0; i<mat->m; i++ ) { 354 v = mat->v + i; 355 for ( j=0; j<mat->n; j++ ) { 356 #if defined(PETSC_COMPLEX) 357 printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 358 #else 359 printf("%6.4e ",*v); v += mat->m; 360 #endif 361 } 362 printf("\n"); 363 } 364 } 365 return 0; 366 } 367 368 369 static int MatiSDdestroy(PetscObject obj) 370 { 371 Mat mat = (Mat) obj; 372 MatiSD *l = (MatiSD *) mat->data; 373 if (l->pivots) FREE(l->pivots); 374 FREE(l); 375 PETSCHEADERDESTROY(mat); 376 return 0; 377 } 378 379 static int MatiSDtrans(Mat matin) 380 { 381 MatiSD *mat = (MatiSD *) matin->data; 382 int k,j; 383 Scalar *v = mat->v, tmp; 384 if (mat->m != mat->n) { 385 SETERR(1,"Cannot transpose rectangular dense matrix"); 386 } 387 for ( j=0; j<mat->m; j++ ) { 388 for ( k=0; k<j; k++ ) { 389 tmp = v[j + k*mat->n]; 390 v[j + k*mat->n] = v[k + j*mat->n]; 391 v[k + j*mat->n] = tmp; 392 } 393 } 394 return 0; 395 } 396 397 static int MatiSDequal(Mat matin1,Mat matin2) 398 { 399 MatiSD *mat1 = (MatiSD *) matin1->data; 400 MatiSD *mat2 = (MatiSD *) matin2->data; 401 int i; 402 Scalar *v1 = mat1->v, *v2 = mat2->v; 403 if (mat1->m != mat2->m) return 0; 404 if (mat1->n != mat2->n) return 0; 405 for ( i=0; i<mat1->m*mat1->n; i++ ) { 406 if (*v1 != *v2) return 0; 407 v1++; v2++; 408 } 409 return 1; 410 } 411 412 static int MatiSDgetdiag(Mat matin,Vec v) 413 { 414 MatiSD *mat = (MatiSD *) matin->data; 415 int i, n; 416 Scalar *x; 417 CHKTYPE(v,SEQVECTOR); 418 VecGetArray(v,&x); VecGetSize(v,&n); 419 if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 420 for ( i=0; i<mat->m; i++ ) { 421 x[i] = mat->v[i*mat->m + i]; 422 } 423 return 0; 424 } 425 426 static int MatiSDscale(Mat matin,Vec ll,Vec rr) 427 { 428 MatiSD *mat = (MatiSD *) matin->data; 429 Scalar *l,*r,x,*v; 430 int i,j,m = mat->m, n = mat->n; 431 if (ll) { 432 VecGetArray(ll,&l); VecGetSize(ll,&m); 433 if (m != mat->m) SETERR(1,"Left scaling vector wrong length"); 434 for ( i=0; i<m; i++ ) { 435 x = l[i]; 436 v = mat->v + i; 437 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 438 } 439 } 440 if (rr) { 441 VecGetArray(rr,&r); VecGetSize(rr,&n); 442 if (n != mat->n) SETERR(1,"Right scaling vector wrong length"); 443 for ( i=0; i<n; i++ ) { 444 x = r[i]; 445 v = mat->v + i*m; 446 for ( j=0; j<m; j++ ) { (*v++) *= x;} 447 } 448 } 449 return 0; 450 } 451 452 453 static int MatiSDnorm(Mat matin,int type,double *norm) 454 { 455 MatiSD *mat = (MatiSD *) matin->data; 456 Scalar *v = mat->v; 457 double sum = 0.0; 458 int i, j; 459 if (type == NORM_FROBENIUS) { 460 for (i=0; i<mat->n*mat->m; i++ ) { 461 #if defined(PETSC_COMPLEX) 462 sum += real(conj(*v)*(*v)); v++; 463 #else 464 sum += (*v)*(*v); v++; 465 #endif 466 } 467 *norm = sqrt(sum); 468 } 469 else if (type == NORM_1) { 470 *norm = 0.0; 471 for ( j=0; j<mat->n; j++ ) { 472 sum = 0.0; 473 for ( i=0; i<mat->m; i++ ) { 474 #if defined(PETSC_COMPLEX) 475 sum += abs(*v++); 476 #else 477 sum += fabs(*v++); 478 #endif 479 } 480 if (sum > *norm) *norm = sum; 481 } 482 } 483 else if (type == NORM_INFINITY) { 484 *norm = 0.0; 485 for ( j=0; j<mat->m; j++ ) { 486 v = mat->v + j; 487 sum = 0.0; 488 for ( i=0; i<mat->n; i++ ) { 489 #if defined(PETSC_COMPLEX) 490 sum += abs(*v); v += mat->m; 491 #else 492 sum += fabs(*v); v += mat->m; 493 #endif 494 } 495 if (sum > *norm) *norm = sum; 496 } 497 } 498 else { 499 SETERR(1,"No support for the two norm yet"); 500 } 501 return 0; 502 } 503 504 static int MatiDenseinsopt(Mat aijin,int op) 505 { 506 MatiSD *aij = (MatiSD *) aijin->data; 507 if (op == ROW_ORIENTED) aij->roworiented = 1; 508 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 509 /* doesn't care about sorted rows or columns */ 510 return 0; 511 } 512 513 static int MatiZero(Mat A) 514 { 515 MatiSD *l = (MatiSD *) A->data; 516 MEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 517 return 0; 518 } 519 520 static int MatiZerorows(Mat A,IS is,Scalar *diag) 521 { 522 MatiSD *l = (MatiSD *) A->data; 523 int n = l->n, i, j,ierr,N, *rows; 524 Scalar *slot; 525 ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 526 ierr = ISGetIndices(is,&rows); CHKERR(ierr); 527 for ( i=0; i<N; i++ ) { 528 slot = l->v + rows[i]; 529 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 530 } 531 if (diag) { 532 for ( i=0; i<N; i++ ) { 533 slot = l->v + (n+1)*rows[i]; 534 *slot = *diag; 535 } 536 } 537 ISRestoreIndices(is,&rows); 538 return 0; 539 } 540 /* -------------------------------------------------------------------*/ 541 static struct _MatOps MatOps = {MatiSDinsert, 542 MatiSDgetrow, MatiSDrestorerow, 543 MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd, 544 MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd, 545 MatiSDlufactor,MatiSDchfactor, 546 MatiSDrelax, 547 MatiSDtrans, 548 MatiSDnz,MatiSDmemory,MatiSDequal, 549 MatiSDcopy, 550 MatiSDgetdiag,MatiSDscale,MatiSDnorm, 551 0,0, 552 0, MatiDenseinsopt,MatiZero,MatiZerorows,0, 553 MatiSDlufactorsymbolic,MatiSDlufactornumeric, 554 MatiSDchfactorsymbolic,MatiSDchfactornumeric 555 }; 556 /*@ 557 MatCreateSequentialDense - Creates a sequential dense matrix that 558 is stored in the usual Fortran 77 manner. Many of the matrix 559 operations use the BLAS and LAPACK routines. 560 561 Input Parameters: 562 . m, n - the number of rows and columns in the matrix. 563 564 Output Parameter: 565 . newmat - the matrix created. 566 567 Keywords: dense matrix, lapack, blas 568 @*/ 569 int MatCreateSequentialDense(int m,int n,Mat *newmat) 570 { 571 int size = sizeof(MatiSD) + m*n*sizeof(Scalar); 572 Mat mat; 573 MatiSD *l; 574 *newmat = 0; 575 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSESEQ,MPI_COMM_SELF); 576 l = (MatiSD *) MALLOC(size); CHKPTR(l); 577 mat->ops = &MatOps; 578 mat->destroy = MatiSDdestroy; 579 mat->view = MatiSDview; 580 mat->data = (void *) l; 581 mat->factor = 0; 582 mat->col = 0; 583 mat->row = 0; 584 585 l->m = m; 586 l->n = n; 587 l->v = (Scalar *) (l + 1); 588 l->pivots = 0; 589 l->roworiented = 1; 590 591 MEMSET(l->v,0,m*n*sizeof(Scalar)); 592 *newmat = mat; 593 return 0; 594 } 595 596 int MatiSDCreate(Mat matin,Mat *newmat) 597 { 598 MatiSD *m = (MatiSD *) matin->data; 599 return MatCreateSequentialDense(m->m,m->n,newmat); 600 } 601