1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.48 1995/08/17 13:36:24 curfman Exp bsmith $"; 3 #endif 4 5 /* 6 Standard Fortran style matrices 7 */ 8 #include "petsc.h" 9 #include "plapack.h" 10 #include "matimpl.h" 11 #include "math.h" 12 #include "vec/vecimpl.h" 13 #include "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; 367 PetscObject vobj = (PetscObject) ptr; 368 369 if (ptr == 0) { 370 ptr = STDOUT_VIEWER; 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 = ViewerFileGetPointer_Private(ptr); 377 for ( i=0; i<mat->m; i++ ) { 378 v = mat->v + i; 379 for ( j=0; j<mat->n; j++ ) { 380 #if defined(PETSC_COMPLEX) 381 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 382 #else 383 fprintf(fd,"%6.4e ",*v); v += mat->m; 384 #endif 385 } 386 fprintf(fd,"\n"); 387 } 388 } 389 return 0; 390 } 391 392 393 static int MatDestroy_Dense(PetscObject obj) 394 { 395 Mat mat = (Mat) obj; 396 Mat_Dense *l = (Mat_Dense *) mat->data; 397 #if defined(PETSC_LOG) 398 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 399 #endif 400 if (l->pivots) PETSCFREE(l->pivots); 401 PETSCFREE(l); 402 PLogObjectDestroy(mat); 403 PETSCHEADERDESTROY(mat); 404 return 0; 405 } 406 407 static int MatTranspose_Dense(Mat matin,Mat *matout) 408 { 409 Mat_Dense *mat = (Mat_Dense *) matin->data; 410 int k,j; 411 Scalar *v = mat->v, tmp; 412 413 if (!matout) { /* in place transpose */ 414 if (mat->m != mat->n) { 415 SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix"); 416 } 417 for ( j=0; j<mat->m; j++ ) { 418 for ( k=0; k<j; k++ ) { 419 tmp = v[j + k*mat->n]; 420 v[j + k*mat->n] = v[k + j*mat->n]; 421 v[k + j*mat->n] = tmp; 422 } 423 } 424 } 425 else { 426 SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet"); 427 } 428 return 0; 429 } 430 431 static int MatEqual_Dense(Mat matin1,Mat matin2) 432 { 433 Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 434 Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 435 int i; 436 Scalar *v1 = mat1->v, *v2 = mat2->v; 437 if (mat1->m != mat2->m) return 0; 438 if (mat1->n != mat2->n) return 0; 439 for ( i=0; i<mat1->m*mat1->n; i++ ) { 440 if (*v1 != *v2) return 0; 441 v1++; v2++; 442 } 443 return 1; 444 } 445 446 static int MatGetDiagonal_Dense(Mat matin,Vec v) 447 { 448 Mat_Dense *mat = (Mat_Dense *) matin->data; 449 int i, n; 450 Scalar *x; 451 VecGetArray(v,&x); VecGetSize(v,&n); 452 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec"); 453 for ( i=0; i<mat->m; i++ ) { 454 x[i] = mat->v[i*mat->m + i]; 455 } 456 return 0; 457 } 458 459 static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 460 { 461 Mat_Dense *mat = (Mat_Dense *) matin->data; 462 Scalar *l,*r,x,*v; 463 int i,j,m = mat->m, n = mat->n; 464 if (ll) { 465 VecGetArray(ll,&l); VecGetSize(ll,&m); 466 if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size"); 467 for ( i=0; i<m; i++ ) { 468 x = l[i]; 469 v = mat->v + i; 470 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 471 } 472 } 473 if (rr) { 474 VecGetArray(rr,&r); VecGetSize(rr,&n); 475 if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size"); 476 for ( i=0; i<n; i++ ) { 477 x = r[i]; 478 v = mat->v + i*m; 479 for ( j=0; j<m; j++ ) { (*v++) *= x;} 480 } 481 } 482 return 0; 483 } 484 485 static int MatNorm_Dense(Mat matin,MatNormType type,double *norm) 486 { 487 Mat_Dense *mat = (Mat_Dense *) matin->data; 488 Scalar *v = mat->v; 489 double sum = 0.0; 490 int i, j; 491 if (type == NORM_FROBENIUS) { 492 for (i=0; i<mat->n*mat->m; i++ ) { 493 #if defined(PETSC_COMPLEX) 494 sum += real(conj(*v)*(*v)); v++; 495 #else 496 sum += (*v)*(*v); v++; 497 #endif 498 } 499 *norm = sqrt(sum); 500 } 501 else if (type == NORM_1) { 502 *norm = 0.0; 503 for ( j=0; j<mat->n; j++ ) { 504 sum = 0.0; 505 for ( i=0; i<mat->m; i++ ) { 506 #if defined(PETSC_COMPLEX) 507 sum += abs(*v++); 508 #else 509 sum += fabs(*v++); 510 #endif 511 } 512 if (sum > *norm) *norm = sum; 513 } 514 } 515 else if (type == NORM_INFINITY) { 516 *norm = 0.0; 517 for ( j=0; j<mat->m; j++ ) { 518 v = mat->v + j; 519 sum = 0.0; 520 for ( i=0; i<mat->n; i++ ) { 521 #if defined(PETSC_COMPLEX) 522 sum += abs(*v); v += mat->m; 523 #else 524 sum += fabs(*v); v += mat->m; 525 #endif 526 } 527 if (sum > *norm) *norm = sum; 528 } 529 } 530 else { 531 SETERRQ(1,"MatNorm_Dense:No two norm yet"); 532 } 533 return 0; 534 } 535 536 static int MatSetOption_Dense(Mat aijin,MatOption op) 537 { 538 Mat_Dense *aij = (Mat_Dense *) aijin->data; 539 if (op == ROW_ORIENTED) aij->roworiented = 1; 540 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 541 /* doesn't care about sorted rows or columns */ 542 return 0; 543 } 544 545 static int MatZero_Dense(Mat A) 546 { 547 Mat_Dense *l = (Mat_Dense *) A->data; 548 PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 549 return 0; 550 } 551 552 static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 553 { 554 Mat_Dense *l = (Mat_Dense *) A->data; 555 int n = l->n, i, j,ierr,N, *rows; 556 Scalar *slot; 557 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 558 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 559 for ( i=0; i<N; i++ ) { 560 slot = l->v + rows[i]; 561 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 562 } 563 if (diag) { 564 for ( i=0; i<N; i++ ) { 565 slot = l->v + (n+1)*rows[i]; 566 *slot = *diag; 567 } 568 } 569 ISRestoreIndices(is,&rows); 570 return 0; 571 } 572 573 static int MatGetSize_Dense(Mat matin,int *m,int *n) 574 { 575 Mat_Dense *mat = (Mat_Dense *) matin->data; 576 *m = mat->m; *n = mat->n; 577 return 0; 578 } 579 580 static int MatGetArray_Dense(Mat matin,Scalar **array) 581 { 582 Mat_Dense *mat = (Mat_Dense *) matin->data; 583 *array = mat->v; 584 return 0; 585 } 586 587 588 static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 589 { 590 SETERRQ(1,"MatGetSubMatrixInPlace_Dense: not done"); 591 } 592 593 static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 594 { 595 Mat_Dense *mat = (Mat_Dense *) matin->data; 596 int nznew, *smap, i, j, ierr, oldcols = mat->n; 597 int *irow, *icol, nrows, ncols, *cwork; 598 Scalar *vwork, *val; 599 Mat newmat; 600 601 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 602 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 603 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 604 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 605 606 smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 607 cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 608 vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 609 PETSCMEMSET((char*)smap,0,oldcols*sizeof(int)); 610 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 611 612 /* Create and fill new matrix */ 613 ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 614 CHKERRQ(ierr); 615 for (i=0; i<nrows; i++) { 616 nznew = 0; 617 val = mat->v + irow[i]; 618 for (j=0; j<oldcols; j++) { 619 if (smap[j]) { 620 cwork[nznew] = smap[j] - 1; 621 vwork[nznew++] = val[j * mat->m]; 622 } 623 } 624 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 625 CHKERRQ(ierr); 626 } 627 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 628 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 629 630 /* Free work space */ 631 PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 632 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 633 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 634 *submat = newmat; 635 return 0; 636 } 637 638 /* -------------------------------------------------------------------*/ 639 static struct _MatOps MatOps = {MatInsert_Dense, 640 MatGetRow_Dense, MatRestoreRow_Dense, 641 MatMult_Dense, MatMultAdd_Dense, 642 MatMultTrans_Dense, MatMultTransAdd_Dense, 643 MatSolve_Dense,MatSolveAdd_Dense, 644 MatSolveTrans_Dense,MatSolveTransAdd_Dense, 645 MatLUFactor_Dense,MatCholeskyFactor_Dense, 646 MatRelax_Dense, 647 MatTranspose_Dense, 648 MatGetInfo_Dense,MatEqual_Dense, 649 MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 650 0,0, 651 0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0, 652 MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 653 MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 654 MatGetSize_Dense,MatGetSize_Dense,0, 655 0,0,MatGetArray_Dense,0,0, 656 MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense, 657 MatCopyPrivate_Dense}; 658 659 /*@ 660 MatCreateSequentialDense - Creates a sequential dense matrix that 661 is stored in column major order (the usual Fortran 77 manner). Many 662 of the matrix operations use the BLAS and LAPACK routines. 663 664 Input Parameters: 665 . comm - MPI communicator, set to MPI_COMM_SELF 666 . m - number of rows 667 . n - number of column 668 669 Output Parameter: 670 . newmat - the matrix 671 672 .keywords: dense, matrix, LAPACK, BLAS 673 674 .seealso: MatCreate(), MatSetValues() 675 @*/ 676 int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 677 { 678 int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 679 Mat mat; 680 Mat_Dense *l; 681 *newmat = 0; 682 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 683 PLogObjectCreate(mat); 684 l = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l); 685 mat->ops = &MatOps; 686 mat->destroy = MatDestroy_Dense; 687 mat->view = MatView_Dense; 688 mat->data = (void *) l; 689 mat->factor = 0; 690 PLogObjectMemory(mat,sizeof(struct _Mat) + size); 691 692 l->m = m; 693 l->n = n; 694 l->v = (Scalar *) (l + 1); 695 l->pivots = 0; 696 l->roworiented = 1; 697 698 PETSCMEMSET(l->v,0,m*n*sizeof(Scalar)); 699 *newmat = mat; 700 return 0; 701 } 702 703 int MatCreate_Dense(Mat matin,Mat *newmat) 704 { 705 Mat_Dense *m = (Mat_Dense *) matin->data; 706 return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 707 } 708