1 2 #ifndef lint 3 static char vcid[] = "$Id: baij.c,v 1.48 1996/05/17 21:49:04 curfman Exp bsmith $"; 4 #endif 5 6 /* 7 Defines the basic matrix operations for the BAIJ (compressed row) 8 matrix storage format. 9 */ 10 #include "baij.h" 11 #include "src/vec/vecimpl.h" 12 #include "src/inline/spops.h" 13 #include "petsc.h" 14 15 extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 16 17 static int MatGetReordering_SeqBAIJ(Mat A,MatOrdering type,IS *rperm,IS *cperm) 18 { 19 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 20 int ierr, *ia, *ja,n = a->mbs,*idx,i,ishift,oshift; 21 22 /* 23 this is tacky: In the future when we have written special factorization 24 and solve routines for the identity permutation we should use a 25 stride index set instead of the general one. 26 */ 27 if (type == ORDER_NATURAL) { 28 idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 29 for ( i=0; i<n; i++ ) idx[i] = i; 30 ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 31 ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 32 PetscFree(idx); 33 ISSetPermutation(*rperm); 34 ISSetPermutation(*cperm); 35 ISSetIdentity(*rperm); 36 ISSetIdentity(*cperm); 37 return 0; 38 } 39 40 MatReorderingRegisterAll(); 41 ishift = 0; 42 oshift = -MatReorderingIndexShift[(int)type]; 43 if (MatReorderingRequiresSymmetric[(int)type]) { 44 ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,ishift,oshift,&ia,&ja);CHKERRQ(ierr); 45 ierr = MatGetReordering_IJ(n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 46 PetscFree(ia); PetscFree(ja); 47 } else { 48 if (ishift == oshift) { 49 ierr = MatGetReordering_IJ(n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 50 } 51 else if (ishift == -1) { 52 /* temporarily subtract 1 from i and j indices */ 53 int nz = a->i[n] - 1; 54 for ( i=0; i<nz; i++ ) a->j[i]--; 55 for ( i=0; i<n+1; i++ ) a->i[i]--; 56 ierr = MatGetReordering_IJ(n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 57 for ( i=0; i<nz; i++ ) a->j[i]++; 58 for ( i=0; i<n+1; i++ ) a->i[i]++; 59 } else { 60 /* temporarily add 1 to i and j indices */ 61 int nz = a->i[n] - 1; 62 for ( i=0; i<nz; i++ ) a->j[i]++; 63 for ( i=0; i<n+1; i++ ) a->i[i]++; 64 ierr = MatGetReordering_IJ(n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 65 for ( i=0; i<nz; i++ ) a->j[i]--; 66 for ( i=0; i<n+1; i++ ) a->i[i]--; 67 } 68 } 69 return 0; 70 } 71 72 /* 73 Adds diagonal pointers to sparse matrix structure. 74 */ 75 76 int MatMarkDiag_SeqBAIJ(Mat A) 77 { 78 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 79 int i,j, *diag, m = a->mbs; 80 81 diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 82 PLogObjectMemory(A,(m+1)*sizeof(int)); 83 for ( i=0; i<m; i++ ) { 84 for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 85 if (a->j[j] == i) { 86 diag[i] = j; 87 break; 88 } 89 } 90 } 91 a->diag = diag; 92 return 0; 93 } 94 95 #include "draw.h" 96 #include "pinclude/pviewer.h" 97 #include "sys.h" 98 99 static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 100 { 101 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 102 int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2; 103 Scalar *aa; 104 105 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 106 col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 107 col_lens[0] = MAT_COOKIE; 108 col_lens[1] = a->m; 109 col_lens[2] = a->n; 110 col_lens[3] = a->nz*bs2; 111 112 /* store lengths of each row and write (including header) to file */ 113 count = 0; 114 for ( i=0; i<a->mbs; i++ ) { 115 for ( j=0; j<bs; j++ ) { 116 col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 117 } 118 } 119 ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 120 PetscFree(col_lens); 121 122 /* store column indices (zero start index) */ 123 jj = (int *) PetscMalloc( a->nz*bs2*sizeof(int) ); CHKPTRQ(jj); 124 count = 0; 125 for ( i=0; i<a->mbs; i++ ) { 126 for ( j=0; j<bs; j++ ) { 127 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 128 for ( l=0; l<bs; l++ ) { 129 jj[count++] = bs*a->j[k] + l; 130 } 131 } 132 } 133 } 134 ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,BINARY_INT,0); CHKERRQ(ierr); 135 PetscFree(jj); 136 137 /* store nonzero values */ 138 aa = (Scalar *) PetscMalloc(a->nz*bs2*sizeof(Scalar)); CHKPTRQ(aa); 139 count = 0; 140 for ( i=0; i<a->mbs; i++ ) { 141 for ( j=0; j<bs; j++ ) { 142 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 143 for ( l=0; l<bs; l++ ) { 144 aa[count++] = a->a[bs2*k + l*bs + j]; 145 } 146 } 147 } 148 } 149 ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 150 PetscFree(aa); 151 return 0; 152 } 153 154 static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 155 { 156 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 157 int ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2; 158 FILE *fd; 159 char *outputname; 160 161 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 162 ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 163 ierr = ViewerGetFormat(viewer,&format); 164 if (format == ASCII_FORMAT_INFO || format == ASCII_FORMAT_INFO_DETAILED) { 165 fprintf(fd," block size is %d\n",bs); 166 } 167 else if (format == ASCII_FORMAT_MATLAB) { 168 SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported"); 169 } 170 else if (format == ASCII_FORMAT_COMMON) { 171 for ( i=0; i<a->mbs; i++ ) { 172 for ( j=0; j<bs; j++ ) { 173 fprintf(fd,"row %d:",i*bs+j); 174 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 175 for ( l=0; l<bs; l++ ) { 176 #if defined(PETSC_COMPLEX) 177 if (imag(a->a[bs2*k + l*bs + j]) != 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0) 178 fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 179 real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 180 else if (real(a->a[bs2*k + l*bs + j]) != 0.0) 181 fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 182 #else 183 if (a->a[bs2*k + l*bs + j] != 0.0) 184 fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 185 #endif 186 } 187 } 188 fprintf(fd,"\n"); 189 } 190 } 191 } 192 else { 193 for ( i=0; i<a->mbs; i++ ) { 194 for ( j=0; j<bs; j++ ) { 195 fprintf(fd,"row %d:",i*bs+j); 196 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 197 for ( l=0; l<bs; l++ ) { 198 #if defined(PETSC_COMPLEX) 199 if (imag(a->a[bs2*k + l*bs + j]) != 0.0) { 200 fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 201 real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j])); 202 } 203 else { 204 fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j])); 205 } 206 #else 207 fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]); 208 #endif 209 } 210 } 211 fprintf(fd,"\n"); 212 } 213 } 214 } 215 fflush(fd); 216 return 0; 217 } 218 219 static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer) 220 { 221 Mat A = (Mat) obj; 222 ViewerType vtype; 223 int ierr; 224 225 if (!viewer) { 226 viewer = STDOUT_VIEWER_SELF; 227 } 228 229 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 230 if (vtype == MATLAB_VIEWER) { 231 SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported"); 232 } 233 else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 234 return MatView_SeqBAIJ_ASCII(A,viewer); 235 } 236 else if (vtype == BINARY_FILE_VIEWER) { 237 return MatView_SeqBAIJ_Binary(A,viewer); 238 } 239 else if (vtype == DRAW_VIEWER) { 240 SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported"); 241 } 242 return 0; 243 } 244 245 #define CHUNKSIZE 10 246 247 /* This version has row oriented v */ 248 static int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 249 { 250 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 251 int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 252 int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 253 int *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol; 254 int ridx,cidx,bs2=a->bs2; 255 Scalar *ap,value,*aa=a->a,*bap; 256 257 for ( k=0; k<m; k++ ) { /* loop over added rows */ 258 row = im[k]; brow = row/bs; 259 if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row"); 260 if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large"); 261 rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 262 rmax = imax[brow]; nrow = ailen[brow]; 263 low = 0; 264 for ( l=0; l<n; l++ ) { /* loop over added columns */ 265 if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column"); 266 if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large"); 267 col = in[l]; bcol = col/bs; 268 ridx = row % bs; cidx = col % bs; 269 if (roworiented) { 270 value = *v++; 271 } 272 else { 273 value = v[k + l*m]; 274 } 275 if (!sorted) low = 0; high = nrow; 276 while (high-low > 5) { 277 t = (low+high)/2; 278 if (rp[t] > bcol) high = t; 279 else low = t; 280 } 281 for ( i=low; i<high; i++ ) { 282 if (rp[i] > bcol) break; 283 if (rp[i] == bcol) { 284 bap = ap + bs2*i + bs*cidx + ridx; 285 if (is == ADD_VALUES) *bap += value; 286 else *bap = value; 287 goto noinsert; 288 } 289 } 290 if (nonew) goto noinsert; 291 if (nrow >= rmax) { 292 /* there is no extra room in row, therefore enlarge */ 293 int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 294 Scalar *new_a; 295 296 /* malloc new storage space */ 297 len = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 298 new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 299 new_j = (int *) (new_a + bs2*new_nz); 300 new_i = new_j + new_nz; 301 302 /* copy over old data into new slots */ 303 for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 304 for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 305 PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int)); 306 len = (new_nz - CHUNKSIZE - ai[brow] - nrow); 307 PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow, 308 len*sizeof(int)); 309 PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar)); 310 PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar)); 311 PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE), 312 aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar)); 313 /* free up old matrix storage */ 314 PetscFree(a->a); 315 if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 316 aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 317 a->singlemalloc = 1; 318 319 rp = aj + ai[brow]; ap = aa + bs2*ai[brow]; 320 rmax = imax[brow] = imax[brow] + CHUNKSIZE; 321 PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar))); 322 a->maxnz += bs2*CHUNKSIZE; 323 a->reallocs++; 324 a->nz++; 325 } 326 N = nrow++ - 1; 327 /* shift up all the later entries in this row */ 328 for ( ii=N; ii>=i; ii-- ) { 329 rp[ii+1] = rp[ii]; 330 PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar)); 331 } 332 if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar)); 333 rp[i] = bcol; 334 ap[bs2*i + bs*cidx + ridx] = value; 335 noinsert:; 336 low = i; 337 } 338 ailen[brow] = nrow; 339 } 340 return 0; 341 } 342 343 static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 344 { 345 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 346 *m = a->m; *n = a->n; 347 return 0; 348 } 349 350 static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 351 { 352 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 353 *m = 0; *n = a->m; 354 return 0; 355 } 356 357 int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 358 { 359 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 360 int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2; 361 Scalar *aa,*v_i,*aa_i; 362 363 bs = a->bs; 364 ai = a->i; 365 aj = a->j; 366 aa = a->a; 367 bs2 = a->bs2; 368 369 if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range"); 370 371 bn = row/bs; /* Block number */ 372 bp = row % bs; /* Block Position */ 373 M = ai[bn+1] - ai[bn]; 374 *nz = bs*M; 375 376 if (v) { 377 *v = 0; 378 if (*nz) { 379 *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 380 for ( i=0; i<M; i++ ) { /* for each block in the block row */ 381 v_i = *v + i*bs; 382 aa_i = aa + bs2*(ai[bn] + i); 383 for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];} 384 } 385 } 386 } 387 388 if (idx) { 389 *idx = 0; 390 if (*nz) { 391 *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 392 for ( i=0; i<M; i++ ) { /* for each block in the block row */ 393 idx_i = *idx + i*bs; 394 itmp = bs*aj[ai[bn] + i]; 395 for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 396 } 397 } 398 } 399 return 0; 400 } 401 402 int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 403 { 404 if (idx) {if (*idx) PetscFree(*idx);} 405 if (v) {if (*v) PetscFree(*v);} 406 return 0; 407 } 408 409 static int MatTranspose_SeqBAIJ(Mat A,Mat *B) 410 { 411 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 412 Mat C; 413 int i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col; 414 int *rows,*cols,bs2=a->bs2; 415 Scalar *array=a->a; 416 417 if (B==PETSC_NULL && mbs!=nbs) 418 SETERRQ(1,"MatTranspose_SeqBAIJ:Square matrix only for in-place"); 419 col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col); 420 PetscMemzero(col,(1+nbs)*sizeof(int)); 421 422 for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1; 423 ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr); 424 PetscFree(col); 425 rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows); 426 cols = rows + bs; 427 for ( i=0; i<mbs; i++ ) { 428 cols[0] = i*bs; 429 for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1; 430 len = ai[i+1] - ai[i]; 431 for ( j=0; j<len; j++ ) { 432 rows[0] = (*aj++)*bs; 433 for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1; 434 ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr); 435 array += bs2; 436 } 437 } 438 PetscFree(rows); 439 440 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 441 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 442 443 if (B != PETSC_NULL) { 444 *B = C; 445 } else { 446 /* This isn't really an in-place transpose */ 447 PetscFree(a->a); 448 if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 449 if (a->diag) PetscFree(a->diag); 450 if (a->ilen) PetscFree(a->ilen); 451 if (a->imax) PetscFree(a->imax); 452 if (a->solve_work) PetscFree(a->solve_work); 453 PetscFree(a); 454 PetscMemcpy(A,C,sizeof(struct _Mat)); 455 PetscHeaderDestroy(C); 456 } 457 return 0; 458 } 459 460 461 static int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode) 462 { 463 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 464 int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 465 int m = a->m,*ip, N, *ailen = a->ilen; 466 int mbs = a->mbs, bs2 = a->bs2; 467 Scalar *aa = a->a, *ap; 468 469 if (mode == FLUSH_ASSEMBLY) return 0; 470 471 for ( i=1; i<mbs; i++ ) { 472 /* move each row back by the amount of empty slots (fshift) before it*/ 473 fshift += imax[i-1] - ailen[i-1]; 474 if (fshift) { 475 ip = aj + ai[i]; ap = aa + bs2*ai[i]; 476 N = ailen[i]; 477 for ( j=0; j<N; j++ ) { 478 ip[j-fshift] = ip[j]; 479 PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar)); 480 } 481 } 482 ai[i] = ai[i-1] + ailen[i-1]; 483 } 484 if (mbs) { 485 fshift += imax[mbs-1] - ailen[mbs-1]; 486 ai[mbs] = ai[mbs-1] + ailen[mbs-1]; 487 } 488 /* reset ilen and imax for each row */ 489 for ( i=0; i<mbs; i++ ) { 490 ailen[i] = imax[i] = ai[i+1] - ai[i]; 491 } 492 a->nz = ai[mbs]; 493 494 /* diagonals may have moved, so kill the diagonal pointers */ 495 if (fshift && a->diag) { 496 PetscFree(a->diag); 497 PLogObjectMemory(A,-(m+1)*sizeof(int)); 498 a->diag = 0; 499 } 500 PLogInfo(A,"MatAssemblyEnd_SeqBAIJ: Unneed storage space (blocks) %d used %d, rows %d, block size %d\n", fshift*bs2,a->nz*bs2,m,a->bs); 501 PLogInfo(A,"MatAssemblyEnd_SeqBAIJ: Number of mallocs during MatSetValues %d\n", 502 a->reallocs); 503 return 0; 504 } 505 506 static int MatZeroEntries_SeqBAIJ(Mat A) 507 { 508 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 509 PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(Scalar)); 510 return 0; 511 } 512 513 int MatDestroy_SeqBAIJ(PetscObject obj) 514 { 515 Mat A = (Mat) obj; 516 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 517 518 #if defined(PETSC_LOG) 519 PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 520 #endif 521 PetscFree(a->a); 522 if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 523 if (a->diag) PetscFree(a->diag); 524 if (a->ilen) PetscFree(a->ilen); 525 if (a->imax) PetscFree(a->imax); 526 if (a->solve_work) PetscFree(a->solve_work); 527 if (a->mult_work) PetscFree(a->mult_work); 528 PetscFree(a); 529 PLogObjectDestroy(A); 530 PetscHeaderDestroy(A); 531 return 0; 532 } 533 534 static int MatSetOption_SeqBAIJ(Mat A,MatOption op) 535 { 536 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 537 if (op == ROW_ORIENTED) a->roworiented = 1; 538 else if (op == COLUMN_ORIENTED) a->roworiented = 0; 539 else if (op == COLUMNS_SORTED) a->sorted = 1; 540 else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 541 else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 542 else if (op == ROWS_SORTED || 543 op == SYMMETRIC_MATRIX || 544 op == STRUCTURALLY_SYMMETRIC_MATRIX || 545 op == YES_NEW_DIAGONALS) 546 PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n"); 547 else if (op == NO_NEW_DIAGONALS) 548 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");} 549 else 550 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");} 551 return 0; 552 } 553 554 555 /* -------------------------------------------------------*/ 556 /* Should check that shapes of vectors and matrices match */ 557 /* -------------------------------------------------------*/ 558 #include "pinclude/plapack.h" 559 560 static int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 561 { 562 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 563 Scalar *xg,*zg; 564 register Scalar *x,*z,*v,sum; 565 int mbs=a->mbs,i,*idx,*ii,n,ierr; 566 567 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 568 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 569 570 idx = a->j; 571 v = a->a; 572 ii = a->i; 573 574 for ( i=0; i<mbs; i++ ) { 575 n = ii[1] - ii[0]; ii++; 576 sum = 0.0; 577 while (n--) sum += *v++ * x[*idx++]; 578 z[i] = sum; 579 } 580 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 581 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 582 PLogFlops(2*a->nz - a->m); 583 return 0; 584 } 585 586 static int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 587 { 588 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 589 Scalar *xg,*zg; 590 register Scalar *x,*z,*v,*xb,sum1,sum2; 591 register Scalar x1,x2; 592 int mbs=a->mbs,i,*idx,*ii; 593 int j,n,ierr; 594 595 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 596 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 597 598 idx = a->j; 599 v = a->a; 600 ii = a->i; 601 602 for ( i=0; i<mbs; i++ ) { 603 n = ii[1] - ii[0]; ii++; 604 sum1 = 0.0; sum2 = 0.0; 605 for ( j=0; j<n; j++ ) { 606 xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 607 sum1 += v[0]*x1 + v[2]*x2; 608 sum2 += v[1]*x1 + v[3]*x2; 609 v += 4; 610 } 611 z[0] = sum1; z[1] = sum2; 612 z += 2; 613 } 614 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 615 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 616 PLogFlops(4*a->nz - a->m); 617 return 0; 618 } 619 620 static int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 621 { 622 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 623 Scalar *xg,*zg; 624 register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,x1,x2,x3; 625 int mbs=a->mbs,i,*idx,*ii,j,n,ierr; 626 627 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 628 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 629 630 idx = a->j; 631 v = a->a; 632 ii = a->i; 633 634 for ( i=0; i<mbs; i++ ) { 635 n = ii[1] - ii[0]; ii++; 636 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 637 for ( j=0; j<n; j++ ) { 638 xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 639 sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 640 sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 641 sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 642 v += 9; 643 } 644 z[0] = sum1; z[1] = sum2; z[2] = sum3; 645 z += 3; 646 } 647 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 648 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 649 PLogFlops(18*a->nz - a->m); 650 return 0; 651 } 652 653 static int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) 654 { 655 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 656 Scalar *xg,*zg; 657 register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,sum4; 658 register Scalar x1,x2,x3,x4; 659 int mbs=a->mbs,i,*idx,*ii; 660 int j,n,ierr; 661 662 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 663 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 664 665 idx = a->j; 666 v = a->a; 667 ii = a->i; 668 669 for ( i=0; i<mbs; i++ ) { 670 n = ii[1] - ii[0]; ii++; 671 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 672 for ( j=0; j<n; j++ ) { 673 xb = x + 4*(*idx++); 674 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 675 sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 676 sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 677 sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 678 sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 679 v += 16; 680 } 681 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 682 z += 4; 683 } 684 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 685 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 686 PLogFlops(32*a->nz - a->m); 687 return 0; 688 } 689 690 static int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 691 { 692 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 693 Scalar *xg,*zg; 694 register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,sum4,sum5; 695 register Scalar x1,x2,x3,x4,x5; 696 int mbs=a->mbs,i,*idx,*ii,j,n,ierr; 697 698 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 699 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 700 701 idx = a->j; 702 v = a->a; 703 ii = a->i; 704 705 for ( i=0; i<mbs; i++ ) { 706 n = ii[1] - ii[0]; ii++; 707 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 708 for ( j=0; j<n; j++ ) { 709 xb = x + 5*(*idx++); 710 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 711 sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 712 sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 713 sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 714 sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 715 sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 716 v += 25; 717 } 718 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 719 z += 5; 720 } 721 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 722 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 723 PLogFlops(50*a->nz - a->m); 724 return 0; 725 } 726 727 static int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 728 { 729 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 730 Scalar *xg,*zg; 731 register Scalar *x,*z,*v,*xb; 732 int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr; 733 int _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt, _DZero = 0.0; 734 735 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 736 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 737 738 idx = a->j; 739 v = a->a; 740 ii = a->i; 741 742 743 if (!a->mult_work) { 744 k = PetscMax(a->m,a->n); 745 a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));CHKPTRQ(a->mult_work); 746 } 747 work = a->mult_work; 748 for ( i=0; i<mbs; i++ ) { 749 n = ii[1] - ii[0]; ii++; 750 ncols = n*bs; 751 workt = work; 752 for ( j=0; j<n; j++ ) { 753 xb = x + bs*(*idx++); 754 for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 755 workt += bs; 756 } 757 LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); 758 v += n*bs2; 759 z += bs; 760 } 761 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 762 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 763 PLogFlops(2*a->nz*bs2 - a->m); 764 return 0; 765 } 766 767 static int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz) 768 { 769 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 770 Scalar *xg,*zg,*zb; 771 register Scalar *x,*z,*v,*xb,x1,x2,x3,x4,x5; 772 int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n; 773 int bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 774 775 776 ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg; 777 ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg; 778 PetscMemzero(z,N*sizeof(Scalar)); 779 780 idx = a->j; 781 v = a->a; 782 ii = a->i; 783 784 switch (bs) { 785 case 1: 786 for ( i=0; i<mbs; i++ ) { 787 n = ii[1] - ii[0]; ii++; 788 xb = x + i; x1 = xb[0]; 789 ib = idx + ai[i]; 790 for ( j=0; j<n; j++ ) { 791 rval = ib[j]; 792 z[rval] += *v++ * x1; 793 } 794 } 795 break; 796 case 2: 797 for ( i=0; i<mbs; i++ ) { 798 n = ii[1] - ii[0]; ii++; 799 xb = x + 2*i; x1 = xb[0]; x2 = xb[1]; 800 ib = idx + ai[i]; 801 for ( j=0; j<n; j++ ) { 802 rval = ib[j]*2; 803 z[rval++] += v[0]*x1 + v[1]*x2; 804 z[rval++] += v[2]*x1 + v[3]*x2; 805 v += 4; 806 } 807 } 808 break; 809 case 3: 810 for ( i=0; i<mbs; i++ ) { 811 n = ii[1] - ii[0]; ii++; 812 xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 813 ib = idx + ai[i]; 814 for ( j=0; j<n; j++ ) { 815 rval = ib[j]*3; 816 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 817 z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 818 z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 819 v += 9; 820 } 821 } 822 break; 823 case 4: 824 for ( i=0; i<mbs; i++ ) { 825 n = ii[1] - ii[0]; ii++; 826 xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 827 ib = idx + ai[i]; 828 for ( j=0; j<n; j++ ) { 829 rval = ib[j]*4; 830 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 831 z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 832 z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 833 z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 834 v += 16; 835 } 836 } 837 break; 838 case 5: 839 for ( i=0; i<mbs; i++ ) { 840 n = ii[1] - ii[0]; ii++; 841 xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 842 x4 = xb[3]; x5 = xb[4]; 843 ib = idx + ai[i]; 844 for ( j=0; j<n; j++ ) { 845 rval = ib[j]*5; 846 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 847 z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 848 z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 849 z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 850 z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 851 v += 25; 852 } 853 } 854 break; 855 /* block sizes larger then 5 by 5 are handled by BLAS */ 856 default: { 857 int _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt; 858 if (!a->mult_work) { 859 k = PetscMax(a->m,a->n); 860 a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar)); 861 CHKPTRQ(a->mult_work); 862 } 863 work = a->mult_work; 864 for ( i=0; i<mbs; i++ ) { 865 n = ii[1] - ii[0]; ii++; 866 ncols = n*bs; 867 PetscMemzero(work,ncols*sizeof(Scalar)); 868 LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); 869 v += n*bs2; 870 x += bs; 871 workt = work; 872 for ( j=0; j<n; j++ ) { 873 zb = z + bs*(*idx++); 874 for ( k=0; k<bs; k++ ) zb[k] += workt[k] ; 875 workt += bs; 876 } 877 } 878 } 879 } 880 ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr); 881 ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr); 882 return 0; 883 } 884 885 static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem) 886 { 887 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 888 if (nz) *nz = a->bs2*a->nz; 889 if (nza) *nza = a->maxnz; 890 if (mem) *mem = (int)A->mem; 891 return 0; 892 } 893 894 static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg) 895 { 896 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 897 898 if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type"); 899 900 /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 901 if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| 902 (a->nz != b->nz)) { 903 *flg = PETSC_FALSE; return 0; 904 } 905 906 /* if the a->i are the same */ 907 if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) { 908 *flg = PETSC_FALSE; return 0; 909 } 910 911 /* if a->j are the same */ 912 if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) { 913 *flg = PETSC_FALSE; return 0; 914 } 915 916 /* if a->a are the same */ 917 if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) { 918 *flg = PETSC_FALSE; return 0; 919 } 920 *flg = PETSC_TRUE; 921 return 0; 922 923 } 924 925 static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 926 { 927 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 928 int i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 929 Scalar *x, zero = 0.0,*aa,*aa_j; 930 931 bs = a->bs; 932 aa = a->a; 933 ai = a->i; 934 aj = a->j; 935 ambs = a->mbs; 936 bs2 = a->bs2; 937 938 VecSet(&zero,v); 939 VecGetArray(v,&x); VecGetLocalSize(v,&n); 940 if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector"); 941 for ( i=0; i<ambs; i++ ) { 942 for ( j=ai[i]; j<ai[i+1]; j++ ) { 943 if (aj[j] == i) { 944 row = i*bs; 945 aa_j = aa+j*bs2; 946 for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 947 break; 948 } 949 } 950 } 951 return 0; 952 } 953 954 static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 955 { 956 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 957 Scalar *l,*r,x,*v,*aa,*li,*ri; 958 int i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; 959 960 ai = a->i; 961 aj = a->j; 962 aa = a->a; 963 m = a->m; 964 n = a->n; 965 bs = a->bs; 966 mbs = a->mbs; 967 bs2 = a->bs2; 968 if (ll) { 969 VecGetArray(ll,&l); VecGetSize(ll,&lm); 970 if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length"); 971 for ( i=0; i<mbs; i++ ) { /* for each block row */ 972 M = ai[i+1] - ai[i]; 973 li = l + i*bs; 974 v = aa + bs2*ai[i]; 975 for ( j=0; j<M; j++ ) { /* for each block */ 976 for ( k=0; k<bs2; k++ ) { 977 (*v++) *= li[k%bs]; 978 } 979 } 980 } 981 } 982 983 if (rr) { 984 VecGetArray(rr,&r); VecGetSize(rr,&rn); 985 if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length"); 986 for ( i=0; i<mbs; i++ ) { /* for each block row */ 987 M = ai[i+1] - ai[i]; 988 v = aa + bs2*ai[i]; 989 for ( j=0; j<M; j++ ) { /* for each block */ 990 ri = r + bs*aj[ai[i]+j]; 991 for ( k=0; k<bs; k++ ) { 992 x = ri[k]; 993 for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x; 994 } 995 } 996 } 997 } 998 return 0; 999 } 1000 1001 1002 extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 1003 extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 1004 extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int); 1005 extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,MatGetSubMatrixCall,Mat*); 1006 extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**); 1007 1008 extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec); 1009 extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec); 1010 extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec); 1011 extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec); 1012 extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec); 1013 extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec); 1014 1015 extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 1016 extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 1017 extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 1018 extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 1019 extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 1020 extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 1021 1022 static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 1023 { 1024 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1025 Scalar *v = a->a; 1026 double sum = 0.0; 1027 int i,nz=a->nz,bs2=a->bs2; 1028 1029 if (type == NORM_FROBENIUS) { 1030 for (i=0; i< bs2*nz; i++ ) { 1031 #if defined(PETSC_COMPLEX) 1032 sum += real(conj(*v)*(*v)); v++; 1033 #else 1034 sum += (*v)*(*v); v++; 1035 #endif 1036 } 1037 *norm = sqrt(sum); 1038 } 1039 else { 1040 SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet"); 1041 } 1042 return 0; 1043 } 1044 1045 /* 1046 note: This can only work for identity for row and col. It would 1047 be good to check this and otherwise generate an error. 1048 */ 1049 static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 1050 { 1051 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 1052 Mat outA; 1053 int ierr; 1054 1055 if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported"); 1056 1057 outA = inA; 1058 inA->factor = FACTOR_LU; 1059 a->row = row; 1060 a->col = col; 1061 1062 a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 1063 1064 if (!a->diag) { 1065 ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 1066 } 1067 ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 1068 return 0; 1069 } 1070 1071 static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 1072 { 1073 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 1074 int one = 1, totalnz = a->bs2*a->nz; 1075 BLscal_( &totalnz, alpha, a->a, &one ); 1076 PLogFlops(totalnz); 1077 return 0; 1078 } 1079 1080 static int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 1081 { 1082 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1083 int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 1084 int *ai = a->i, *ailen = a->ilen; 1085 int brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2; 1086 Scalar *ap, *aa = a->a, zero = 0.0; 1087 1088 for ( k=0; k<m; k++ ) { /* loop over rows */ 1089 row = im[k]; brow = row/bs; 1090 if (row < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative row"); 1091 if (row >= a->m) SETERRQ(1,"MatGetValues_SeqBAIJ:Row too large"); 1092 rp = aj + ai[brow] ; ap = aa + bs2*ai[brow] ; 1093 nrow = ailen[brow]; 1094 for ( l=0; l<n; l++ ) { /* loop over columns */ 1095 if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative column"); 1096 if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqBAIJ:Column too large"); 1097 col = in[l] ; 1098 bcol = col/bs; 1099 cidx = col%bs; 1100 ridx = row%bs; 1101 high = nrow; 1102 low = 0; /* assume unsorted */ 1103 while (high-low > 5) { 1104 t = (low+high)/2; 1105 if (rp[t] > bcol) high = t; 1106 else low = t; 1107 } 1108 for ( i=low; i<high; i++ ) { 1109 if (rp[i] > bcol) break; 1110 if (rp[i] == bcol) { 1111 *v++ = ap[bs2*i+bs*cidx+ridx]; 1112 goto finished; 1113 } 1114 } 1115 *v++ = zero; 1116 finished:; 1117 } 1118 } 1119 return 0; 1120 } 1121 1122 /* -------------------------------------------------------------------*/ 1123 static struct _MatOps MatOps = {MatSetValues_SeqBAIJ, 1124 MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ, 1125 MatMult_SeqBAIJ_N,0, 1126 MatMultTrans_SeqBAIJ,0, 1127 MatSolve_SeqBAIJ_N,0, 1128 0,0, 1129 MatLUFactor_SeqBAIJ,0, 1130 0, 1131 MatTranspose_SeqBAIJ, 1132 MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ, 1133 MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ, 1134 0,MatAssemblyEnd_SeqBAIJ, 1135 0, 1136 MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0, 1137 MatGetReordering_SeqBAIJ, 1138 MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 1139 MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 1140 MatILUFactorSymbolic_SeqBAIJ,0, 1141 0,0,/* MatConvert_SeqBAIJ */ 0, 1142 MatGetSubMatrix_SeqBAIJ,0, 1143 MatConvertSameType_SeqBAIJ,0,0, 1144 MatILUFactor_SeqBAIJ,0,0, 1145 MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ, 1146 MatGetValues_SeqBAIJ,0, 1147 0,MatScale_SeqBAIJ, 1148 0}; 1149 1150 /*@C 1151 MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block 1152 compressed row) format. For good matrix assembly performance the 1153 user should preallocate the matrix storage by setting the parameter nz 1154 (or nzz). By setting these parameters accurately, performance can be 1155 increased by more than a factor of 50. 1156 1157 Input Parameters: 1158 . comm - MPI communicator, set to MPI_COMM_SELF 1159 . bs - size of block 1160 . m - number of rows 1161 . n - number of columns 1162 . nz - number of block nonzeros per block row (same for all rows) 1163 . nzz - number of block nonzeros per block row or PETSC_NULL 1164 (possibly different for each row) 1165 1166 Output Parameter: 1167 . A - the matrix 1168 1169 Notes: 1170 The block AIJ format is fully compatible with standard Fortran 77 1171 storage. That is, the stored row and column indices can begin at 1172 either one (as in Fortran) or zero. See the users' manual for details. 1173 1174 Specify the preallocated storage with either nz or nnz (not both). 1175 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1176 allocation. For additional details, see the users manual chapter on 1177 matrices and the file $(PETSC_DIR)/Performance. 1178 1179 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues() 1180 @*/ 1181 int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 1182 { 1183 Mat B; 1184 Mat_SeqBAIJ *b; 1185 int i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs; 1186 1187 if (mbs*bs!=m || nbs*bs!=n) 1188 SETERRQ(1,"MatCreateSeqBAIJ:Number rows, cols must be divisible by blocksize"); 1189 1190 *A = 0; 1191 PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm); 1192 PLogObjectCreate(B); 1193 B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 1194 PetscMemzero(b,sizeof(Mat_SeqBAIJ)); 1195 PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1196 ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr); 1197 if (!flg) { 1198 switch (bs) { 1199 case 1: 1200 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 1201 B->ops.solve = MatSolve_SeqBAIJ_1; 1202 B->ops.mult = MatMult_SeqBAIJ_1; 1203 break; 1204 case 2: 1205 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 1206 B->ops.solve = MatSolve_SeqBAIJ_2; 1207 B->ops.mult = MatMult_SeqBAIJ_2; 1208 break; 1209 case 3: 1210 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 1211 B->ops.solve = MatSolve_SeqBAIJ_3; 1212 B->ops.mult = MatMult_SeqBAIJ_3; 1213 break; 1214 case 4: 1215 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 1216 B->ops.solve = MatSolve_SeqBAIJ_4; 1217 B->ops.mult = MatMult_SeqBAIJ_4; 1218 break; 1219 case 5: 1220 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 1221 B->ops.solve = MatSolve_SeqBAIJ_5; 1222 B->ops.mult = MatMult_SeqBAIJ_5; 1223 break; 1224 } 1225 } 1226 B->destroy = MatDestroy_SeqBAIJ; 1227 B->view = MatView_SeqBAIJ; 1228 B->factor = 0; 1229 B->lupivotthreshold = 1.0; 1230 b->row = 0; 1231 b->col = 0; 1232 b->reallocs = 0; 1233 1234 b->m = m; B->m = m; B->M = m; 1235 b->n = n; B->n = n; B->N = n; 1236 b->mbs = mbs; 1237 b->nbs = nbs; 1238 b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 1239 if (nnz == PETSC_NULL) { 1240 if (nz == PETSC_DEFAULT) nz = 5; 1241 else if (nz <= 0) nz = 1; 1242 for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 1243 nz = nz*mbs; 1244 } 1245 else { 1246 nz = 0; 1247 for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 1248 } 1249 1250 /* allocate the matrix space */ 1251 len = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int); 1252 b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1253 PetscMemzero(b->a,nz*bs2*sizeof(Scalar)); 1254 b->j = (int *) (b->a + nz*bs2); 1255 PetscMemzero(b->j,nz*sizeof(int)); 1256 b->i = b->j + nz; 1257 b->singlemalloc = 1; 1258 1259 b->i[0] = 0; 1260 for (i=1; i<mbs+1; i++) { 1261 b->i[i] = b->i[i-1] + b->imax[i-1]; 1262 } 1263 1264 /* b->ilen will count nonzeros in each block row so far. */ 1265 b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 1266 PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 1267 for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 1268 1269 b->bs = bs; 1270 b->bs2 = bs2; 1271 b->mbs = mbs; 1272 b->nz = 0; 1273 b->maxnz = nz*bs2; 1274 b->sorted = 0; 1275 b->roworiented = 1; 1276 b->nonew = 0; 1277 b->diag = 0; 1278 b->solve_work = 0; 1279 b->mult_work = 0; 1280 b->spptr = 0; 1281 *A = B; 1282 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1283 if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 1284 return 0; 1285 } 1286 1287 int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 1288 { 1289 Mat C; 1290 Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 1291 int i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2; 1292 1293 if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix"); 1294 1295 *B = 0; 1296 PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm); 1297 PLogObjectCreate(C); 1298 C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 1299 PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1300 C->destroy = MatDestroy_SeqBAIJ; 1301 C->view = MatView_SeqBAIJ; 1302 C->factor = A->factor; 1303 c->row = 0; 1304 c->col = 0; 1305 C->assembled = PETSC_TRUE; 1306 1307 c->m = C->m = a->m; 1308 c->n = C->n = a->n; 1309 C->M = a->m; 1310 C->N = a->n; 1311 1312 c->bs = a->bs; 1313 c->bs2 = a->bs2; 1314 c->mbs = a->mbs; 1315 c->nbs = a->nbs; 1316 1317 c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 1318 c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 1319 for ( i=0; i<mbs; i++ ) { 1320 c->imax[i] = a->imax[i]; 1321 c->ilen[i] = a->ilen[i]; 1322 } 1323 1324 /* allocate the matrix space */ 1325 c->singlemalloc = 1; 1326 len = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int)); 1327 c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1328 c->j = (int *) (c->a + nz*bs2); 1329 c->i = c->j + nz; 1330 PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 1331 if (mbs > 0) { 1332 PetscMemcpy(c->j,a->j,nz*sizeof(int)); 1333 if (cpvalues == COPY_VALUES) { 1334 PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar)); 1335 } 1336 } 1337 1338 PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 1339 c->sorted = a->sorted; 1340 c->roworiented = a->roworiented; 1341 c->nonew = a->nonew; 1342 1343 if (a->diag) { 1344 c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 1345 PLogObjectMemory(C,(mbs+1)*sizeof(int)); 1346 for ( i=0; i<mbs; i++ ) { 1347 c->diag[i] = a->diag[i]; 1348 } 1349 } 1350 else c->diag = 0; 1351 c->nz = a->nz; 1352 c->maxnz = a->maxnz; 1353 c->solve_work = 0; 1354 c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1355 c->mult_work = 0; 1356 *B = C; 1357 return 0; 1358 } 1359 1360 int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 1361 { 1362 Mat_SeqBAIJ *a; 1363 Mat B; 1364 int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 1365 int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 1366 int kmax,jcount,block,idx,point,nzcountb,extra_rows; 1367 int *masked, nmask,tmp,bs2,ishift; 1368 Scalar *aa; 1369 MPI_Comm comm = ((PetscObject) viewer)->comm; 1370 1371 ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr); 1372 bs2 = bs*bs; 1373 1374 MPI_Comm_size(comm,&size); 1375 if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor"); 1376 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 1377 ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 1378 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object"); 1379 M = header[1]; N = header[2]; nz = header[3]; 1380 1381 if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices"); 1382 1383 /* 1384 This code adds extra rows to make sure the number of rows is 1385 divisible by the blocksize 1386 */ 1387 mbs = M/bs; 1388 extra_rows = bs - M + bs*(mbs); 1389 if (extra_rows == bs) extra_rows = 0; 1390 else mbs++; 1391 if (extra_rows) { 1392 PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize"); 1393 } 1394 1395 /* read in row lengths */ 1396 rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 1397 ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 1398 for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 1399 1400 /* read in column indices */ 1401 jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 1402 ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr); 1403 for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1404 1405 /* loop over row lengths determining block row lengths */ 1406 browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1407 PetscMemzero(browlengths,mbs*sizeof(int)); 1408 mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 1409 PetscMemzero(mask,mbs*sizeof(int)); 1410 masked = mask + mbs; 1411 rowcount = 0; nzcount = 0; 1412 for ( i=0; i<mbs; i++ ) { 1413 nmask = 0; 1414 for ( j=0; j<bs; j++ ) { 1415 kmax = rowlengths[rowcount]; 1416 for ( k=0; k<kmax; k++ ) { 1417 tmp = jj[nzcount++]/bs; 1418 if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1419 } 1420 rowcount++; 1421 } 1422 browlengths[i] += nmask; 1423 /* zero out the mask elements we set */ 1424 for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1425 } 1426 1427 /* create our matrix */ 1428 ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A); 1429 CHKERRQ(ierr); 1430 B = *A; 1431 a = (Mat_SeqBAIJ *) B->data; 1432 1433 /* set matrix "i" values */ 1434 a->i[0] = 0; 1435 for ( i=1; i<= mbs; i++ ) { 1436 a->i[i] = a->i[i-1] + browlengths[i-1]; 1437 a->ilen[i-1] = browlengths[i-1]; 1438 } 1439 a->nz = 0; 1440 for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 1441 1442 /* read in nonzero values */ 1443 aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 1444 ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr); 1445 for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1446 1447 /* set "a" and "j" values into matrix */ 1448 nzcount = 0; jcount = 0; 1449 for ( i=0; i<mbs; i++ ) { 1450 nzcountb = nzcount; 1451 nmask = 0; 1452 for ( j=0; j<bs; j++ ) { 1453 kmax = rowlengths[i*bs+j]; 1454 for ( k=0; k<kmax; k++ ) { 1455 tmp = jj[nzcount++]/bs; 1456 if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1457 } 1458 rowcount++; 1459 } 1460 /* sort the masked values */ 1461 PetscSortInt(nmask,masked); 1462 1463 /* set "j" values into matrix */ 1464 maskcount = 1; 1465 for ( j=0; j<nmask; j++ ) { 1466 a->j[jcount++] = masked[j]; 1467 mask[masked[j]] = maskcount++; 1468 } 1469 /* set "a" values into matrix */ 1470 ishift = bs2*a->i[i]; 1471 for ( j=0; j<bs; j++ ) { 1472 kmax = rowlengths[i*bs+j]; 1473 for ( k=0; k<kmax; k++ ) { 1474 tmp = jj[nzcountb]/bs ; 1475 block = mask[tmp] - 1; 1476 point = jj[nzcountb] - bs*tmp; 1477 idx = ishift + bs2*block + j + bs*point; 1478 a->a[idx] = aa[nzcountb++]; 1479 } 1480 } 1481 /* zero out the mask elements we set */ 1482 for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1483 } 1484 if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix"); 1485 1486 PetscFree(rowlengths); 1487 PetscFree(browlengths); 1488 PetscFree(aa); 1489 PetscFree(jj); 1490 PetscFree(mask); 1491 1492 B->assembled = PETSC_TRUE; 1493 1494 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 1495 if (flg) { 1496 Viewer tviewer; 1497 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1498 ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr); 1499 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1500 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1501 } 1502 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 1503 if (flg) { 1504 Viewer tviewer; 1505 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1506 ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr); 1507 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1508 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1509 } 1510 ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr); 1511 if (flg) { 1512 Viewer tviewer; 1513 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1514 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1515 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1516 } 1517 ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr); 1518 if (flg) { 1519 Viewer tviewer; 1520 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1521 ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr); 1522 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1523 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1524 } 1525 ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 1526 if (flg) { 1527 Viewer tviewer; 1528 ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr); 1529 ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr); 1530 ierr = ViewerFlush(tviewer); CHKERRQ(ierr); 1531 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1532 } 1533 return 0; 1534 } 1535 1536 1537 1538