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