1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.49 1995/08/17 13:29:40 curfman Exp $"; 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 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 if (flag & SOR_BACKWARD_SWEEP) { 200 for ( i=m-1; i>=0; i-- ) { 201 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 202 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 203 } 204 } 205 } 206 return 0; 207 } 208 209 /* -----------------------------------------------------------------*/ 210 static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 211 { 212 Mat_Dense *mat = (Mat_Dense *) matin->data; 213 Scalar *v = mat->v, *x, *y; 214 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 215 VecGetArray(xx,&x), VecGetArray(yy,&y); 216 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 217 x, &_One, &_DZero, y, &_One ); 218 return 0; 219 } 220 static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 221 { 222 Mat_Dense *mat = (Mat_Dense *) matin->data; 223 Scalar *v = mat->v, *x, *y; 224 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 225 VecGetArray(xx,&x); VecGetArray(yy,&y); 226 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 227 x, &_One, &_DZero, y, &_One ); 228 return 0; 229 } 230 static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 231 { 232 Mat_Dense *mat = (Mat_Dense *) matin->data; 233 Scalar *v = mat->v, *x, *y, *z; 234 int _One=1; Scalar _DOne=1.0; 235 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 236 if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 237 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 238 x, &_One, &_DOne, y, &_One ); 239 return 0; 240 } 241 static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 242 { 243 Mat_Dense *mat = (Mat_Dense *) matin->data; 244 Scalar *v = mat->v, *x, *y, *z; 245 int _One=1; 246 Scalar _DOne=1.0; 247 VecGetArray(xx,&x); VecGetArray(yy,&y); 248 VecGetArray(zz,&z); 249 if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 250 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 251 x, &_One, &_DOne, y, &_One ); 252 return 0; 253 } 254 255 /* -----------------------------------------------------------------*/ 256 static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 257 Scalar **vals) 258 { 259 Mat_Dense *mat = (Mat_Dense *) matin->data; 260 Scalar *v; 261 int i; 262 *ncols = mat->n; 263 if (cols) { 264 *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 265 for ( i=0; i<mat->n; i++ ) *cols[i] = i; 266 } 267 if (vals) { 268 *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 269 v = mat->v + row; 270 for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 271 } 272 return 0; 273 } 274 static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 275 Scalar **vals) 276 { 277 if (cols) { PETSCFREE(*cols); } 278 if (vals) { PETSCFREE(*vals); } 279 return 0; 280 } 281 /* ----------------------------------------------------------------*/ 282 static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 283 int *indexn,Scalar *v,InsertMode addv) 284 { 285 Mat_Dense *mat = (Mat_Dense *) matin->data; 286 int i,j; 287 288 if (!mat->roworiented) { 289 if (addv == INSERTVALUES) { 290 for ( j=0; j<n; j++ ) { 291 for ( i=0; i<m; i++ ) { 292 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 293 } 294 } 295 } 296 else { 297 for ( j=0; j<n; j++ ) { 298 for ( i=0; i<m; i++ ) { 299 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 300 } 301 } 302 } 303 } 304 else { 305 if (addv == INSERTVALUES) { 306 for ( i=0; i<m; i++ ) { 307 for ( j=0; j<n; j++ ) { 308 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 309 } 310 } 311 } 312 else { 313 for ( i=0; i<m; i++ ) { 314 for ( j=0; j<n; j++ ) { 315 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 316 } 317 } 318 } 319 } 320 return 0; 321 } 322 323 /* -----------------------------------------------------------------*/ 324 static int MatCopyPrivate_Dense(Mat matin,Mat *newmat) 325 { 326 Mat_Dense *mat = (Mat_Dense *) matin->data; 327 int ierr; 328 Mat newi; 329 Mat_Dense *l; 330 ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi); 331 CHKERRQ(ierr); 332 l = (Mat_Dense *) newi->data; 333 PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 334 *newmat = newi; 335 return 0; 336 } 337 #include "viewer.h" 338 339 int MatView_Dense(PetscObject obj,Viewer ptr) 340 { 341 Mat matin = (Mat) obj; 342 Mat_Dense *mat = (Mat_Dense *) matin->data; 343 Scalar *v; 344 int i,j; 345 PetscObject vobj = (PetscObject) ptr; 346 347 if (ptr == 0) { 348 ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 349 } 350 if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 351 return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 352 } 353 else { 354 FILE *fd = ViewerFileGetPointer_Private(ptr); 355 for ( i=0; i<mat->m; i++ ) { 356 v = mat->v + i; 357 for ( j=0; j<mat->n; j++ ) { 358 #if defined(PETSC_COMPLEX) 359 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 360 #else 361 fprintf(fd,"%6.4e ",*v); v += mat->m; 362 #endif 363 } 364 fprintf(fd,"\n"); 365 } 366 } 367 return 0; 368 } 369 370 371 static int MatDestroy_Dense(PetscObject obj) 372 { 373 Mat mat = (Mat) obj; 374 Mat_Dense *l = (Mat_Dense *) mat->data; 375 #if defined(PETSC_LOG) 376 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 377 #endif 378 if (l->pivots) PETSCFREE(l->pivots); 379 PETSCFREE(l); 380 PLogObjectDestroy(mat); 381 PETSCHEADERDESTROY(mat); 382 return 0; 383 } 384 385 static int MatTranspose_Dense(Mat matin,Mat *matout) 386 { 387 Mat_Dense *mat = (Mat_Dense *) matin->data; 388 int k,j; 389 Scalar *v = mat->v, tmp; 390 391 if (!matout) { /* in place transpose */ 392 if (mat->m != mat->n) { 393 SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix"); 394 } 395 for ( j=0; j<mat->m; j++ ) { 396 for ( k=0; k<j; k++ ) { 397 tmp = v[j + k*mat->n]; 398 v[j + k*mat->n] = v[k + j*mat->n]; 399 v[k + j*mat->n] = tmp; 400 } 401 } 402 } 403 else { 404 SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet"); 405 } 406 return 0; 407 } 408 409 static int MatEqual_Dense(Mat matin1,Mat matin2) 410 { 411 Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 412 Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 413 int i; 414 Scalar *v1 = mat1->v, *v2 = mat2->v; 415 if (mat1->m != mat2->m) return 0; 416 if (mat1->n != mat2->n) return 0; 417 for ( i=0; i<mat1->m*mat1->n; i++ ) { 418 if (*v1 != *v2) return 0; 419 v1++; v2++; 420 } 421 return 1; 422 } 423 424 static int MatGetDiagonal_Dense(Mat matin,Vec v) 425 { 426 Mat_Dense *mat = (Mat_Dense *) matin->data; 427 int i, n; 428 Scalar *x; 429 VecGetArray(v,&x); VecGetSize(v,&n); 430 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec"); 431 for ( i=0; i<mat->m; i++ ) { 432 x[i] = mat->v[i*mat->m + i]; 433 } 434 return 0; 435 } 436 437 static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 438 { 439 Mat_Dense *mat = (Mat_Dense *) matin->data; 440 Scalar *l,*r,x,*v; 441 int i,j,m = mat->m, n = mat->n; 442 if (ll) { 443 VecGetArray(ll,&l); VecGetSize(ll,&m); 444 if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size"); 445 for ( i=0; i<m; i++ ) { 446 x = l[i]; 447 v = mat->v + i; 448 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 449 } 450 } 451 if (rr) { 452 VecGetArray(rr,&r); VecGetSize(rr,&n); 453 if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size"); 454 for ( i=0; i<n; i++ ) { 455 x = r[i]; 456 v = mat->v + i*m; 457 for ( j=0; j<m; j++ ) { (*v++) *= x;} 458 } 459 } 460 return 0; 461 } 462 463 static int MatNorm_Dense(Mat matin,MatNormType type,double *norm) 464 { 465 Mat_Dense *mat = (Mat_Dense *) matin->data; 466 Scalar *v = mat->v; 467 double sum = 0.0; 468 int i, j; 469 if (type == NORM_FROBENIUS) { 470 for (i=0; i<mat->n*mat->m; i++ ) { 471 #if defined(PETSC_COMPLEX) 472 sum += real(conj(*v)*(*v)); v++; 473 #else 474 sum += (*v)*(*v); v++; 475 #endif 476 } 477 *norm = sqrt(sum); 478 } 479 else if (type == NORM_1) { 480 *norm = 0.0; 481 for ( j=0; j<mat->n; j++ ) { 482 sum = 0.0; 483 for ( i=0; i<mat->m; i++ ) { 484 #if defined(PETSC_COMPLEX) 485 sum += abs(*v++); 486 #else 487 sum += fabs(*v++); 488 #endif 489 } 490 if (sum > *norm) *norm = sum; 491 } 492 } 493 else if (type == NORM_INFINITY) { 494 *norm = 0.0; 495 for ( j=0; j<mat->m; j++ ) { 496 v = mat->v + j; 497 sum = 0.0; 498 for ( i=0; i<mat->n; i++ ) { 499 #if defined(PETSC_COMPLEX) 500 sum += abs(*v); v += mat->m; 501 #else 502 sum += fabs(*v); v += mat->m; 503 #endif 504 } 505 if (sum > *norm) *norm = sum; 506 } 507 } 508 else { 509 SETERRQ(1,"MatNorm_Dense:No two norm yet"); 510 } 511 return 0; 512 } 513 514 static int MatSetOption_Dense(Mat aijin,MatOption op) 515 { 516 Mat_Dense *aij = (Mat_Dense *) aijin->data; 517 if (op == ROW_ORIENTED) aij->roworiented = 1; 518 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 519 /* doesn't care about sorted rows or columns */ 520 return 0; 521 } 522 523 static int MatZero_Dense(Mat A) 524 { 525 Mat_Dense *l = (Mat_Dense *) A->data; 526 PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 527 return 0; 528 } 529 530 static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 531 { 532 Mat_Dense *l = (Mat_Dense *) A->data; 533 int n = l->n, i, j,ierr,N, *rows; 534 Scalar *slot; 535 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 536 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 537 for ( i=0; i<N; i++ ) { 538 slot = l->v + rows[i]; 539 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 540 } 541 if (diag) { 542 for ( i=0; i<N; i++ ) { 543 slot = l->v + (n+1)*rows[i]; 544 *slot = *diag; 545 } 546 } 547 ISRestoreIndices(is,&rows); 548 return 0; 549 } 550 551 static int MatGetSize_Dense(Mat matin,int *m,int *n) 552 { 553 Mat_Dense *mat = (Mat_Dense *) matin->data; 554 *m = mat->m; *n = mat->n; 555 return 0; 556 } 557 558 static int MatGetArray_Dense(Mat matin,Scalar **array) 559 { 560 Mat_Dense *mat = (Mat_Dense *) matin->data; 561 *array = mat->v; 562 return 0; 563 } 564 565 566 static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 567 { 568 SETERRQ(1,"MatGetSubMatrixInPlace_Dense: not done"); 569 } 570 571 static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 572 { 573 Mat_Dense *mat = (Mat_Dense *) matin->data; 574 int nznew, *smap, i, j, ierr, oldcols = mat->n; 575 int *irow, *icol, nrows, ncols, *cwork; 576 Scalar *vwork, *val; 577 Mat newmat; 578 579 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 580 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 581 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 582 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 583 584 smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 585 cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 586 vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 587 PETSCMEMSET((char*)smap,0,oldcols*sizeof(int)); 588 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 589 590 /* Create and fill new matrix */ 591 ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 592 CHKERRQ(ierr); 593 for (i=0; i<nrows; i++) { 594 nznew = 0; 595 val = mat->v + irow[i]; 596 for (j=0; j<oldcols; j++) { 597 if (smap[j]) { 598 cwork[nznew] = smap[j] - 1; 599 vwork[nznew++] = val[j * mat->m]; 600 } 601 } 602 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 603 CHKERRQ(ierr); 604 } 605 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 606 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 607 608 /* Free work space */ 609 PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 610 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 611 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 612 *submat = newmat; 613 return 0; 614 } 615 616 /* -------------------------------------------------------------------*/ 617 static struct _MatOps MatOps = {MatInsert_Dense, 618 MatGetRow_Dense, MatRestoreRow_Dense, 619 MatMult_Dense, MatMultAdd_Dense, 620 MatMultTrans_Dense, MatMultTransAdd_Dense, 621 MatSolve_Dense,MatSolveAdd_Dense, 622 MatSolveTrans_Dense,MatSolveTransAdd_Dense, 623 MatLUFactor_Dense,MatCholeskyFactor_Dense, 624 MatRelax_Dense, 625 MatTranspose_Dense, 626 MatGetInfo_Dense,MatEqual_Dense, 627 MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 628 0,0, 629 0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0, 630 MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 631 MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 632 MatGetSize_Dense,MatGetSize_Dense,0, 633 0,0,MatGetArray_Dense,0,0, 634 MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense, 635 MatCopyPrivate_Dense}; 636 637 /*@ 638 MatCreateSequentialDense - Creates a sequential dense matrix that 639 is stored in column major order (the usual Fortran 77 manner). Many 640 of the matrix operations use the BLAS and LAPACK routines. 641 642 Input Parameters: 643 . comm - MPI communicator, set to MPI_COMM_SELF 644 . m - number of rows 645 . n - number of column 646 647 Output Parameter: 648 . newmat - the matrix 649 650 .keywords: dense, matrix, LAPACK, BLAS 651 652 .seealso: MatCreate(), MatSetValues() 653 @*/ 654 int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 655 { 656 int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 657 Mat mat; 658 Mat_Dense *l; 659 *newmat = 0; 660 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 661 PLogObjectCreate(mat); 662 l = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l); 663 mat->ops = &MatOps; 664 mat->destroy = MatDestroy_Dense; 665 mat->view = MatView_Dense; 666 mat->data = (void *) l; 667 mat->factor = 0; 668 PLogObjectMemory(mat,sizeof(struct _Mat) + size); 669 670 l->m = m; 671 l->n = n; 672 l->v = (Scalar *) (l + 1); 673 l->pivots = 0; 674 l->roworiented = 1; 675 676 PETSCMEMSET(l->v,0,m*n*sizeof(Scalar)); 677 *newmat = mat; 678 return 0; 679 } 680 681 int MatCreate_Dense(Mat matin,Mat *newmat) 682 { 683 Mat_Dense *m = (Mat_Dense *) matin->data; 684 return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 685 } 686