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