1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.83 1995/12/23 04:52:52 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 *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 MatConvertSameType_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,PETSC_NULL,&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,PETSC_NULL,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 if (!l->user_alloc) PetscFree(l->v); 529 PetscFree(l); 530 PLogObjectDestroy(mat); 531 PetscHeaderDestroy(mat); 532 return 0; 533 } 534 535 static int MatTranspose_SeqDense(Mat A,Mat *matout) 536 { 537 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 538 int k, j, m, n; 539 Scalar *v, tmp; 540 541 v = mat->v; m = mat->m; n = mat->n; 542 if (!matout) { /* in place transpose */ 543 if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 544 for ( j=0; j<m; j++ ) { 545 for ( k=0; k<j; k++ ) { 546 tmp = v[j + k*n]; 547 v[j + k*n] = v[k + j*n]; 548 v[k + j*n] = tmp; 549 } 550 } 551 } 552 else { /* out-of-place transpose */ 553 int ierr; 554 Mat tmat; 555 Mat_SeqDense *tmatd; 556 Scalar *v2; 557 ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 558 tmatd = (Mat_SeqDense *) tmat->data; 559 v = mat->v; v2 = tmatd->v; 560 for ( j=0; j<n; j++ ) { 561 for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 562 } 563 ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 564 ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 565 *matout = tmat; 566 } 567 return 0; 568 } 569 570 static int MatEqual_SeqDense(Mat A1,Mat A2) 571 { 572 Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 573 Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 574 int i; 575 Scalar *v1 = mat1->v, *v2 = mat2->v; 576 if (mat1->m != mat2->m) return 0; 577 if (mat1->n != mat2->n) return 0; 578 for ( i=0; i<mat1->m*mat1->n; i++ ) { 579 if (*v1 != *v2) return 0; 580 v1++; v2++; 581 } 582 return 1; 583 } 584 585 static int MatGetDiagonal_SeqDense(Mat A,Vec v) 586 { 587 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 588 int i, n; 589 Scalar *x; 590 VecGetArray(v,&x); VecGetSize(v,&n); 591 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 592 for ( i=0; i<mat->m; i++ ) { 593 x[i] = mat->v[i*mat->m + i]; 594 } 595 return 0; 596 } 597 598 static int MatScale_SeqDense(Mat A,Vec ll,Vec rr) 599 { 600 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 601 Scalar *l,*r,x,*v; 602 int i,j,m = mat->m, n = mat->n; 603 604 if (ll) { 605 VecGetArray(ll,&l); VecGetSize(ll,&m); 606 if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size"); 607 PLogFlops(n*m); 608 for ( i=0; i<m; i++ ) { 609 x = l[i]; 610 v = mat->v + i; 611 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 612 } 613 } 614 if (rr) { 615 VecGetArray(rr,&r); VecGetSize(rr,&n); 616 if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size"); 617 PLogFlops(n*m); 618 for ( i=0; i<n; i++ ) { 619 x = r[i]; 620 v = mat->v + i*m; 621 for ( j=0; j<m; j++ ) { (*v++) *= x;} 622 } 623 } 624 return 0; 625 } 626 627 static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 628 { 629 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 630 Scalar *v = mat->v; 631 double sum = 0.0; 632 int i, j; 633 634 if (type == NORM_FROBENIUS) { 635 for (i=0; i<mat->n*mat->m; i++ ) { 636 #if defined(PETSC_COMPLEX) 637 sum += real(conj(*v)*(*v)); v++; 638 #else 639 sum += (*v)*(*v); v++; 640 #endif 641 } 642 *norm = sqrt(sum); 643 PLogFlops(2*mat->n*mat->m); 644 } 645 else if (type == NORM_1) { 646 *norm = 0.0; 647 for ( j=0; j<mat->n; j++ ) { 648 sum = 0.0; 649 for ( i=0; i<mat->m; i++ ) { 650 sum += PetscAbsScalar(*v); v++; 651 } 652 if (sum > *norm) *norm = sum; 653 } 654 PLogFlops(mat->n*mat->m); 655 } 656 else if (type == NORM_INFINITY) { 657 *norm = 0.0; 658 for ( j=0; j<mat->m; j++ ) { 659 v = mat->v + j; 660 sum = 0.0; 661 for ( i=0; i<mat->n; i++ ) { 662 sum += PetscAbsScalar(*v); v += mat->m; 663 } 664 if (sum > *norm) *norm = sum; 665 } 666 PLogFlops(mat->n*mat->m); 667 } 668 else { 669 SETERRQ(1,"MatNorm_SeqDense:No two norm"); 670 } 671 return 0; 672 } 673 674 static int MatSetOption_SeqDense(Mat A,MatOption op) 675 { 676 Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 677 678 if (op == ROW_ORIENTED) aij->roworiented = 1; 679 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 680 else if (op == ROWS_SORTED || 681 op == COLUMNS_SORTED || 682 op == SYMMETRIC_MATRIX || 683 op == STRUCTURALLY_SYMMETRIC_MATRIX || 684 op == NO_NEW_NONZERO_LOCATIONS || 685 op == YES_NEW_NONZERO_LOCATIONS || 686 op == NO_NEW_DIAGONALS || 687 op == YES_NEW_DIAGONALS) 688 PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n"); 689 else 690 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 691 return 0; 692 } 693 694 static int MatZeroEntries_SeqDense(Mat A) 695 { 696 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 697 PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 698 return 0; 699 } 700 701 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 702 { 703 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 704 int n = l->n, i, j,ierr,N, *rows; 705 Scalar *slot; 706 707 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 708 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 709 for ( i=0; i<N; i++ ) { 710 slot = l->v + rows[i]; 711 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 712 } 713 if (diag) { 714 for ( i=0; i<N; i++ ) { 715 slot = l->v + (n+1)*rows[i]; 716 *slot = *diag; 717 } 718 } 719 ISRestoreIndices(is,&rows); 720 return 0; 721 } 722 723 static int MatGetSize_SeqDense(Mat A,int *m,int *n) 724 { 725 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 726 *m = mat->m; *n = mat->n; 727 return 0; 728 } 729 730 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 731 { 732 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 733 *m = 0; *n = mat->m; 734 return 0; 735 } 736 737 static int MatGetArray_SeqDense(Mat A,Scalar **array) 738 { 739 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 740 *array = mat->v; 741 return 0; 742 } 743 744 745 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 746 { 747 SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 748 } 749 750 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 751 Mat *submat) 752 { 753 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 754 int nznew, *smap, i, j, ierr, oldcols = mat->n; 755 int *irow, *icol, nrows, ncols, *cwork; 756 Scalar *vwork, *val; 757 Mat newmat; 758 759 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 760 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 761 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 762 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 763 764 smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 765 cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 766 vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 767 PetscMemzero((char*)smap,oldcols*sizeof(int)); 768 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 769 770 /* Create and fill new matrix */ 771 ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 772 for (i=0; i<nrows; i++) { 773 nznew = 0; 774 val = mat->v + irow[i]; 775 for (j=0; j<oldcols; j++) { 776 if (smap[j]) { 777 cwork[nznew] = smap[j] - 1; 778 vwork[nznew++] = val[j * mat->m]; 779 } 780 } 781 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 782 CHKERRQ(ierr); 783 } 784 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 785 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 786 787 /* Free work space */ 788 PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 789 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 790 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 791 *submat = newmat; 792 return 0; 793 } 794 795 static int MatCopy_SeqDense(Mat A, Mat B) 796 { 797 Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 798 if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B); 799 if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)"); 800 PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 801 return 0; 802 } 803 804 /* -------------------------------------------------------------------*/ 805 static struct _MatOps MatOps = {MatSetValues_SeqDense, 806 MatGetRow_SeqDense, MatRestoreRow_SeqDense, 807 MatMult_SeqDense, MatMultAdd_SeqDense, 808 MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 809 MatSolve_SeqDense,MatSolveAdd_SeqDense, 810 MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 811 MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 812 MatRelax_SeqDense, 813 MatTranspose_SeqDense, 814 MatGetInfo_SeqDense,MatEqual_SeqDense, 815 MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense, 816 0,0, 817 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 818 MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 819 MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 820 MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 821 0,0,MatGetArray_SeqDense,0,0, 822 MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 823 MatConvertSameType_SeqDense,0,0,0,0, 824 MatAXPY_SeqDense,0,0, 825 MatGetValues_SeqDense, 826 MatCopy_SeqDense}; 827 828 /*@C 829 MatCreateSeqDense - Creates a sequential dense matrix that 830 is stored in column major order (the usual Fortran 77 manner). Many 831 of the matrix operations use the BLAS and LAPACK routines. 832 833 Input Parameters: 834 . comm - MPI communicator, set to MPI_COMM_SELF 835 . m - number of rows 836 . n - number of columns 837 . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 838 to control all matrix memory allocation. 839 840 Output Parameter: 841 . newmat - the matrix 842 843 Notes: 844 The data input variable is intended primarily for Fortran programmers 845 who wish to allocate their own matrix memory space. Most users should 846 set data=PETSC_NULL. 847 848 .keywords: dense, matrix, LAPACK, BLAS 849 850 .seealso: MatCreate(), MatSetValues() 851 @*/ 852 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat) 853 { 854 Mat mat; 855 Mat_SeqDense *l; 856 857 *newmat = 0; 858 859 PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 860 PLogObjectCreate(mat); 861 l = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l); 862 PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 863 mat->destroy = MatDestroy_SeqDense; 864 mat->view = MatView_SeqDense; 865 mat->factor = 0; 866 PLogObjectMemory(mat,sizeof(struct _Mat)); 867 mat->data = (void *) l; 868 869 l->m = m; 870 l->n = n; 871 l->pivots = 0; 872 l->roworiented = 1; 873 if (data == PETSC_NULL) { 874 l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v); 875 PetscMemzero(l->v,m*n*sizeof(Scalar)); 876 l->user_alloc = 0; 877 PLogObjectMemory(mat,n*m*sizeof(Scalar)); 878 } 879 else { /* user-allocated storage */ 880 l->v = data; 881 l->user_alloc = 1; 882 } 883 884 *newmat = mat; 885 return 0; 886 } 887 888 int MatCreate_SeqDense(Mat A,Mat *newmat) 889 { 890 Mat_SeqDense *m = (Mat_SeqDense *) A->data; 891 return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat); 892 } 893