1 2 3 #ifndef lint 4 static char vcid[] = "$Id: dense.c,v 1.41 1995/06/14 17:24:00 bsmith Exp bsmith $"; 5 #endif 6 7 /* 8 Standard Fortran style matrices 9 */ 10 #include "ptscimpl.h" 11 #include "plapack.h" 12 #include "matimpl.h" 13 #include "math.h" 14 #include "vec/vecimpl.h" 15 #include "pviewer.h" 16 17 typedef struct { 18 Scalar *v; 19 int roworiented; 20 int m,n,pad; 21 int *pivots; /* pivots in LU factorization */ 22 } Mat_Dense; 23 24 static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz, 25 int *nzalloc,int *mem) 26 { 27 Mat_Dense *mat = (Mat_Dense *) matin->data; 28 int i,N = mat->m*mat->n,count = 0; 29 Scalar *v = mat->v; 30 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 31 *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar); 32 return 0; 33 } 34 35 /* ---------------------------------------------------------------*/ 36 /* COMMENT: I have chosen to hide column permutation in the pivots, 37 rather than put it in the Mat->col slot.*/ 38 static int MatLUFactor_Dense(Mat matin,IS row,IS col) 39 { 40 Mat_Dense *mat = (Mat_Dense *) matin->data; 41 int info; 42 if (!mat->pivots) { 43 mat->pivots = (int *) PETSCMALLOC( mat->m*sizeof(int) ); 44 CHKPTRQ(mat->pivots); 45 } 46 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 47 if (info) SETERRQ(1,"Bad LU factorization"); 48 matin->factor = FACTOR_LU; 49 return 0; 50 } 51 static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact) 52 { 53 int ierr; 54 if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0); 55 return 0; 56 } 57 static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 58 { 59 return MatLUFactor(*fact,0,0); 60 } 61 static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,Mat *fact) 62 { 63 int ierr; 64 if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0); 65 return 0; 66 } 67 static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact) 68 { 69 return MatCholeskyFactor(*fact,0); 70 } 71 static int MatCholeskyFactor_Dense(Mat matin,IS perm) 72 { 73 Mat_Dense *mat = (Mat_Dense *) matin->data; 74 int info; 75 if (mat->pivots) {PETSCFREE(mat->pivots); mat->pivots = 0;} 76 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 77 if (info) SETERRQ(1,"Bad Cholesky factorization"); 78 matin->factor = FACTOR_CHOLESKY; 79 return 0; 80 } 81 82 static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 83 { 84 Mat_Dense *mat = (Mat_Dense *) matin->data; 85 int one = 1, info; 86 Scalar *x, *y; 87 VecGetArray(xx,&x); VecGetArray(yy,&y); 88 PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 89 if (matin->factor == FACTOR_LU) { 90 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 91 y, &mat->m, &info ); 92 } 93 else if (matin->factor == FACTOR_CHOLESKY){ 94 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 95 y, &mat->m, &info ); 96 } 97 else SETERRQ(1,"Matrix must be factored to solve"); 98 if (info) SETERRQ(1,"Bad solve"); 99 return 0; 100 } 101 static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 102 { 103 Mat_Dense *mat = (Mat_Dense *) matin->data; 104 int one = 1, info; 105 Scalar *x, *y; 106 VecGetArray(xx,&x); VecGetArray(yy,&y); 107 PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 108 /* assume if pivots exist then LU else Cholesky */ 109 if (mat->pivots) { 110 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 111 y, &mat->m, &info ); 112 } 113 else { 114 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 115 y, &mat->m, &info ); 116 } 117 if (info) SETERRQ(1,"Bad solve"); 118 return 0; 119 } 120 static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 121 { 122 Mat_Dense *mat = (Mat_Dense *) matin->data; 123 int one = 1, info,ierr; 124 Scalar *x, *y, sone = 1.0; 125 Vec tmp = 0; 126 VecGetArray(xx,&x); VecGetArray(yy,&y); 127 if (yy == zz) { 128 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 129 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 130 } 131 PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 132 /* assume if pivots exist then LU else Cholesky */ 133 if (mat->pivots) { 134 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 135 y, &mat->m, &info ); 136 } 137 else { 138 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 139 y, &mat->m, &info ); 140 } 141 if (info) SETERRQ(1,"Bad solve"); 142 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 143 else VecAXPY(&sone,zz,yy); 144 return 0; 145 } 146 static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 147 { 148 Mat_Dense *mat = (Mat_Dense *) matin->data; 149 int one = 1, info,ierr; 150 Scalar *x, *y, sone = 1.0; 151 Vec tmp; 152 VecGetArray(xx,&x); VecGetArray(yy,&y); 153 if (yy == zz) { 154 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 155 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 156 } 157 PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 158 /* assume if pivots exist then LU else Cholesky */ 159 if (mat->pivots) { 160 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 161 y, &mat->m, &info ); 162 } 163 else { 164 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 165 y, &mat->m, &info ); 166 } 167 if (info) SETERRQ(1,"Bad solve"); 168 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 169 else VecAXPY(&sone,zz,yy); 170 return 0; 171 } 172 /* ------------------------------------------------------------------*/ 173 static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag, 174 double shift,int its,Vec xx) 175 { 176 Mat_Dense *mat = (Mat_Dense *) matin->data; 177 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 178 int o = 1,ierr, m = mat->m, i; 179 180 if (flag & SOR_ZERO_INITIAL_GUESS) { 181 /* this is a hack fix, should have another version without 182 the second BLdot */ 183 if ((ierr = VecSet(&zero,xx))) SETERRQ(ierr,0); 184 } 185 VecGetArray(xx,&x); VecGetArray(bb,&b); 186 while (its--) { 187 if (flag & SOR_FORWARD_SWEEP){ 188 for ( i=0; i<m; i++ ) { 189 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 190 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 191 } 192 } 193 if (flag & SOR_BACKWARD_SWEEP) { 194 for ( i=m-1; i>=0; i-- ) { 195 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 196 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 197 } 198 } 199 } 200 return 0; 201 } 202 203 /* -----------------------------------------------------------------*/ 204 static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 205 { 206 Mat_Dense *mat = (Mat_Dense *) matin->data; 207 Scalar *v = mat->v, *x, *y; 208 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 209 VecGetArray(xx,&x), VecGetArray(yy,&y); 210 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 211 x, &_One, &_DZero, y, &_One ); 212 return 0; 213 } 214 static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 215 { 216 Mat_Dense *mat = (Mat_Dense *) matin->data; 217 Scalar *v = mat->v, *x, *y; 218 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 219 VecGetArray(xx,&x); VecGetArray(yy,&y); 220 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 221 x, &_One, &_DZero, y, &_One ); 222 return 0; 223 } 224 static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 225 { 226 Mat_Dense *mat = (Mat_Dense *) matin->data; 227 Scalar *v = mat->v, *x, *y, *z; 228 int _One=1; Scalar _DOne=1.0; 229 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 230 if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 231 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 232 x, &_One, &_DOne, y, &_One ); 233 return 0; 234 } 235 static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 236 { 237 Mat_Dense *mat = (Mat_Dense *) matin->data; 238 Scalar *v = mat->v, *x, *y, *z; 239 int _One=1; 240 Scalar _DOne=1.0; 241 VecGetArray(xx,&x); VecGetArray(yy,&y); 242 VecGetArray(zz,&z); 243 if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 244 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 245 x, &_One, &_DOne, y, &_One ); 246 return 0; 247 } 248 249 /* -----------------------------------------------------------------*/ 250 static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 251 Scalar **vals) 252 { 253 Mat_Dense *mat = (Mat_Dense *) matin->data; 254 Scalar *v; 255 int i; 256 *ncols = mat->n; 257 if (cols) { 258 *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 259 for ( i=0; i<mat->n; i++ ) *cols[i] = i; 260 } 261 if (vals) { 262 *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 263 v = mat->v + row; 264 for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 265 } 266 return 0; 267 } 268 static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 269 Scalar **vals) 270 { 271 if (cols) { PETSCFREE(*cols); } 272 if (vals) { PETSCFREE(*vals); } 273 return 0; 274 } 275 /* ----------------------------------------------------------------*/ 276 static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 277 int *indexn,Scalar *v,InsertMode addv) 278 { 279 Mat_Dense *mat = (Mat_Dense *) matin->data; 280 int i,j; 281 282 if (!mat->roworiented) { 283 if (addv == INSERTVALUES) { 284 for ( j=0; j<n; j++ ) { 285 for ( i=0; i<m; i++ ) { 286 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 287 } 288 } 289 } 290 else { 291 for ( j=0; j<n; j++ ) { 292 for ( i=0; i<m; i++ ) { 293 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 294 } 295 } 296 } 297 } 298 else { 299 if (addv == INSERTVALUES) { 300 for ( i=0; i<m; i++ ) { 301 for ( j=0; j<n; j++ ) { 302 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 303 } 304 } 305 } 306 else { 307 for ( i=0; i<m; i++ ) { 308 for ( j=0; j<n; j++ ) { 309 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 310 } 311 } 312 } 313 } 314 return 0; 315 } 316 317 /* -----------------------------------------------------------------*/ 318 static int MatCopyPrivate_Dense(Mat matin,Mat *newmat) 319 { 320 Mat_Dense *mat = (Mat_Dense *) matin->data; 321 int ierr; 322 Mat newi; 323 Mat_Dense *l; 324 if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi))) 325 SETERRQ(ierr,0); 326 l = (Mat_Dense *) newi->data; 327 PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 328 *newmat = newi; 329 return 0; 330 } 331 #include "viewer.h" 332 333 int MatView_Dense(PetscObject obj,Viewer ptr) 334 { 335 Mat matin = (Mat) obj; 336 Mat_Dense *mat = (Mat_Dense *) matin->data; 337 Scalar *v; 338 int i,j; 339 PetscObject vobj = (PetscObject) ptr; 340 341 if (ptr == 0) { 342 ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 343 } 344 if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 345 return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 346 } 347 else { 348 FILE *fd = ViewerFileGetPointer_Private(ptr); 349 for ( i=0; i<mat->m; i++ ) { 350 v = mat->v + i; 351 for ( j=0; j<mat->n; j++ ) { 352 #if defined(PETSC_COMPLEX) 353 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 354 #else 355 fprintf(fd,"%6.4e ",*v); v += mat->m; 356 #endif 357 } 358 fprintf(fd,"\n"); 359 } 360 } 361 return 0; 362 } 363 364 365 static int MatDestroy_Dense(PetscObject obj) 366 { 367 Mat mat = (Mat) obj; 368 Mat_Dense *l = (Mat_Dense *) mat->data; 369 #if defined(PETSC_LOG) 370 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 371 #endif 372 if (l->pivots) PETSCFREE(l->pivots); 373 PETSCFREE(l); 374 PLogObjectDestroy(mat); 375 PETSCHEADERDESTROY(mat); 376 return 0; 377 } 378 379 static int MatTranspose_Dense(Mat matin,Mat *matout) 380 { 381 Mat_Dense *mat = (Mat_Dense *) matin->data; 382 int k,j; 383 Scalar *v = mat->v, tmp; 384 385 if (!matout) { /* in place transpose */ 386 if (mat->m != mat->n) { 387 SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix"); 388 } 389 for ( j=0; j<mat->m; j++ ) { 390 for ( k=0; k<j; k++ ) { 391 tmp = v[j + k*mat->n]; 392 v[j + k*mat->n] = v[k + j*mat->n]; 393 v[k + j*mat->n] = tmp; 394 } 395 } 396 } 397 else { 398 SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet"); 399 } 400 return 0; 401 } 402 403 static int MatEqual_Dense(Mat matin1,Mat matin2) 404 { 405 Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 406 Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 407 int i; 408 Scalar *v1 = mat1->v, *v2 = mat2->v; 409 if (mat1->m != mat2->m) return 0; 410 if (mat1->n != mat2->n) return 0; 411 for ( i=0; i<mat1->m*mat1->n; i++ ) { 412 if (*v1 != *v2) return 0; 413 v1++; v2++; 414 } 415 return 1; 416 } 417 418 static int MatGetDiagonal_Dense(Mat matin,Vec v) 419 { 420 Mat_Dense *mat = (Mat_Dense *) matin->data; 421 int i, n; 422 Scalar *x; 423 VecGetArray(v,&x); VecGetSize(v,&n); 424 if (n != mat->m) SETERRQ(1,"Nonconforming matrix and vector"); 425 for ( i=0; i<mat->m; i++ ) { 426 x[i] = mat->v[i*mat->m + i]; 427 } 428 return 0; 429 } 430 431 static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 432 { 433 Mat_Dense *mat = (Mat_Dense *) matin->data; 434 Scalar *l,*r,x,*v; 435 int i,j,m = mat->m, n = mat->n; 436 if (ll) { 437 VecGetArray(ll,&l); VecGetSize(ll,&m); 438 if (m != mat->m) SETERRQ(1,"Left scaling vector wrong length"); 439 for ( i=0; i<m; i++ ) { 440 x = l[i]; 441 v = mat->v + i; 442 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 443 } 444 } 445 if (rr) { 446 VecGetArray(rr,&r); VecGetSize(rr,&n); 447 if (n != mat->n) SETERRQ(1,"Right scaling vector wrong length"); 448 for ( i=0; i<n; i++ ) { 449 x = r[i]; 450 v = mat->v + i*m; 451 for ( j=0; j<m; j++ ) { (*v++) *= x;} 452 } 453 } 454 return 0; 455 } 456 457 458 static int MatNorm_Dense(Mat matin,int type,double *norm) 459 { 460 Mat_Dense *mat = (Mat_Dense *) matin->data; 461 Scalar *v = mat->v; 462 double sum = 0.0; 463 int i, j; 464 if (type == NORM_FROBENIUS) { 465 for (i=0; i<mat->n*mat->m; i++ ) { 466 #if defined(PETSC_COMPLEX) 467 sum += real(conj(*v)*(*v)); v++; 468 #else 469 sum += (*v)*(*v); v++; 470 #endif 471 } 472 *norm = sqrt(sum); 473 } 474 else if (type == NORM_1) { 475 *norm = 0.0; 476 for ( j=0; j<mat->n; j++ ) { 477 sum = 0.0; 478 for ( i=0; i<mat->m; i++ ) { 479 #if defined(PETSC_COMPLEX) 480 sum += abs(*v++); 481 #else 482 sum += fabs(*v++); 483 #endif 484 } 485 if (sum > *norm) *norm = sum; 486 } 487 } 488 else if (type == NORM_INFINITY) { 489 *norm = 0.0; 490 for ( j=0; j<mat->m; j++ ) { 491 v = mat->v + j; 492 sum = 0.0; 493 for ( i=0; i<mat->n; i++ ) { 494 #if defined(PETSC_COMPLEX) 495 sum += abs(*v); v += mat->m; 496 #else 497 sum += fabs(*v); v += mat->m; 498 #endif 499 } 500 if (sum > *norm) *norm = sum; 501 } 502 } 503 else { 504 SETERRQ(1,"No support for the two norm yet"); 505 } 506 return 0; 507 } 508 509 static int MatSetOption_Dense(Mat aijin,MatOption op) 510 { 511 Mat_Dense *aij = (Mat_Dense *) aijin->data; 512 if (op == ROW_ORIENTED) aij->roworiented = 1; 513 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 514 /* doesn't care about sorted rows or columns */ 515 return 0; 516 } 517 518 static int MatZero_Dense(Mat A) 519 { 520 Mat_Dense *l = (Mat_Dense *) A->data; 521 PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 522 return 0; 523 } 524 525 static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 526 { 527 Mat_Dense *l = (Mat_Dense *) A->data; 528 int n = l->n, i, j,ierr,N, *rows; 529 Scalar *slot; 530 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 531 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 532 for ( i=0; i<N; i++ ) { 533 slot = l->v + rows[i]; 534 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 535 } 536 if (diag) { 537 for ( i=0; i<N; i++ ) { 538 slot = l->v + (n+1)*rows[i]; 539 *slot = *diag; 540 } 541 } 542 ISRestoreIndices(is,&rows); 543 return 0; 544 } 545 546 static int MatGetSize_Dense(Mat matin,int *m,int *n) 547 { 548 Mat_Dense *mat = (Mat_Dense *) matin->data; 549 *m = mat->m; *n = mat->n; 550 return 0; 551 } 552 553 static int MatGetArray_Dense(Mat matin,Scalar **array) 554 { 555 Mat_Dense *mat = (Mat_Dense *) matin->data; 556 *array = mat->v; 557 return 0; 558 } 559 560 561 static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 562 { 563 SETERRQ(1,"MatGetSubMatrixInPlace_Dense not yet done."); 564 } 565 566 static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 567 { 568 Mat_Dense *mat = (Mat_Dense *) matin->data; 569 int nznew, *smap, i, j, ierr, oldcols = mat->n; 570 int *irow, *icol, nrows, ncols, *cwork; 571 Scalar *vwork, *val; 572 Mat newmat; 573 574 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 575 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 576 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 577 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 578 579 smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 580 cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 581 vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 582 PETSCMEMSET((char*)smap,0,oldcols*sizeof(int)); 583 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 584 585 /* Create and fill new matrix */ 586 ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 587 CHKERRQ(ierr); 588 for (i=0; i<nrows; i++) { 589 nznew = 0; 590 val = mat->v + irow[i]; 591 for (j=0; j<oldcols; j++) { 592 if (smap[j]) { 593 cwork[nznew] = smap[j] - 1; 594 vwork[nznew++] = val[j * mat->m]; 595 } 596 } 597 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 598 CHKERRQ(ierr); 599 } 600 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 601 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 602 603 /* Free work space */ 604 PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 605 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 606 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 607 *submat = newmat; 608 return 0; 609 } 610 611 /* -------------------------------------------------------------------*/ 612 static struct _MatOps MatOps = {MatInsert_Dense, 613 MatGetRow_Dense, MatRestoreRow_Dense, 614 MatMult_Dense, MatMultAdd_Dense, 615 MatMultTrans_Dense, MatMultTransAdd_Dense, 616 MatSolve_Dense,MatSolveAdd_Dense, 617 MatSolveTrans_Dense,MatSolveTransAdd_Dense, 618 MatLUFactor_Dense,MatCholeskyFactor_Dense, 619 MatRelax_Dense, 620 MatTranspose_Dense, 621 MatGetInfo_Dense,MatEqual_Dense, 622 MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 623 0,0, 624 0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0, 625 MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 626 MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 627 MatGetSize_Dense,MatGetSize_Dense,0, 628 0,0,MatGetArray_Dense,0,0, 629 MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense, 630 MatCopyPrivate_Dense}; 631 632 /*@ 633 MatCreateSequentialDense - Creates a sequential dense matrix that 634 is stored in column major order (the usual Fortran 77 manner). Many 635 of the matrix operations use the BLAS and LAPACK routines. 636 637 Input Parameters: 638 . comm - MPI communicator, set to MPI_COMM_SELF 639 . m - number of rows 640 . n - number of column 641 642 Output Parameter: 643 . newmat - the matrix 644 645 .keywords: dense, matrix, LAPACK, BLAS 646 647 .seealso: MatCreate(), MatSetValues() 648 @*/ 649 int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 650 { 651 int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 652 Mat mat; 653 Mat_Dense *l; 654 *newmat = 0; 655 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 656 PLogObjectCreate(mat); 657 l = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l); 658 mat->ops = &MatOps; 659 mat->destroy = MatDestroy_Dense; 660 mat->view = MatView_Dense; 661 mat->data = (void *) l; 662 mat->factor = 0; 663 664 l->m = m; 665 l->n = n; 666 l->v = (Scalar *) (l + 1); 667 l->pivots = 0; 668 l->roworiented = 1; 669 670 PETSCMEMSET(l->v,0,m*n*sizeof(Scalar)); 671 *newmat = mat; 672 return 0; 673 } 674 675 int MatCreate_Dense(Mat matin,Mat *newmat) 676 { 677 Mat_Dense *m = (Mat_Dense *) matin->data; 678 return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 679 } 680