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