1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.98 1996/03/21 00:44:51 bsmith Exp bsmith $"; 3 #endif 4 /* 5 Defines the basic matrix operations for sequential dense. 6 */ 7 8 #include "dense.h" 9 #include "pinclude/plapack.h" 10 #include "pinclude/pviewer.h" 11 12 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 13 { 14 Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 15 int N = x->m*x->n, one = 1; 16 BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 17 PLogFlops(2*N-1); 18 return 0; 19 } 20 21 static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 22 { 23 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 24 int i,N = mat->m*mat->n,count = 0; 25 Scalar *v = mat->v; 26 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 27 if (nz) *nz = count; 28 if (nzalloc) *nzalloc = N; 29 if (mem) *mem = (int)A->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_SeqDense(Mat A,IS row,IS col,double f) 37 { 38 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 39 int info; 40 41 if (!mat->pivots) { 42 mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 43 PLogObjectMemory(A,mat->m*sizeof(int)); 44 } 45 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 46 if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 47 A->factor = FACTOR_LU; 48 PLogFlops((2*mat->n*mat->n*mat->n)/3); 49 return 0; 50 } 51 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 52 { 53 return MatConvert(A,MATSAME,fact); 54 } 55 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 56 { 57 return MatLUFactor(*fact,0,0,1.0); 58 } 59 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 60 { 61 return MatConvert(A,MATSAME,fact); 62 } 63 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 64 { 65 return MatCholeskyFactor(*fact,0,1.0); 66 } 67 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 68 { 69 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 70 int info; 71 72 if (mat->pivots) { 73 PetscFree(mat->pivots); 74 PLogObjectMemory(A,-mat->m*sizeof(int)); 75 mat->pivots = 0; 76 } 77 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 78 if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 79 A->factor = FACTOR_CHOLESKY; 80 PLogFlops((mat->n*mat->n*mat->n)/3); 81 return 0; 82 } 83 84 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 85 { 86 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 87 int one = 1, info, ierr; 88 Scalar *x, *y; 89 90 ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 91 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 92 if (A->factor == FACTOR_LU) { 93 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 94 } 95 else if (A->factor == FACTOR_CHOLESKY){ 96 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 97 } 98 else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 99 if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 100 PLogFlops(mat->n*mat->n - mat->n); 101 return 0; 102 } 103 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 104 { 105 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 106 int one = 1, info; 107 Scalar *x, *y; 108 109 VecGetArray(xx,&x); VecGetArray(yy,&y); 110 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 111 /* assume if pivots exist then use LU; else Cholesky */ 112 if (mat->pivots) { 113 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 114 } 115 else { 116 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 117 } 118 if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 119 PLogFlops(mat->n*mat->n - mat->n); 120 return 0; 121 } 122 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 123 { 124 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 125 int one = 1, info,ierr; 126 Scalar *x, *y, sone = 1.0; 127 Vec tmp = 0; 128 129 VecGetArray(xx,&x); VecGetArray(yy,&y); 130 if (yy == zz) { 131 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 132 PLogObjectParent(A,tmp); 133 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 134 } 135 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 136 /* assume if pivots exist then use LU; else Cholesky */ 137 if (mat->pivots) { 138 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 139 } 140 else { 141 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 142 } 143 if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 144 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 145 else VecAXPY(&sone,zz,yy); 146 PLogFlops(mat->n*mat->n - mat->n); 147 return 0; 148 } 149 150 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 151 { 152 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 153 int one = 1, info,ierr; 154 Scalar *x, *y, sone = 1.0; 155 Vec tmp; 156 157 VecGetArray(xx,&x); VecGetArray(yy,&y); 158 if (yy == zz) { 159 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 160 PLogObjectParent(A,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,y, &mat->m, &info ); 167 } 168 else { 169 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 170 } 171 if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 172 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 173 else VecAXPY(&sone,zz,yy); 174 PLogFlops(mat->n*mat->n - mat->n); 175 return 0; 176 } 177 /* ------------------------------------------------------------------*/ 178 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 179 double shift,int its,Vec xx) 180 { 181 Mat_SeqDense *mat = (Mat_SeqDense *) A->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 A,Vec xx,Vec yy) 232 { 233 Mat_SeqDense *mat = (Mat_SeqDense *) A->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 PLogFlops(2*mat->m*mat->n - mat->n); 239 return 0; 240 } 241 static int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 242 { 243 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 244 Scalar *v = mat->v, *x, *y; 245 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 246 VecGetArray(xx,&x); VecGetArray(yy,&y); 247 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 248 PLogFlops(2*mat->m*mat->n - mat->m); 249 return 0; 250 } 251 static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 252 { 253 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 254 Scalar *v = mat->v, *x, *y, *z; 255 int _One=1; Scalar _DOne=1.0; 256 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 257 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 258 LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 259 PLogFlops(2*mat->m*mat->n); 260 return 0; 261 } 262 static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 263 { 264 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 265 Scalar *v = mat->v, *x, *y, *z; 266 int _One=1; 267 Scalar _DOne=1.0; 268 VecGetArray(xx,&x); VecGetArray(yy,&y); 269 VecGetArray(zz,&z); 270 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 271 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 272 PLogFlops(2*mat->m*mat->n); 273 return 0; 274 } 275 276 /* -----------------------------------------------------------------*/ 277 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 278 { 279 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 280 Scalar *v; 281 int i; 282 283 *ncols = mat->n; 284 if (cols) { 285 *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols); 286 for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 287 } 288 if (vals) { 289 *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 290 v = mat->v + row; 291 for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 292 } 293 return 0; 294 } 295 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 296 { 297 if (cols) { PetscFree(*cols); } 298 if (vals) { PetscFree(*vals); } 299 return 0; 300 } 301 /* ----------------------------------------------------------------*/ 302 static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 303 int *indexn,Scalar *v,InsertMode addv) 304 { 305 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 306 int i,j; 307 308 if (!mat->roworiented) { 309 if (addv == INSERT_VALUES) { 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 else { 317 for ( j=0; j<n; j++ ) { 318 for ( i=0; i<m; i++ ) { 319 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 320 } 321 } 322 } 323 } 324 else { 325 if (addv == INSERT_VALUES) { 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 else { 333 for ( i=0; i<m; i++ ) { 334 for ( j=0; j<n; j++ ) { 335 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 336 } 337 } 338 } 339 } 340 return 0; 341 } 342 343 static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 344 { 345 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 346 int i, j; 347 Scalar *vpt = v; 348 349 /* row-oriented output */ 350 for ( i=0; i<m; i++ ) { 351 for ( j=0; j<n; j++ ) { 352 *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 353 } 354 } 355 return 0; 356 } 357 358 /* -----------------------------------------------------------------*/ 359 static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues) 360 { 361 Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 362 int ierr; 363 Mat newi; 364 365 ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr); 366 l = (Mat_SeqDense *) newi->data; 367 if (cpvalues == COPY_VALUES) { 368 PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 369 } 370 newi->assembled = PETSC_TRUE; 371 *newmat = newi; 372 return 0; 373 } 374 375 #include "sys.h" 376 377 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 378 { 379 Mat_SeqDense *a; 380 Mat B; 381 int *scols, i, j, nz, ierr, fd, header[4], size; 382 int *rowlengths = 0, M, N, *cols; 383 Scalar *vals, *svals, *v; 384 MPI_Comm comm = ((PetscObject)viewer)->comm; 385 386 MPI_Comm_size(comm,&size); 387 if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor"); 388 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 389 ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 390 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object"); 391 M = header[1]; N = header[2]; nz = header[3]; 392 393 if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 394 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 395 B = *A; 396 a = (Mat_SeqDense *) B->data; 397 398 /* read in nonzero values */ 399 ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr); 400 401 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 402 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 403 ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); 404 } else { 405 /* read row lengths */ 406 rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 407 ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 408 409 /* create our matrix */ 410 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 411 B = *A; 412 a = (Mat_SeqDense *) B->data; 413 v = a->v; 414 415 /* read column indices and nonzeros */ 416 cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols); 417 ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 418 vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 419 ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 420 421 /* insert into matrix */ 422 for ( i=0; i<M; i++ ) { 423 for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 424 svals += rowlengths[i]; scols += rowlengths[i]; 425 } 426 PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 427 428 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 429 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 430 } 431 return 0; 432 } 433 434 #include "pinclude/pviewer.h" 435 #include "sys.h" 436 437 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 438 { 439 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 440 int ierr, i, j, format; 441 FILE *fd; 442 char *outputname; 443 Scalar *v; 444 445 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 446 ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 447 ierr = ViewerGetFormat(viewer,&format); 448 if (format == ASCII_FORMAT_INFO) { 449 ; /* do nothing for now */ 450 } 451 else { 452 #if defined(PETSC_COMPLEX) 453 int allreal = 1; 454 /* determine if matrix has all real values */ 455 v = a->v; 456 for ( i=0; i<a->m*a->n; i++ ) { 457 if (imag(v[i])) { allreal = 0; break ;} 458 } 459 #endif 460 for ( i=0; i<a->m; i++ ) { 461 v = a->v + i; 462 for ( j=0; j<a->n; j++ ) { 463 #if defined(PETSC_COMPLEX) 464 if (allreal) { 465 fprintf(fd,"%6.4e ",real(*v)); v += a->m; 466 } else { 467 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 468 } 469 #else 470 fprintf(fd,"%6.4e ",*v); v += a->m; 471 #endif 472 } 473 fprintf(fd,"\n"); 474 } 475 } 476 fflush(fd); 477 return 0; 478 } 479 480 static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 481 { 482 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 483 int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 484 int format; 485 Scalar *v, *anonz,*vals; 486 487 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 488 489 ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 490 if (format == BINARY_FORMAT_NATIVE) { 491 /* store the matrix as a dense matrix */ 492 col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens); 493 col_lens[0] = MAT_COOKIE; 494 col_lens[1] = m; 495 col_lens[2] = n; 496 col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 497 ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr); 498 PetscFree(col_lens); 499 500 /* write out matrix, by rows */ 501 vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals); 502 v = a->v; 503 for ( i=0; i<m; i++ ) { 504 for ( j=0; j<n; j++ ) { 505 vals[i + j*m] = *v++; 506 } 507 } 508 ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr); 509 PetscFree(vals); 510 } else { 511 col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 512 col_lens[0] = MAT_COOKIE; 513 col_lens[1] = m; 514 col_lens[2] = n; 515 col_lens[3] = nz; 516 517 /* store lengths of each row and write (including header) to file */ 518 for ( i=0; i<m; i++ ) col_lens[4+i] = n; 519 ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr); 520 521 /* Possibly should write in smaller increments, not whole matrix at once? */ 522 /* store column indices (zero start index) */ 523 ict = 0; 524 for ( i=0; i<m; i++ ) { 525 for ( j=0; j<n; j++ ) col_lens[ict++] = j; 526 } 527 ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr); 528 PetscFree(col_lens); 529 530 /* store nonzero values */ 531 anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 532 ict = 0; 533 for ( i=0; i<m; i++ ) { 534 v = a->v + i; 535 for ( j=0; j<n; j++ ) { 536 anonz[ict++] = *v; v += a->m; 537 } 538 } 539 ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr); 540 PetscFree(anonz); 541 } 542 return 0; 543 } 544 545 static int MatView_SeqDense(PetscObject obj,Viewer viewer) 546 { 547 Mat A = (Mat) obj; 548 Mat_SeqDense *a = (Mat_SeqDense*) A->data; 549 ViewerType vtype; 550 int ierr; 551 552 if (!viewer) { 553 viewer = STDOUT_VIEWER_SELF; 554 } 555 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 556 557 if (vtype == MATLAB_VIEWER) { 558 return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 559 } 560 else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 561 return MatView_SeqDense_ASCII(A,viewer); 562 } 563 else if (vtype == BINARY_FILE_VIEWER) { 564 return MatView_SeqDense_Binary(A,viewer); 565 } 566 return 0; 567 } 568 569 static int MatDestroy_SeqDense(PetscObject obj) 570 { 571 Mat mat = (Mat) obj; 572 Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 573 #if defined(PETSC_LOG) 574 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 575 #endif 576 if (l->pivots) PetscFree(l->pivots); 577 if (!l->user_alloc) PetscFree(l->v); 578 PetscFree(l); 579 PLogObjectDestroy(mat); 580 PetscHeaderDestroy(mat); 581 return 0; 582 } 583 584 static int MatTranspose_SeqDense(Mat A,Mat *matout) 585 { 586 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 587 int k, j, m, n; 588 Scalar *v, tmp; 589 590 v = mat->v; m = mat->m; n = mat->n; 591 if (matout == PETSC_NULL) { /* in place transpose */ 592 if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 593 for ( j=0; j<m; j++ ) { 594 for ( k=0; k<j; k++ ) { 595 tmp = v[j + k*n]; 596 v[j + k*n] = v[k + j*n]; 597 v[k + j*n] = tmp; 598 } 599 } 600 } 601 else { /* out-of-place transpose */ 602 int ierr; 603 Mat tmat; 604 Mat_SeqDense *tmatd; 605 Scalar *v2; 606 ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 607 tmatd = (Mat_SeqDense *) tmat->data; 608 v = mat->v; v2 = tmatd->v; 609 for ( j=0; j<n; j++ ) { 610 for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 611 } 612 ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 613 ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 614 *matout = tmat; 615 } 616 return 0; 617 } 618 619 static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 620 { 621 Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 622 Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 623 int i; 624 Scalar *v1 = mat1->v, *v2 = mat2->v; 625 626 if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type MATSEQDENSE"); 627 if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;} 628 if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;} 629 for ( i=0; i<mat1->m*mat1->n; i++ ) { 630 if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;} 631 v1++; v2++; 632 } 633 *flg = PETSC_TRUE; 634 return 0; 635 } 636 637 static int MatGetDiagonal_SeqDense(Mat A,Vec v) 638 { 639 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 640 int i, n; 641 Scalar *x; 642 VecGetArray(v,&x); VecGetSize(v,&n); 643 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 644 for ( i=0; i<mat->m; i++ ) { 645 x[i] = mat->v[i*mat->m + i]; 646 } 647 return 0; 648 } 649 650 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 651 { 652 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 653 Scalar *l,*r,x,*v; 654 int i,j,m = mat->m, n = mat->n; 655 656 if (ll) { 657 VecGetArray(ll,&l); VecGetSize(ll,&m); 658 if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size"); 659 PLogFlops(n*m); 660 for ( i=0; i<m; i++ ) { 661 x = l[i]; 662 v = mat->v + i; 663 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 664 } 665 } 666 if (rr) { 667 VecGetArray(rr,&r); VecGetSize(rr,&n); 668 if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size"); 669 PLogFlops(n*m); 670 for ( i=0; i<n; i++ ) { 671 x = r[i]; 672 v = mat->v + i*m; 673 for ( j=0; j<m; j++ ) { (*v++) *= x;} 674 } 675 } 676 return 0; 677 } 678 679 static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 680 { 681 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 682 Scalar *v = mat->v; 683 double sum = 0.0; 684 int i, j; 685 686 if (type == NORM_FROBENIUS) { 687 for (i=0; i<mat->n*mat->m; i++ ) { 688 #if defined(PETSC_COMPLEX) 689 sum += real(conj(*v)*(*v)); v++; 690 #else 691 sum += (*v)*(*v); v++; 692 #endif 693 } 694 *norm = sqrt(sum); 695 PLogFlops(2*mat->n*mat->m); 696 } 697 else if (type == NORM_1) { 698 *norm = 0.0; 699 for ( j=0; j<mat->n; j++ ) { 700 sum = 0.0; 701 for ( i=0; i<mat->m; i++ ) { 702 sum += PetscAbsScalar(*v); v++; 703 } 704 if (sum > *norm) *norm = sum; 705 } 706 PLogFlops(mat->n*mat->m); 707 } 708 else if (type == NORM_INFINITY) { 709 *norm = 0.0; 710 for ( j=0; j<mat->m; j++ ) { 711 v = mat->v + j; 712 sum = 0.0; 713 for ( i=0; i<mat->n; i++ ) { 714 sum += PetscAbsScalar(*v); v += mat->m; 715 } 716 if (sum > *norm) *norm = sum; 717 } 718 PLogFlops(mat->n*mat->m); 719 } 720 else { 721 SETERRQ(1,"MatNorm_SeqDense:No two norm"); 722 } 723 return 0; 724 } 725 726 static int MatSetOption_SeqDense(Mat A,MatOption op) 727 { 728 Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 729 730 if (op == ROW_ORIENTED) aij->roworiented = 1; 731 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 732 else if (op == ROWS_SORTED || 733 op == COLUMNS_SORTED || 734 op == SYMMETRIC_MATRIX || 735 op == STRUCTURALLY_SYMMETRIC_MATRIX || 736 op == NO_NEW_NONZERO_LOCATIONS || 737 op == YES_NEW_NONZERO_LOCATIONS || 738 op == NO_NEW_DIAGONALS || 739 op == YES_NEW_DIAGONALS) 740 PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n"); 741 else 742 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 743 return 0; 744 } 745 746 static int MatZeroEntries_SeqDense(Mat A) 747 { 748 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 749 PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 750 return 0; 751 } 752 753 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 754 { 755 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 756 int n = l->n, i, j,ierr,N, *rows; 757 Scalar *slot; 758 759 ierr = ISGetSize(is,&N); CHKERRQ(ierr); 760 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 761 for ( i=0; i<N; i++ ) { 762 slot = l->v + rows[i]; 763 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 764 } 765 if (diag) { 766 for ( i=0; i<N; i++ ) { 767 slot = l->v + (n+1)*rows[i]; 768 *slot = *diag; 769 } 770 } 771 ISRestoreIndices(is,&rows); 772 return 0; 773 } 774 775 static int MatGetSize_SeqDense(Mat A,int *m,int *n) 776 { 777 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 778 *m = mat->m; *n = mat->n; 779 return 0; 780 } 781 782 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 783 { 784 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 785 *m = 0; *n = mat->m; 786 return 0; 787 } 788 789 static int MatGetArray_SeqDense(Mat A,Scalar **array) 790 { 791 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 792 *array = mat->v; 793 return 0; 794 } 795 796 static int MatRestoreArray_SeqDense(Mat A,Scalar **array) 797 { 798 return 0; 799 } 800 801 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 802 { 803 SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 804 } 805 806 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 807 Mat *submat) 808 { 809 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 810 int nznew, *smap, i, j, ierr, oldcols = mat->n; 811 int *irow, *icol, nrows, ncols, *cwork; 812 Scalar *vwork, *val; 813 Mat newmat; 814 815 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 816 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 817 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 818 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 819 820 smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 821 cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 822 vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 823 PetscMemzero((char*)smap,oldcols*sizeof(int)); 824 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 825 826 /* Create and fill new matrix */ 827 ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 828 for (i=0; i<nrows; i++) { 829 nznew = 0; 830 val = mat->v + irow[i]; 831 for (j=0; j<oldcols; j++) { 832 if (smap[j]) { 833 cwork[nznew] = smap[j] - 1; 834 vwork[nznew++] = val[j * mat->m]; 835 } 836 } 837 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 838 CHKERRQ(ierr); 839 } 840 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 841 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 842 843 /* Free work space */ 844 PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 845 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 846 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 847 *submat = newmat; 848 return 0; 849 } 850 851 static int MatCopy_SeqDense(Mat A, Mat B) 852 { 853 Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 854 if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B); 855 if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)"); 856 PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 857 return 0; 858 } 859 860 /* -------------------------------------------------------------------*/ 861 static struct _MatOps MatOps = {MatSetValues_SeqDense, 862 MatGetRow_SeqDense, MatRestoreRow_SeqDense, 863 MatMult_SeqDense, MatMultAdd_SeqDense, 864 MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 865 MatSolve_SeqDense,MatSolveAdd_SeqDense, 866 MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 867 MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 868 MatRelax_SeqDense, 869 MatTranspose_SeqDense, 870 MatGetInfo_SeqDense,MatEqual_SeqDense, 871 MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense, 872 0,0, 873 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 874 MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 875 MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 876 MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 877 0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0, 878 MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 879 MatConvertSameType_SeqDense,0,0,0,0, 880 MatAXPY_SeqDense,0,0, 881 MatGetValues_SeqDense, 882 MatCopy_SeqDense}; 883 884 885 /*@C 886 MatCreateSeqDense - Creates a sequential dense matrix that 887 is stored in column major order (the usual Fortran 77 manner). Many 888 of the matrix operations use the BLAS and LAPACK routines. 889 890 Input Parameters: 891 . comm - MPI communicator, set to MPI_COMM_SELF 892 . m - number of rows 893 . n - number of columns 894 . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 895 to control all matrix memory allocation. 896 897 Output Parameter: 898 . newmat - the matrix 899 900 Notes: 901 The data input variable is intended primarily for Fortran programmers 902 who wish to allocate their own matrix memory space. Most users should 903 set data=PETSC_NULL. 904 905 .keywords: dense, matrix, LAPACK, BLAS 906 907 .seealso: MatCreate(), MatSetValues() 908 @*/ 909 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat) 910 { 911 Mat mat; 912 Mat_SeqDense *l; 913 int ierr,flg; 914 915 *newmat = 0; 916 917 PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 918 PLogObjectCreate(mat); 919 l = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l); 920 PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 921 mat->destroy = MatDestroy_SeqDense; 922 mat->view = MatView_SeqDense; 923 mat->factor = 0; 924 PLogObjectMemory(mat,sizeof(struct _Mat)); 925 mat->data = (void *) l; 926 927 l->m = m; 928 l->n = n; 929 l->pivots = 0; 930 l->roworiented = 1; 931 if (data == PETSC_NULL) { 932 l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v); 933 PetscMemzero(l->v,m*n*sizeof(Scalar)); 934 l->user_alloc = 0; 935 PLogObjectMemory(mat,n*m*sizeof(Scalar)); 936 } 937 else { /* user-allocated storage */ 938 l->v = data; 939 l->user_alloc = 1; 940 } 941 942 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 943 if (flg) { 944 ierr = MatPrintHelp(mat); CHKERRQ(ierr); 945 } 946 *newmat = mat; 947 return 0; 948 } 949 950 int MatCreate_SeqDense(Mat A,Mat *newmat) 951 { 952 Mat_SeqDense *m = (Mat_SeqDense *) A->data; 953 return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat); 954 } 955