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