1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.71 1995/11/01 19:10:00 bsmith Exp bsmith $"; 3 #endif 4 5 #include "dense.h" 6 #include "pinclude/plapack.h" 7 #include "pinclude/pviewer.h" 8 9 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 10 { 11 Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 12 int N = x->m*x->n, one = 1; 13 BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 14 PLogFlops(2*N-1); 15 return 0; 16 } 17 18 static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 19 { 20 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 21 int i,N = mat->m*mat->n,count = 0; 22 Scalar *v = mat->v; 23 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 24 *nz = count; *nzalloc = N; *mem = (int)A->mem; 25 return 0; 26 } 27 28 /* ---------------------------------------------------------------*/ 29 /* COMMENT: I have chosen to hide column permutation in the pivots, 30 rather than put it in the Mat->col slot.*/ 31 static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 32 { 33 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 34 int info; 35 36 if (!mat->pivots) { 37 mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 38 PLogObjectMemory(A,mat->m*sizeof(int)); 39 } 40 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 41 if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 42 A->factor = FACTOR_LU; 43 PLogFlops((2*mat->n*mat->n*mat->n)/3); 44 return 0; 45 } 46 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 47 { 48 int ierr; 49 ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 50 return 0; 51 } 52 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 53 { 54 return MatLUFactor(*fact,0,0,1.0); 55 } 56 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 57 { 58 int ierr; 59 ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 60 return 0; 61 } 62 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 63 { 64 return MatCholeskyFactor(*fact,0,1.0); 65 } 66 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 67 { 68 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 69 int info; 70 71 if (mat->pivots) { 72 PetscFree(mat->pivots); 73 PLogObjectMemory(A,-mat->m*sizeof(int)); 74 mat->pivots = 0; 75 } 76 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 77 if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 78 A->factor = FACTOR_CHOLESKY; 79 PLogFlops((mat->n*mat->n*mat->n)/3); 80 return 0; 81 } 82 83 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 84 { 85 Mat_SeqDense *mat = (Mat_SeqDense *) A->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 (A->factor == FACTOR_LU) { 91 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 92 } 93 else if (A->factor == FACTOR_CHOLESKY){ 94 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 95 } 96 else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 97 if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 98 PLogFlops(mat->n*mat->n - mat->n); 99 return 0; 100 } 101 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 102 { 103 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 104 int one = 1, info; 105 Scalar *x, *y; 106 VecGetArray(xx,&x); VecGetArray(yy,&y); 107 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 108 /* assume if pivots exist then use LU; else Cholesky */ 109 if (mat->pivots) { 110 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 111 } 112 else { 113 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 114 } 115 if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 116 PLogFlops(mat->n*mat->n - mat->n); 117 return 0; 118 } 119 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 120 { 121 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 122 int one = 1, info,ierr; 123 Scalar *x, *y, sone = 1.0; 124 Vec tmp = 0; 125 VecGetArray(xx,&x); VecGetArray(yy,&y); 126 if (yy == zz) { 127 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 128 PLogObjectParent(A,tmp); 129 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 130 } 131 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 132 /* assume if pivots exist then use LU; else Cholesky */ 133 if (mat->pivots) { 134 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 135 } 136 else { 137 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 138 } 139 if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 140 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 141 else VecAXPY(&sone,zz,yy); 142 PLogFlops(mat->n*mat->n - mat->n); 143 return 0; 144 } 145 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 146 { 147 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 148 int one = 1, info,ierr; 149 Scalar *x, *y, sone = 1.0; 150 Vec tmp; 151 VecGetArray(xx,&x); VecGetArray(yy,&y); 152 if (yy == zz) { 153 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 154 PLogObjectParent(A,tmp); 155 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 156 } 157 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 158 /* assume if pivots exist then use LU; else Cholesky */ 159 if (mat->pivots) { 160 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 161 } 162 else { 163 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 164 } 165 if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 166 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 167 else VecAXPY(&sone,zz,yy); 168 PLogFlops(mat->n*mat->n - mat->n); 169 return 0; 170 } 171 /* ------------------------------------------------------------------*/ 172 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 173 double shift,int its,Vec xx) 174 { 175 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 176 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 177 int o = 1,ierr, m = mat->m, i; 178 179 if (flag & SOR_ZERO_INITIAL_GUESS) { 180 /* this is a hack fix, should have another version without 181 the second BLdot */ 182 ierr = VecSet(&zero,xx); CHKERRQ(ierr); 183 } 184 VecGetArray(xx,&x); VecGetArray(bb,&b); 185 while (its--) { 186 if (flag & SOR_FORWARD_SWEEP){ 187 for ( i=0; i<m; i++ ) { 188 #if defined(PETSC_COMPLEX) 189 /* cannot use BLAS dot for complex because compiler/linker is 190 not happy about returning a double complex */ 191 int _i; 192 Scalar sum = b[i]; 193 for ( _i=0; _i<m; _i++ ) { 194 sum -= conj(v[i+_i*m])*x[_i]; 195 } 196 xt = sum; 197 #else 198 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 199 #endif 200 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 201 } 202 } 203 if (flag & SOR_BACKWARD_SWEEP) { 204 for ( i=m-1; i>=0; i-- ) { 205 #if defined(PETSC_COMPLEX) 206 /* cannot use BLAS dot for complex because compiler/linker is 207 not happy about returning a double complex */ 208 int _i; 209 Scalar sum = b[i]; 210 for ( _i=0; _i<m; _i++ ) { 211 sum -= conj(v[i+_i*m])*x[_i]; 212 } 213 xt = sum; 214 #else 215 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 216 #endif 217 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 218 } 219 } 220 } 221 return 0; 222 } 223 224 /* -----------------------------------------------------------------*/ 225 static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 226 { 227 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 228 Scalar *v = mat->v, *x, *y; 229 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 230 VecGetArray(xx,&x), VecGetArray(yy,&y); 231 LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 232 PLogFlops(2*mat->m*mat->n - mat->n); 233 return 0; 234 } 235 static int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 236 { 237 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 238 Scalar *v = mat->v, *x, *y; 239 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 240 VecGetArray(xx,&x); VecGetArray(yy,&y); 241 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 242 PLogFlops(2*mat->m*mat->n - mat->m); 243 return 0; 244 } 245 static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 246 { 247 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 248 Scalar *v = mat->v, *x, *y, *z; 249 int _One=1; Scalar _DOne=1.0; 250 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 251 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 252 LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 253 PLogFlops(2*mat->m*mat->n); 254 return 0; 255 } 256 static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 257 { 258 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 259 Scalar *v = mat->v, *x, *y, *z; 260 int _One=1; 261 Scalar _DOne=1.0; 262 VecGetArray(xx,&x); VecGetArray(yy,&y); 263 VecGetArray(zz,&z); 264 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 265 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 266 PLogFlops(2*mat->m*mat->n); 267 return 0; 268 } 269 270 /* -----------------------------------------------------------------*/ 271 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 272 { 273 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 274 Scalar *v; 275 int i; 276 *ncols = mat->n; 277 if (cols) { 278 *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols); 279 for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 280 } 281 if (vals) { 282 *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 283 v = mat->v + row; 284 for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 285 } 286 return 0; 287 } 288 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 289 { 290 if (cols) { PetscFree(*cols); } 291 if (vals) { PetscFree(*vals); } 292 return 0; 293 } 294 /* ----------------------------------------------------------------*/ 295 static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n, 296 int *indexn,Scalar *v,InsertMode addv) 297 { 298 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 299 int i,j; 300 301 if (!mat->roworiented) { 302 if (addv == INSERT_VALUES) { 303 for ( j=0; j<n; j++ ) { 304 for ( i=0; i<m; i++ ) { 305 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 306 } 307 } 308 } 309 else { 310 for ( j=0; j<n; j++ ) { 311 for ( i=0; i<m; i++ ) { 312 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 313 } 314 } 315 } 316 } 317 else { 318 if (addv == INSERT_VALUES) { 319 for ( i=0; i<m; i++ ) { 320 for ( j=0; j<n; j++ ) { 321 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 322 } 323 } 324 } 325 else { 326 for ( i=0; i<m; i++ ) { 327 for ( j=0; j<n; j++ ) { 328 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 329 } 330 } 331 } 332 } 333 return 0; 334 } 335 336 /* -----------------------------------------------------------------*/ 337 static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat,int cpvalues) 338 { 339 Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 340 int ierr; 341 Mat newi; 342 343 ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,&newi); CHKERRQ(ierr); 344 l = (Mat_SeqDense *) newi->data; 345 if (cpvalues == COPY_VALUES) { 346 PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 347 } 348 *newmat = newi; 349 return 0; 350 } 351 352 #include "sysio.h" 353 354 int MatLoad_SeqDense(Viewer bview,MatType type,Mat *A) 355 { 356 Mat_SeqDense *a; 357 Mat B; 358 int *scols, i, j, nz, ierr, fd, header[4], size; 359 int *rowlengths = 0, M, N, *cols; 360 Scalar *vals, *svals, *v; 361 PetscObject vobj = (PetscObject) bview; 362 MPI_Comm comm = vobj->comm; 363 364 MPI_Comm_size(comm,&size); 365 if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor"); 366 ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 367 ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 368 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object"); 369 M = header[1]; N = header[2]; nz = header[3]; 370 371 /* read row lengths */ 372 rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 373 ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 374 375 /* create our matrix */ 376 ierr = MatCreateSeqDense(comm,M,N,A); CHKERRQ(ierr); 377 B = *A; 378 a = (Mat_SeqDense *) B->data; 379 v = a->v; 380 381 /* read column indices and nonzeros */ 382 cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols); 383 ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr); 384 vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 385 ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 386 387 /* insert into matrix */ 388 for ( i=0; i<M; i++ ) { 389 for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 390 svals += rowlengths[i]; scols += rowlengths[i]; 391 } 392 PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 393 394 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 395 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 396 return 0; 397 } 398 399 #include "pinclude/pviewer.h" 400 #include "sysio.h" 401 402 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 403 { 404 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 405 int ierr, i, j, format; 406 FILE *fd; 407 char *outputname; 408 Scalar *v; 409 410 ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 411 ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 412 ierr = ViewerFileGetFormat_Private(viewer,&format); 413 if (format == FILE_FORMAT_INFO) { 414 ; /* do nothing for now */ 415 } 416 else { 417 for ( i=0; i<a->m; i++ ) { 418 v = a->v + i; 419 for ( j=0; j<a->n; j++ ) { 420 #if defined(PETSC_COMPLEX) 421 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 422 #else 423 fprintf(fd,"%6.4e ",*v); v += a->m; 424 #endif 425 } 426 fprintf(fd,"\n"); 427 } 428 } 429 fflush(fd); 430 return 0; 431 } 432 433 static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 434 { 435 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 436 int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 437 Scalar *v, *anonz; 438 439 ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 440 col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 441 col_lens[0] = MAT_COOKIE; 442 col_lens[1] = m; 443 col_lens[2] = n; 444 col_lens[3] = nz; 445 446 /* store lengths of each row and write (including header) to file */ 447 for ( i=0; i<m; i++ ) col_lens[4+i] = n; 448 ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr); 449 450 /* Possibly should write in smaller increments, not whole matrix at once? */ 451 /* store column indices (zero start index) */ 452 ict = 0; 453 for ( i=0; i<m; i++ ) { 454 for ( j=0; j<n; j++ ) col_lens[ict++] = j; 455 } 456 ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr); 457 PetscFree(col_lens); 458 459 /* store nonzero values */ 460 anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 461 ict = 0; 462 for ( i=0; i<m; i++ ) { 463 v = a->v + i; 464 for ( j=0; j<n; j++ ) { 465 anonz[ict++] = *v; v += a->m; 466 } 467 } 468 ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr); 469 PetscFree(anonz); 470 return 0; 471 } 472 473 static int MatView_SeqDense(PetscObject obj,Viewer viewer) 474 { 475 Mat A = (Mat) obj; 476 Mat_SeqDense *a = (Mat_SeqDense*) A->data; 477 PetscObject vobj = (PetscObject) viewer; 478 479 if (!viewer) { 480 viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 481 } 482 if (vobj->cookie == VIEWER_COOKIE) { 483 if (vobj->type == MATLAB_VIEWER) { 484 return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 485 } 486 else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 487 return MatView_SeqDense_ASCII(A,viewer); 488 } 489 else if (vobj->type == BINARY_FILE_VIEWER) { 490 return MatView_SeqDense_Binary(A,viewer); 491 } 492 } 493 return 0; 494 } 495 496 static int MatDestroy_SeqDense(PetscObject obj) 497 { 498 Mat mat = (Mat) obj; 499 Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 500 #if defined(PETSC_LOG) 501 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 502 #endif 503 if (l->pivots) PetscFree(l->pivots); 504 PetscFree(l); 505 PLogObjectDestroy(mat); 506 PetscHeaderDestroy(mat); 507 return 0; 508 } 509 510 static int MatTranspose_SeqDense(Mat A,Mat *matout) 511 { 512 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 513 int k, j, m, n; 514 Scalar *v, tmp; 515 516 v = mat->v; m = mat->m; n = mat->n; 517 if (!matout) { /* in place transpose */ 518 if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 519 for ( j=0; j<m; j++ ) { 520 for ( k=0; k<j; k++ ) { 521 tmp = v[j + k*n]; 522 v[j + k*n] = v[k + j*n]; 523 v[k + j*n] = tmp; 524 } 525 } 526 } 527 else { /* out-of-place transpose */ 528 int ierr; 529 Mat tmat; 530 Mat_SeqDense *tmatd; 531 Scalar *v2; 532 ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr); 533 tmatd = (Mat_SeqDense *) tmat->data; 534 v = mat->v; v2 = tmatd->v; 535 for ( j=0; j<n; j++ ) { 536 for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 537 } 538 ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 539 ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 540 *matout = tmat; 541 } 542 return 0; 543 } 544 545 static int MatEqual_SeqDense(Mat A1,Mat A2) 546 { 547 Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 548 Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 549 int i; 550 Scalar *v1 = mat1->v, *v2 = mat2->v; 551 if (mat1->m != mat2->m) return 0; 552 if (mat1->n != mat2->n) return 0; 553 for ( i=0; i<mat1->m*mat1->n; i++ ) { 554 if (*v1 != *v2) return 0; 555 v1++; v2++; 556 } 557 return 1; 558 } 559 560 static int MatGetDiagonal_SeqDense(Mat A,Vec v) 561 { 562 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 563 int i, n; 564 Scalar *x; 565 VecGetArray(v,&x); VecGetSize(v,&n); 566 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 567 for ( i=0; i<mat->m; i++ ) { 568 x[i] = mat->v[i*mat->m + i]; 569 } 570 return 0; 571 } 572 573 static int MatScale_SeqDense(Mat A,Vec ll,Vec rr) 574 { 575 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 576 Scalar *l,*r,x,*v; 577 int i,j,m = mat->m, n = mat->n; 578 579 if (ll) { 580 VecGetArray(ll,&l); VecGetSize(ll,&m); 581 if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size"); 582 PLogFlops(n*m); 583 for ( i=0; i<m; i++ ) { 584 x = l[i]; 585 v = mat->v + i; 586 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 587 } 588 } 589 if (rr) { 590 VecGetArray(rr,&r); VecGetSize(rr,&n); 591 if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size"); 592 PLogFlops(n*m); 593 for ( i=0; i<n; i++ ) { 594 x = r[i]; 595 v = mat->v + i*m; 596 for ( j=0; j<m; j++ ) { (*v++) *= x;} 597 } 598 } 599 return 0; 600 } 601 602 static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 603 { 604 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 605 Scalar *v = mat->v; 606 double sum = 0.0; 607 int i, j; 608 609 if (type == NORM_FROBENIUS) { 610 for (i=0; i<mat->n*mat->m; i++ ) { 611 #if defined(PETSC_COMPLEX) 612 sum += real(conj(*v)*(*v)); v++; 613 #else 614 sum += (*v)*(*v); v++; 615 #endif 616 } 617 *norm = sqrt(sum); 618 PLogFlops(2*mat->n*mat->m); 619 } 620 else if (type == NORM_1) { 621 *norm = 0.0; 622 for ( j=0; j<mat->n; j++ ) { 623 sum = 0.0; 624 for ( i=0; i<mat->m; i++ ) { 625 sum += PetscAbsScalar(*v++); 626 } 627 if (sum > *norm) *norm = sum; 628 } 629 PLogFlops(mat->n*mat->m); 630 } 631 else if (type == NORM_INFINITY) { 632 *norm = 0.0; 633 for ( j=0; j<mat->m; j++ ) { 634 v = mat->v + j; 635 sum = 0.0; 636 for ( i=0; i<mat->n; i++ ) { 637 sum += PetscAbsScalar(*v); v += mat->m; 638 } 639 if (sum > *norm) *norm = sum; 640 } 641 PLogFlops(mat->n*mat->m); 642 } 643 else { 644 SETERRQ(1,"MatNorm_SeqDense:No two norm"); 645 } 646 return 0; 647 } 648 649 static int MatSetOption_SeqDense(Mat A,MatOption op) 650 { 651 Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 652 if (op == ROW_ORIENTED) aij->roworiented = 1; 653 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 654 else if (op == ROWS_SORTED || 655 op == COLUMNS_SORTED || 656 op == SYMMETRIC_MATRIX || 657 op == STRUCTURALLY_SYMMETRIC_MATRIX || 658 op == NO_NEW_NONZERO_LOCATIONS || 659 op == YES_NEW_NONZERO_LOCATIONS || 660 op == NO_NEW_DIAGONALS || 661 op == YES_NEW_DIAGONALS) 662 PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n"); 663 else 664 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 665 return 0; 666 } 667 668 static int MatZeroEntries_SeqDense(Mat A) 669 { 670 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 671 PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 672 return 0; 673 } 674 675 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 676 { 677 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 678 int n = l->n, i, j,ierr,N, *rows; 679 Scalar *slot; 680 681 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 682 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 683 for ( i=0; i<N; i++ ) { 684 slot = l->v + rows[i]; 685 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 686 } 687 if (diag) { 688 for ( i=0; i<N; i++ ) { 689 slot = l->v + (n+1)*rows[i]; 690 *slot = *diag; 691 } 692 } 693 ISRestoreIndices(is,&rows); 694 return 0; 695 } 696 697 static int MatGetSize_SeqDense(Mat A,int *m,int *n) 698 { 699 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 700 *m = mat->m; *n = mat->n; 701 return 0; 702 } 703 704 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 705 { 706 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 707 *m = 0; *n = mat->m; 708 return 0; 709 } 710 711 static int MatGetArray_SeqDense(Mat A,Scalar **array) 712 { 713 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 714 *array = mat->v; 715 return 0; 716 } 717 718 719 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 720 { 721 SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 722 } 723 724 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 725 Mat *submat) 726 { 727 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 728 int nznew, *smap, i, j, ierr, oldcols = mat->n; 729 int *irow, *icol, nrows, ncols, *cwork; 730 Scalar *vwork, *val; 731 Mat newmat; 732 733 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 734 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 735 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 736 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 737 738 smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 739 cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 740 vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 741 PetscMemzero((char*)smap,oldcols*sizeof(int)); 742 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 743 744 /* Create and fill new matrix */ 745 ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr); 746 for (i=0; i<nrows; i++) { 747 nznew = 0; 748 val = mat->v + irow[i]; 749 for (j=0; j<oldcols; j++) { 750 if (smap[j]) { 751 cwork[nznew] = smap[j] - 1; 752 vwork[nznew++] = val[j * mat->m]; 753 } 754 } 755 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 756 CHKERRQ(ierr); 757 } 758 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 759 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 760 761 /* Free work space */ 762 PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 763 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 764 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 765 *submat = newmat; 766 return 0; 767 } 768 769 /* -------------------------------------------------------------------*/ 770 static struct _MatOps MatOps = {MatInsert_SeqDense, 771 MatGetRow_SeqDense, MatRestoreRow_SeqDense, 772 MatMult_SeqDense, MatMultAdd_SeqDense, 773 MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 774 MatSolve_SeqDense,MatSolveAdd_SeqDense, 775 MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 776 MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 777 MatRelax_SeqDense, 778 MatTranspose_SeqDense, 779 MatGetInfo_SeqDense,MatEqual_SeqDense, 780 MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense, 781 0,0, 782 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 783 MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 784 MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 785 MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 786 0,0,MatGetArray_SeqDense,0,0, 787 MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 788 MatCopyPrivate_SeqDense,0,0,0,0, 789 MatAXPY_SeqDense}; 790 791 /*@C 792 MatCreateSeqDense - Creates a sequential dense matrix that 793 is stored in column major order (the usual Fortran 77 manner). Many 794 of the matrix operations use the BLAS and LAPACK routines. 795 796 Input Parameters: 797 . comm - MPI communicator, set to MPI_COMM_SELF 798 . m - number of rows 799 . n - number of column 800 801 Output Parameter: 802 . newmat - the matrix 803 804 .keywords: dense, matrix, LAPACK, BLAS 805 806 .seealso: MatCreate(), MatSetValues() 807 @*/ 808 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat) 809 { 810 int size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar); 811 Mat mat; 812 Mat_SeqDense *l; 813 814 *newmat = 0; 815 PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 816 PLogObjectCreate(mat); 817 l = (Mat_SeqDense *) PetscMalloc(size); CHKPTRQ(l); 818 PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 819 mat->destroy = MatDestroy_SeqDense; 820 mat->view = MatView_SeqDense; 821 mat->data = (void *) l; 822 mat->factor = 0; 823 PLogObjectMemory(mat,sizeof(struct _Mat) + size); 824 825 l->m = m; 826 l->n = n; 827 l->v = (Scalar *) (l + 1); 828 l->pivots = 0; 829 l->roworiented = 1; 830 831 PetscMemzero(l->v,m*n*sizeof(Scalar)); 832 *newmat = mat; 833 return 0; 834 } 835 836 int MatCreate_SeqDense(Mat A,Mat *newmat) 837 { 838 Mat_SeqDense *m = (Mat_SeqDense *) A->data; 839 return MatCreateSeqDense(A->comm,m->m,m->n,newmat); 840 } 841