1 2 3 #ifndef lint 4 static char vcid[] = "$Id: dense.c,v 1.44 1995/07/13 14:51:42 curfman 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,double f) 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,"MatLUFactor_Dense:Bad LU factorization"); 48 matin->factor = FACTOR_LU; 49 return 0; 50 } 51 static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,double f, 52 Mat *fact) 53 { 54 int ierr; 55 ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr); 56 return 0; 57 } 58 static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 59 { 60 return MatLUFactor(*fact,0,0,1.0); 61 } 62 static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,double f,Mat *fact) 63 { 64 int ierr; 65 ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr); 66 return 0; 67 } 68 static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact) 69 { 70 return MatCholeskyFactor(*fact,0,1.0); 71 } 72 static int MatCholeskyFactor_Dense(Mat matin,IS perm,double f) 73 { 74 Mat_Dense *mat = (Mat_Dense *) matin->data; 75 int info; 76 if (mat->pivots) {PETSCFREE(mat->pivots); mat->pivots = 0;} 77 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 78 if (info) SETERRQ(1,"MatCholeskyFactor_Dense:Bad factorization"); 79 matin->factor = FACTOR_CHOLESKY; 80 return 0; 81 } 82 83 static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 84 { 85 Mat_Dense *mat = (Mat_Dense *) matin->data; 86 int one = 1, info; 87 Scalar *x, *y; 88 VecGetArray(xx,&x); VecGetArray(yy,&y); 89 PETSCMEMCPY(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 SETERRQ(1,"MatSolve_Dense:Matrix must be factored to solve"); 99 if (info) SETERRQ(1,"MatSolve_Dense:Bad solve"); 100 return 0; 101 } 102 static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 103 { 104 Mat_Dense *mat = (Mat_Dense *) matin->data; 105 int one = 1, info; 106 Scalar *x, *y; 107 VecGetArray(xx,&x); VecGetArray(yy,&y); 108 PETSCMEMCPY(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) SETERRQ(1,"MatSolveTrans_Dense:Bad solve"); 119 return 0; 120 } 121 static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 122 { 123 Mat_Dense *mat = (Mat_Dense *) 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 = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 130 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 131 } 132 PETSCMEMCPY(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) SETERRQ(1,"MatSolveAdd_Dense:Bad solve"); 143 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 144 else VecAXPY(&sone,zz,yy); 145 return 0; 146 } 147 static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 148 { 149 Mat_Dense *mat = (Mat_Dense *) 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 = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 156 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 157 } 158 PETSCMEMCPY(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) SETERRQ(1,"MatSolveTransAdd_Dense: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 MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag, 175 double shift,int its,Vec xx) 176 { 177 Mat_Dense *mat = (Mat_Dense *) 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 ierr = VecSet(&zero,xx); CHKERRQ(ierr); 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 MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 206 { 207 Mat_Dense *mat = (Mat_Dense *) 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 MatMult_Dense(Mat matin,Vec xx,Vec yy) 216 { 217 Mat_Dense *mat = (Mat_Dense *) 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 MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 226 { 227 Mat_Dense *mat = (Mat_Dense *) 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) PETSCMEMCPY(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 MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 237 { 238 Mat_Dense *mat = (Mat_Dense *) 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) PETSCMEMCPY(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 MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 252 Scalar **vals) 253 { 254 Mat_Dense *mat = (Mat_Dense *) matin->data; 255 Scalar *v; 256 int i; 257 *ncols = mat->n; 258 if (cols) { 259 *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 260 for ( i=0; i<mat->n; i++ ) *cols[i] = i; 261 } 262 if (vals) { 263 *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*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 MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 270 Scalar **vals) 271 { 272 if (cols) { PETSCFREE(*cols); } 273 if (vals) { PETSCFREE(*vals); } 274 return 0; 275 } 276 /* ----------------------------------------------------------------*/ 277 static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 278 int *indexn,Scalar *v,InsertMode addv) 279 { 280 Mat_Dense *mat = (Mat_Dense *) matin->data; 281 int i,j; 282 283 if (!mat->roworiented) { 284 if (addv == INSERTVALUES) { 285 for ( j=0; j<n; j++ ) { 286 for ( i=0; i<m; i++ ) { 287 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 288 } 289 } 290 } 291 else { 292 for ( j=0; j<n; j++ ) { 293 for ( i=0; i<m; i++ ) { 294 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 295 } 296 } 297 } 298 } 299 else { 300 if (addv == INSERTVALUES) { 301 for ( i=0; i<m; i++ ) { 302 for ( j=0; j<n; j++ ) { 303 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 304 } 305 } 306 } 307 else { 308 for ( i=0; i<m; i++ ) { 309 for ( j=0; j<n; j++ ) { 310 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 311 } 312 } 313 } 314 } 315 return 0; 316 } 317 318 /* -----------------------------------------------------------------*/ 319 static int MatCopyPrivate_Dense(Mat matin,Mat *newmat) 320 { 321 Mat_Dense *mat = (Mat_Dense *) matin->data; 322 int ierr; 323 Mat newi; 324 Mat_Dense *l; 325 ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi); 326 CHKERRQ(ierr); 327 l = (Mat_Dense *) newi->data; 328 PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 329 *newmat = newi; 330 return 0; 331 } 332 #include "viewer.h" 333 334 int MatView_Dense(PetscObject obj,Viewer ptr) 335 { 336 Mat matin = (Mat) obj; 337 Mat_Dense *mat = (Mat_Dense *) matin->data; 338 Scalar *v; 339 int i,j; 340 PetscObject vobj = (PetscObject) ptr; 341 342 if (ptr == 0) { 343 ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 344 } 345 if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 346 return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 347 } 348 else { 349 FILE *fd = ViewerFileGetPointer_Private(ptr); 350 for ( i=0; i<mat->m; i++ ) { 351 v = mat->v + i; 352 for ( j=0; j<mat->n; j++ ) { 353 #if defined(PETSC_COMPLEX) 354 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 355 #else 356 fprintf(fd,"%6.4e ",*v); v += mat->m; 357 #endif 358 } 359 fprintf(fd,"\n"); 360 } 361 } 362 return 0; 363 } 364 365 366 static int MatDestroy_Dense(PetscObject obj) 367 { 368 Mat mat = (Mat) obj; 369 Mat_Dense *l = (Mat_Dense *) mat->data; 370 #if defined(PETSC_LOG) 371 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 372 #endif 373 if (l->pivots) PETSCFREE(l->pivots); 374 PETSCFREE(l); 375 PLogObjectDestroy(mat); 376 PETSCHEADERDESTROY(mat); 377 return 0; 378 } 379 380 static int MatTranspose_Dense(Mat matin,Mat *matout) 381 { 382 Mat_Dense *mat = (Mat_Dense *) matin->data; 383 int k,j; 384 Scalar *v = mat->v, tmp; 385 386 if (!matout) { /* in place transpose */ 387 if (mat->m != mat->n) { 388 SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix"); 389 } 390 for ( j=0; j<mat->m; j++ ) { 391 for ( k=0; k<j; k++ ) { 392 tmp = v[j + k*mat->n]; 393 v[j + k*mat->n] = v[k + j*mat->n]; 394 v[k + j*mat->n] = tmp; 395 } 396 } 397 } 398 else { 399 SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet"); 400 } 401 return 0; 402 } 403 404 static int MatEqual_Dense(Mat matin1,Mat matin2) 405 { 406 Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 407 Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 408 int i; 409 Scalar *v1 = mat1->v, *v2 = mat2->v; 410 if (mat1->m != mat2->m) return 0; 411 if (mat1->n != mat2->n) return 0; 412 for ( i=0; i<mat1->m*mat1->n; i++ ) { 413 if (*v1 != *v2) return 0; 414 v1++; v2++; 415 } 416 return 1; 417 } 418 419 static int MatGetDiagonal_Dense(Mat matin,Vec v) 420 { 421 Mat_Dense *mat = (Mat_Dense *) matin->data; 422 int i, n; 423 Scalar *x; 424 VecGetArray(v,&x); VecGetSize(v,&n); 425 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec"); 426 for ( i=0; i<mat->m; i++ ) { 427 x[i] = mat->v[i*mat->m + i]; 428 } 429 return 0; 430 } 431 432 static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 433 { 434 Mat_Dense *mat = (Mat_Dense *) matin->data; 435 Scalar *l,*r,x,*v; 436 int i,j,m = mat->m, n = mat->n; 437 if (ll) { 438 VecGetArray(ll,&l); VecGetSize(ll,&m); 439 if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size"); 440 for ( i=0; i<m; i++ ) { 441 x = l[i]; 442 v = mat->v + i; 443 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 444 } 445 } 446 if (rr) { 447 VecGetArray(rr,&r); VecGetSize(rr,&n); 448 if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size"); 449 for ( i=0; i<n; i++ ) { 450 x = r[i]; 451 v = mat->v + i*m; 452 for ( j=0; j<m; j++ ) { (*v++) *= x;} 453 } 454 } 455 return 0; 456 } 457 458 static int MatNorm_Dense(Mat matin,MatNormType 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,"MatNorm_Dense:No 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 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