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