1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.12 1995/03/17 04:56:44 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 defined(PETSC_LOG) 374 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 375 #endif 376 if (l->pivots) FREE(l->pivots); 377 FREE(l); 378 PLogObjectDestroy(mat); 379 PETSCHEADERDESTROY(mat); 380 return 0; 381 } 382 383 static int MatiSDtrans(Mat matin) 384 { 385 MatiSD *mat = (MatiSD *) matin->data; 386 int k,j; 387 Scalar *v = mat->v, tmp; 388 if (mat->m != mat->n) { 389 SETERR(1,"Cannot transpose rectangular dense matrix"); 390 } 391 for ( j=0; j<mat->m; j++ ) { 392 for ( k=0; k<j; k++ ) { 393 tmp = v[j + k*mat->n]; 394 v[j + k*mat->n] = v[k + j*mat->n]; 395 v[k + j*mat->n] = tmp; 396 } 397 } 398 return 0; 399 } 400 401 static int MatiSDequal(Mat matin1,Mat matin2) 402 { 403 MatiSD *mat1 = (MatiSD *) matin1->data; 404 MatiSD *mat2 = (MatiSD *) matin2->data; 405 int i; 406 Scalar *v1 = mat1->v, *v2 = mat2->v; 407 if (mat1->m != mat2->m) return 0; 408 if (mat1->n != mat2->n) return 0; 409 for ( i=0; i<mat1->m*mat1->n; i++ ) { 410 if (*v1 != *v2) return 0; 411 v1++; v2++; 412 } 413 return 1; 414 } 415 416 static int MatiSDgetdiag(Mat matin,Vec v) 417 { 418 MatiSD *mat = (MatiSD *) matin->data; 419 int i, n; 420 Scalar *x; 421 CHKTYPE(v,SEQVECTOR); 422 VecGetArray(v,&x); VecGetSize(v,&n); 423 if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 424 for ( i=0; i<mat->m; i++ ) { 425 x[i] = mat->v[i*mat->m + i]; 426 } 427 return 0; 428 } 429 430 static int MatiSDscale(Mat matin,Vec ll,Vec rr) 431 { 432 MatiSD *mat = (MatiSD *) matin->data; 433 Scalar *l,*r,x,*v; 434 int i,j,m = mat->m, n = mat->n; 435 if (ll) { 436 VecGetArray(ll,&l); VecGetSize(ll,&m); 437 if (m != mat->m) SETERR(1,"Left scaling vector wrong length"); 438 for ( i=0; i<m; i++ ) { 439 x = l[i]; 440 v = mat->v + i; 441 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 442 } 443 } 444 if (rr) { 445 VecGetArray(rr,&r); VecGetSize(rr,&n); 446 if (n != mat->n) SETERR(1,"Right scaling vector wrong length"); 447 for ( i=0; i<n; i++ ) { 448 x = r[i]; 449 v = mat->v + i*m; 450 for ( j=0; j<m; j++ ) { (*v++) *= x;} 451 } 452 } 453 return 0; 454 } 455 456 457 static int MatiSDnorm(Mat matin,int type,double *norm) 458 { 459 MatiSD *mat = (MatiSD *) matin->data; 460 Scalar *v = mat->v; 461 double sum = 0.0; 462 int i, j; 463 if (type == NORM_FROBENIUS) { 464 for (i=0; i<mat->n*mat->m; i++ ) { 465 #if defined(PETSC_COMPLEX) 466 sum += real(conj(*v)*(*v)); v++; 467 #else 468 sum += (*v)*(*v); v++; 469 #endif 470 } 471 *norm = sqrt(sum); 472 } 473 else if (type == NORM_1) { 474 *norm = 0.0; 475 for ( j=0; j<mat->n; j++ ) { 476 sum = 0.0; 477 for ( i=0; i<mat->m; i++ ) { 478 #if defined(PETSC_COMPLEX) 479 sum += abs(*v++); 480 #else 481 sum += fabs(*v++); 482 #endif 483 } 484 if (sum > *norm) *norm = sum; 485 } 486 } 487 else if (type == NORM_INFINITY) { 488 *norm = 0.0; 489 for ( j=0; j<mat->m; j++ ) { 490 v = mat->v + j; 491 sum = 0.0; 492 for ( i=0; i<mat->n; i++ ) { 493 #if defined(PETSC_COMPLEX) 494 sum += abs(*v); v += mat->m; 495 #else 496 sum += fabs(*v); v += mat->m; 497 #endif 498 } 499 if (sum > *norm) *norm = sum; 500 } 501 } 502 else { 503 SETERR(1,"No support for the two norm yet"); 504 } 505 return 0; 506 } 507 508 static int MatiDenseinsopt(Mat aijin,int op) 509 { 510 MatiSD *aij = (MatiSD *) aijin->data; 511 if (op == ROW_ORIENTED) aij->roworiented = 1; 512 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 513 /* doesn't care about sorted rows or columns */ 514 return 0; 515 } 516 517 static int MatiZero(Mat A) 518 { 519 MatiSD *l = (MatiSD *) A->data; 520 MEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 521 return 0; 522 } 523 524 static int MatiZerorows(Mat A,IS is,Scalar *diag) 525 { 526 MatiSD *l = (MatiSD *) A->data; 527 int n = l->n, i, j,ierr,N, *rows; 528 Scalar *slot; 529 ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 530 ierr = ISGetIndices(is,&rows); CHKERR(ierr); 531 for ( i=0; i<N; i++ ) { 532 slot = l->v + rows[i]; 533 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 534 } 535 if (diag) { 536 for ( i=0; i<N; i++ ) { 537 slot = l->v + (n+1)*rows[i]; 538 *slot = *diag; 539 } 540 } 541 ISRestoreIndices(is,&rows); 542 return 0; 543 } 544 /* -------------------------------------------------------------------*/ 545 static struct _MatOps MatOps = {MatiSDinsert, 546 MatiSDgetrow, MatiSDrestorerow, 547 MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd, 548 MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd, 549 MatiSDlufactor,MatiSDchfactor, 550 MatiSDrelax, 551 MatiSDtrans, 552 MatiSDnz,MatiSDmemory,MatiSDequal, 553 MatiSDcopy, 554 MatiSDgetdiag,MatiSDscale,MatiSDnorm, 555 0,0, 556 0, MatiDenseinsopt,MatiZero,MatiZerorows,0, 557 MatiSDlufactorsymbolic,MatiSDlufactornumeric, 558 MatiSDchfactorsymbolic,MatiSDchfactornumeric 559 }; 560 /*@ 561 MatCreateSequentialDense - Creates a sequential dense matrix that 562 is stored in the usual Fortran 77 manner. Many of the matrix 563 operations use the BLAS and LAPACK routines. 564 565 Input Parameters: 566 . m, n - the number of rows and columns in the matrix. 567 568 Output Parameter: 569 . newmat - the matrix created. 570 571 Keywords: dense matrix, lapack, blas 572 @*/ 573 int MatCreateSequentialDense(int m,int n,Mat *newmat) 574 { 575 int size = sizeof(MatiSD) + m*n*sizeof(Scalar); 576 Mat mat; 577 MatiSD *l; 578 *newmat = 0; 579 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSESEQ,MPI_COMM_SELF); 580 PLogObjectCreate(mat); 581 l = (MatiSD *) MALLOC(size); CHKPTR(l); 582 mat->ops = &MatOps; 583 mat->destroy = MatiSDdestroy; 584 mat->view = MatiSDview; 585 mat->data = (void *) l; 586 mat->factor = 0; 587 mat->col = 0; 588 mat->row = 0; 589 590 l->m = m; 591 l->n = n; 592 l->v = (Scalar *) (l + 1); 593 l->pivots = 0; 594 l->roworiented = 1; 595 596 MEMSET(l->v,0,m*n*sizeof(Scalar)); 597 *newmat = mat; 598 return 0; 599 } 600 601 int MatiSDCreate(Mat matin,Mat *newmat) 602 { 603 MatiSD *m = (MatiSD *) matin->data; 604 return MatCreateSequentialDense(m->m,m->n,newmat); 605 } 606