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