1 2 #ifndef lint 3 static char vcid[] = "$Id: baij.c,v 1.24 1996/03/31 19:59:04 curfman Exp balay $"; 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 "vec/vecimpl.h" 12 #include "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 = a->indexshift; 42 oshift = -MatReorderingIndexShift[(int)type]; 43 if (MatReorderingRequiresSymmetric[(int)type]) { 44 ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja); 45 CHKERRQ(ierr); 46 ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 47 PetscFree(ia); PetscFree(ja); 48 } else { 49 if (ishift == oshift) { 50 ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 51 } 52 else if (ishift == -1) { 53 /* temporarily subtract 1 from i and j indices */ 54 int nz = a->i[a->n] - 1; 55 for ( i=0; i<nz; i++ ) a->j[i]--; 56 for ( i=0; i<a->n+1; i++ ) a->i[i]--; 57 ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 58 for ( i=0; i<nz; i++ ) a->j[i]++; 59 for ( i=0; i<a->n+1; i++ ) a->i[i]++; 60 } else { 61 /* temporarily add 1 to i and j indices */ 62 int nz = a->i[a->n] - 1; 63 for ( i=0; i<nz; i++ ) a->j[i]++; 64 for ( i=0; i<a->n+1; i++ ) a->i[i]++; 65 ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr); 66 for ( i=0; i<nz; i++ ) a->j[i]--; 67 for ( i=0; i<a->n+1; i++ ) a->i[i]--; 68 } 69 } 70 return 0; 71 } 72 73 /* 74 Adds diagonal pointers to sparse matrix structure. 75 */ 76 77 int MatMarkDiag_SeqBAIJ(Mat A) 78 { 79 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 80 int i,j, *diag, m = a->mbs; 81 82 diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 83 PLogObjectMemory(A,(m+1)*sizeof(int)); 84 for ( i=0; i<m; i++ ) { 85 for ( j=a->i[i]; j<a->i[i+1]; j++ ) { 86 if (a->j[j] == i) { 87 diag[i] = j; 88 break; 89 } 90 } 91 } 92 a->diag = diag; 93 return 0; 94 } 95 96 #include "draw.h" 97 #include "pinclude/pviewer.h" 98 #include "sys.h" 99 100 static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer) 101 { 102 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 103 int i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l; 104 Scalar *aa; 105 106 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 107 col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens); 108 col_lens[0] = MAT_COOKIE; 109 col_lens[1] = a->m; 110 col_lens[2] = a->n; 111 col_lens[3] = a->nz*bs*bs; 112 113 /* store lengths of each row and write (including header) to file */ 114 count = 0; 115 for ( i=0; i<a->mbs; i++ ) { 116 for ( j=0; j<bs; j++ ) { 117 col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]); 118 } 119 } 120 ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 121 PetscFree(col_lens); 122 123 /* store column indices (zero start index) */ 124 jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj); 125 count = 0; 126 for ( i=0; i<a->mbs; i++ ) { 127 for ( j=0; j<bs; j++ ) { 128 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 129 for ( l=0; l<bs; l++ ) { 130 jj[count++] = bs*a->j[k] + l; 131 } 132 } 133 } 134 } 135 ierr = PetscBinaryWrite(fd,jj,bs*bs*a->nz,BINARY_INT,0); CHKERRQ(ierr); 136 PetscFree(jj); 137 138 /* store nonzero values */ 139 aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa); 140 count = 0; 141 for ( i=0; i<a->mbs; i++ ) { 142 for ( j=0; j<bs; j++ ) { 143 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 144 for ( l=0; l<bs; l++ ) { 145 aa[count++] = a->a[bs*bs*k + l*bs + j]; 146 } 147 } 148 } 149 } 150 ierr = PetscBinaryWrite(fd,aa,bs*bs*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 151 PetscFree(aa); 152 return 0; 153 } 154 155 static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer) 156 { 157 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 158 int ierr, i,j,format,bs = a->bs,k,l; 159 FILE *fd; 160 char *outputname; 161 162 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 163 ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr); 164 ierr = ViewerGetFormat(viewer,&format); 165 if (format == ASCII_FORMAT_INFO) { 166 /* no need to print additional information */ ; 167 } 168 else if (format == ASCII_FORMAT_MATLAB) { 169 SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported"); 170 } 171 else { 172 for ( i=0; i<a->mbs; i++ ) { 173 for ( j=0; j<bs; j++ ) { 174 fprintf(fd,"row %d:",i*bs+j); 175 for ( k=a->i[i]; k<a->i[i+1]; k++ ) { 176 for ( l=0; l<bs; l++ ) { 177 #if defined(PETSC_COMPLEX) 178 if (imag(a->a[bs*bs*k + l*bs + j]) != 0.0) { 179 fprintf(fd," %d %g + %g i",bs*a->j[k]+l, 180 real(a->a[bs*bs*k + l*bs + j]),imag(a->a[bs*bs*k + l*bs + j])); 181 } 182 else { 183 fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs*bs*k + l*bs + j])); 184 } 185 #else 186 fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs*bs*k + l*bs + j]); 187 #endif 188 } 189 } 190 fprintf(fd,"\n"); 191 } 192 } 193 } 194 fflush(fd); 195 return 0; 196 } 197 198 static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer) 199 { 200 Mat A = (Mat) obj; 201 ViewerType vtype; 202 int ierr; 203 204 if (!viewer) { 205 viewer = STDOUT_VIEWER_SELF; 206 } 207 208 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 209 if (vtype == MATLAB_VIEWER) { 210 SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported"); 211 } 212 else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 213 return MatView_SeqBAIJ_ASCII(A,viewer); 214 } 215 else if (vtype == BINARY_FILE_VIEWER) { 216 return MatView_SeqBAIJ_Binary(A,viewer); 217 } 218 else if (vtype == DRAW_VIEWER) { 219 SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported"); 220 } 221 return 0; 222 } 223 224 #define CHUNKSIZE 10 225 226 /* This version has row oriented v */ 227 static int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 228 { 229 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 230 int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted; 231 int *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented; 232 int *aj=a->j,nonew=a->nonew,shift=a->indexshift,bs=a->bs,brow,bcol; 233 int ridx,cidx; 234 Scalar *ap,value,*aa=a->a,*bap; 235 236 for ( k=0; k<m; k++ ) { /* loop over added rows */ 237 row = im[k]; brow = row/bs; 238 if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row"); 239 if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large"); 240 rp = aj + ai[brow] + shift; ap = aa + bs*bs*ai[brow] + shift; 241 rmax = imax[brow]; nrow = ailen[brow]; 242 low = 0; 243 for ( l=0; l<n; l++ ) { /* loop over added columns */ 244 if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column"); 245 if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large"); 246 col = in[l] - shift; bcol = col/bs; 247 ridx = row % bs; cidx = col % bs; 248 if (roworiented) { 249 value = *v++; 250 } 251 else { 252 value = v[k + l*m]; 253 } 254 if (!sorted) low = 0; high = nrow; 255 while (high-low > 5) { 256 t = (low+high)/2; 257 if (rp[t] > bcol) high = t; 258 else low = t; 259 } 260 for ( i=low; i<high; i++ ) { 261 if (rp[i] > bcol) break; 262 if (rp[i] == bcol) { 263 bap = ap + bs*bs*i + bs*ridx + cidx; 264 if (is == ADD_VALUES) *bap += value; 265 else *bap = value; 266 goto noinsert; 267 } 268 } 269 if (nonew) goto noinsert; 270 if (nrow >= rmax) { 271 /* there is no extra room in row, therefore enlarge */ 272 int new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j; 273 Scalar *new_a; 274 275 /* malloc new storage space */ 276 len = new_nz*(sizeof(int)+bs*bs*sizeof(Scalar))+(a->mbs+1)*sizeof(int); 277 new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 278 new_j = (int *) (new_a + bs*bs*new_nz); 279 new_i = new_j + new_nz; 280 281 /* copy over old data into new slots */ 282 for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];} 283 for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 284 PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 285 len = (new_nz - CHUNKSIZE - ai[brow] - nrow - shift); 286 PetscMemcpy(new_j+ai[brow]+shift+nrow+CHUNKSIZE,aj+ai[brow]+shift+nrow, 287 len*sizeof(int)); 288 PetscMemcpy(new_a,aa,(ai[brow]+nrow+shift)*bs*bs*sizeof(Scalar)); 289 PetscMemzero(new_a+bs*bs*(ai[brow]+shift+nrow),bs*bs*CHUNKSIZE); 290 PetscMemcpy(new_a+bs*bs*(ai[brow]+shift+nrow+CHUNKSIZE), 291 aa+bs*bs*(ai[brow]+shift+nrow),bs*bs*len*sizeof(Scalar)); 292 /* free up old matrix storage */ 293 PetscFree(a->a); 294 if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 295 aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 296 a->singlemalloc = 1; 297 298 rp = aj + ai[brow] + shift; ap = aa + ai[brow] + shift; 299 rmax = imax[brow] = imax[brow] + CHUNKSIZE; 300 PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs*bs*sizeof(Scalar))); 301 a->maxnz += bs*bs*CHUNKSIZE; 302 a->reallocs++; 303 } 304 N = nrow++ - 1; a->nz++; 305 /* shift up all the later entries in this row */ 306 for ( ii=N; ii>=i; ii-- ) { 307 rp[ii+1] = rp[ii]; 308 PetscMemcpy(ap+bs*bs*(ii+1),ap+bs*bs*(ii),bs*bs); 309 } 310 rp[i] = bcol; 311 ap[bs*bs*i + bs*ridx + cidx] = value; 312 noinsert:; 313 low = i; 314 } 315 ailen[brow] = nrow; 316 } 317 return 0; 318 } 319 320 static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n) 321 { 322 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 323 *m = a->m; *n = a->n; 324 return 0; 325 } 326 327 static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n) 328 { 329 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 330 *m = 0; *n = a->m; 331 return 0; 332 } 333 334 int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 335 { 336 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 337 int itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i; 338 Scalar *aa,*v_i,*aa_i; 339 340 bs = a->bs; 341 ai = a->i; 342 aj = a->j; 343 aa = a->a; 344 345 if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range"); 346 347 bn = row/bs; /* Block number */ 348 bp = row % bs; /* Block Position */ 349 M = ai[bn+1] - ai[bn]; 350 *nz = bs*M; 351 352 if (v) { 353 *v = 0; 354 if (*nz) { 355 *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v); 356 for ( i=0; i<M; i++ ) { /* for each block in the block row */ 357 v_i = *v + i*bs; 358 aa_i = aa + bs*bs*(ai[bn] + i); 359 for ( j=bp,k=0; j<bs*bs; j+=bs,k++ ) {v_i[k] = aa_i[j];} 360 } 361 } 362 } 363 364 if (idx) { 365 *idx = 0; 366 if (*nz) { 367 *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 368 for ( i=0; i<M; i++ ) { /* for each block in the block row */ 369 idx_i = *idx + i*bs; 370 itmp = bs*aj[ai[bn] + i]; 371 for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;} 372 } 373 } 374 } 375 return 0; 376 } 377 378 int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 379 { 380 if (idx) {if (*idx) PetscFree(*idx);} 381 if (v) {if (*v) PetscFree(*v);} 382 return 0; 383 } 384 385 static int MatZeroEntries_SeqBAIJ(Mat A) 386 { 387 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 388 PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar)); 389 return 0; 390 } 391 392 int MatDestroy_SeqBAIJ(PetscObject obj) 393 { 394 Mat A = (Mat) obj; 395 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 396 397 #if defined(PETSC_LOG) 398 PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 399 #endif 400 PetscFree(a->a); 401 if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 402 if (a->diag) PetscFree(a->diag); 403 if (a->ilen) PetscFree(a->ilen); 404 if (a->imax) PetscFree(a->imax); 405 if (a->solve_work) PetscFree(a->solve_work); 406 if (a->mult_work) PetscFree(a->mult_work); 407 PetscFree(a); 408 PLogObjectDestroy(A); 409 PetscHeaderDestroy(A); 410 return 0; 411 } 412 413 static int MatSetOption_SeqBAIJ(Mat A,MatOption op) 414 { 415 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 416 if (op == ROW_ORIENTED) a->roworiented = 1; 417 else if (op == COLUMN_ORIENTED) a->roworiented = 0; 418 else if (op == COLUMNS_SORTED) a->sorted = 1; 419 else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 420 else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 421 else if (op == ROWS_SORTED || 422 op == SYMMETRIC_MATRIX || 423 op == STRUCTURALLY_SYMMETRIC_MATRIX || 424 op == YES_NEW_DIAGONALS) 425 PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n"); 426 else if (op == NO_NEW_DIAGONALS) 427 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");} 428 else 429 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");} 430 return 0; 431 } 432 433 434 /* -------------------------------------------------------*/ 435 /* Should check that shapes of vectors and matrices match */ 436 /* -------------------------------------------------------*/ 437 #include "pinclude/plapack.h" 438 439 static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy) 440 { 441 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 442 Scalar *xg,*yg; 443 register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5; 444 register Scalar x1,x2,x3,x4,x5; 445 int mbs = a->mbs, m = a->m, i, *idx,*ii; 446 int bs = a->bs,j,n,bs2 = bs*bs; 447 448 VecGetArray(xx,&xg); x = xg; VecGetArray(yy,&yg); y = yg; 449 PetscMemzero(y,m*sizeof(Scalar)); 450 x = x; 451 idx = a->j; 452 v = a->a; 453 ii = a->i; 454 455 switch (bs) { 456 case 1: 457 for ( i=0; i<m; i++ ) { 458 n = ii[1] - ii[0]; ii++; 459 sum = 0.0; 460 while (n--) sum += *v++ * x[*idx++]; 461 y[i] = sum; 462 } 463 break; 464 case 2: 465 for ( i=0; i<mbs; i++ ) { 466 n = ii[1] - ii[0]; ii++; 467 sum1 = 0.0; sum2 = 0.0; 468 for ( j=0; j<n; j++ ) { 469 xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 470 sum1 += v[0]*x1 + v[2]*x2; 471 sum2 += v[1]*x1 + v[3]*x2; 472 v += 4; 473 } 474 y[0] += sum1; y[1] += sum2; 475 y += 2; 476 } 477 break; 478 case 3: 479 for ( i=0; i<mbs; i++ ) { 480 n = ii[1] - ii[0]; ii++; 481 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 482 for ( j=0; j<n; j++ ) { 483 xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 484 sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 485 sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 486 sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 487 v += 9; 488 } 489 y[0] += sum1; y[1] += sum2; y[2] += sum3; 490 y += 3; 491 } 492 break; 493 case 4: 494 for ( i=0; i<mbs; i++ ) { 495 n = ii[1] - ii[0]; ii++; 496 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 497 for ( j=0; j<n; j++ ) { 498 xb = x + 4*(*idx++); 499 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 500 sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 501 sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 502 sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 503 sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 504 v += 16; 505 } 506 y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; 507 y += 4; 508 } 509 break; 510 case 5: 511 for ( i=0; i<mbs; i++ ) { 512 n = ii[1] - ii[0]; ii++; 513 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 514 for ( j=0; j<n; j++ ) { 515 xb = x + 5*(*idx++); 516 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 517 sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 518 sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 519 sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 520 sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 521 sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 522 v += 25; 523 } 524 y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5; 525 y += 5; 526 } 527 break; 528 /* block sizes larger then 5 by 5 are handled by BLAS */ 529 default: { 530 int _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt; 531 if (!a->mult_work) { 532 a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar)); 533 CHKPTRQ(a->mult_work); 534 } 535 work = a->mult_work; 536 for ( i=0; i<mbs; i++ ) { 537 n = ii[1] - ii[0]; ii++; 538 ncols = n*bs; 539 workt = work; 540 for ( j=0; j<n; j++ ) { 541 xb = x + bs*(*idx++); 542 for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 543 workt += bs; 544 } 545 LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One); 546 v += n*bs2; 547 y += bs; 548 } 549 } 550 } 551 PLogFlops(2*bs2*a->nz - m); 552 return 0; 553 } 554 555 static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem) 556 { 557 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 558 if (nz) *nz = a->bs*a->bs*a->nz; 559 if (nza) *nza = a->maxnz; 560 if (mem) *mem = (int)A->mem; 561 return 0; 562 } 563 564 static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg) 565 { 566 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 567 568 if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type"); 569 570 /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 571 if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| 572 (a->nz != b->nz)||(a->indexshift != b->indexshift)) { 573 *flg = PETSC_FALSE; return 0; 574 } 575 576 /* if the a->i are the same */ 577 if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) { 578 *flg = PETSC_FALSE; return 0; 579 } 580 581 /* if a->j are the same */ 582 if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) { 583 *flg = PETSC_FALSE; return 0; 584 } 585 586 /* if a->a are the same */ 587 if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) { 588 *flg = PETSC_FALSE; return 0; 589 } 590 *flg = PETSC_TRUE; 591 return 0; 592 593 } 594 595 static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 596 { 597 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 598 int i,j,k,n,row,bs,*ai,*aj,ambs; 599 Scalar *x, zero = 0.0,*aa,*aa_j; 600 601 bs = a->bs; 602 aa = a->a; 603 ai = a->i; 604 aj = a->j; 605 ambs = a->mbs; 606 607 VecSet(&zero,v); 608 VecGetArray(v,&x); VecGetLocalSize(v,&n); 609 if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector"); 610 for ( i=0; i<ambs; i++ ) { 611 for ( j=ai[i]; j<ai[i+1]; j++ ) { 612 if (aj[j] == i) { 613 row = i*bs; 614 aa_j = aa+j*bs*bs; 615 for (k=0; k<bs*bs; k+=(bs+1),row++) x[row] = aa_j[k]; 616 break; 617 } 618 } 619 } 620 return 0; 621 } 622 623 static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 624 { 625 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 626 Scalar *l,*r,x,*v,*aa,*li,*ri; 627 int i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs; 628 629 ai = a->i; 630 aj = a->j; 631 aa = a->a; 632 m = a->m; 633 n = a->n; 634 bs = a->bs; 635 mbs = a->mbs; 636 637 if (ll) { 638 VecGetArray(ll,&l); VecGetSize(ll,&lm); 639 if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length"); 640 for ( i=0; i<mbs; i++ ) { /* for each block row */ 641 M = ai[i+1] - ai[i]; 642 li = l + i*bs; 643 v = aa + bs*bs*ai[i]; 644 for ( j=0; j<M; j++ ) { /* for each block */ 645 for ( k=0; k<bs*bs; k++ ) { 646 (*v++) *= li[k%bs]; 647 } 648 } 649 } 650 } 651 652 if (rr) { 653 VecGetArray(rr,&r); VecGetSize(rr,&rn); 654 if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length"); 655 for ( i=0; i<mbs; i++ ) { /* for each block row */ 656 M = ai[i+1] - ai[i]; 657 v = aa + bs*bs*ai[i]; 658 for ( j=0; j<M; j++ ) { /* for each block */ 659 ri = r + bs*aj[ai[i]+j]; 660 for ( k=0; k<bs; k++ ) { 661 x = ri[k]; 662 for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x; 663 } 664 } 665 } 666 } 667 return 0; 668 } 669 670 671 extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*); 672 extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double); 673 extern int MatSolve_SeqBAIJ(Mat,Vec,Vec); 674 extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*); 675 extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*); 676 extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*); 677 extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*); 678 extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*); 679 extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*); 680 681 static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 682 { 683 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 684 Scalar *v = a->a; 685 double sum = 0.0; 686 int i; 687 688 if (type == NORM_FROBENIUS) { 689 for (i=0; i<a->nz; i++ ) { 690 #if defined(PETSC_COMPLEX) 691 sum += real(conj(*v)*(*v)); v++; 692 #else 693 sum += (*v)*(*v); v++; 694 #endif 695 } 696 *norm = sqrt(sum); 697 } 698 else { 699 SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet"); 700 } 701 return 0; 702 } 703 704 /* 705 note: This can only work for identity for row and col. It would 706 be good to check this and otherwise generate an error. 707 */ 708 static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill) 709 { 710 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 711 Mat outA; 712 int ierr; 713 714 if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported"); 715 716 outA = inA; 717 inA->factor = FACTOR_LU; 718 a->row = row; 719 a->col = col; 720 721 a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work); 722 723 if (!a->diag) { 724 ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr); 725 } 726 ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr); 727 return 0; 728 } 729 730 static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 731 { 732 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 733 int one = 1, totalnz = a->bs*a->bs*a->nz; 734 BLscal_( &totalnz, alpha, a->a, &one ); 735 PLogFlops(totalnz); 736 return 0; 737 } 738 739 int MatPrintHelp_SeqBAIJ(Mat A) 740 { 741 static int called = 0; 742 743 if (called) return 0; else called = 1; 744 return 0; 745 } 746 /* -------------------------------------------------------------------*/ 747 static struct _MatOps MatOps = {MatSetValues_SeqBAIJ, 748 MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ, 749 MatMult_SeqBAIJ,0, 750 0,0, 751 MatSolve_SeqBAIJ,0, 752 0,0, 753 MatLUFactor_SeqBAIJ,0, 754 0, 755 0, 756 MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ, 757 MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ, 758 0,0, 759 0, 760 MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0, 761 MatGetReordering_SeqBAIJ, 762 MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0, 763 MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ, 764 MatILUFactorSymbolic_SeqBAIJ,0, 765 0,0,/* MatConvert_SeqBAIJ */ 0, 766 0,0, 767 MatConvertSameType_SeqBAIJ,0,0, 768 MatILUFactor_SeqBAIJ,0,0, 769 0,0, 770 0,0, 771 MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ, 772 0}; 773 774 /*@C 775 MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format 776 (the default parallel PETSc format). For good matrix assembly performance 777 the user should preallocate the matrix storage by setting the parameter nz 778 (or nzz). By setting these parameters accurately, performance can be 779 increased by more than a factor of 50. 780 781 Input Parameters: 782 . comm - MPI communicator, set to MPI_COMM_SELF 783 . bs - size of block 784 . m - number of rows 785 . n - number of columns 786 . nz - number of block nonzeros per block row (same for all rows) 787 . nzz - number of block nonzeros per block row or PETSC_NULL 788 (possibly different for each row) 789 790 Output Parameter: 791 . A - the matrix 792 793 Notes: 794 The BAIJ format, is fully compatible with standard Fortran 77 795 storage. That is, the stored row and column indices can begin at 796 either one (as in Fortran) or zero. See the users manual for details. 797 798 Specify the preallocated storage with either nz or nnz (not both). 799 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 800 allocation. For additional details, see the users manual chapter on 801 matrices and the file $(PETSC_DIR)/Performance. 802 803 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 804 @*/ 805 int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A) 806 { 807 Mat B; 808 Mat_SeqBAIJ *b; 809 int i,len,ierr, flg,mbs = m/bs; 810 811 if (mbs*bs != m) 812 SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize"); 813 814 *A = 0; 815 PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm); 816 PLogObjectCreate(B); 817 B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b); 818 PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 819 switch (bs) { 820 case 1: 821 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1; 822 break; 823 case 2: 824 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2; 825 break; 826 case 3: 827 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3; 828 break; 829 case 4: 830 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4; 831 break; 832 case 5: 833 B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5; 834 break; 835 } 836 B->destroy = MatDestroy_SeqBAIJ; 837 B->view = MatView_SeqBAIJ; 838 B->factor = 0; 839 B->lupivotthreshold = 1.0; 840 b->row = 0; 841 b->col = 0; 842 b->reallocs = 0; 843 844 b->m = m; 845 b->mbs = mbs; 846 b->n = n; 847 b->imax = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax); 848 if (nnz == PETSC_NULL) { 849 if (nz == PETSC_DEFAULT) nz = 5; 850 else if (nz <= 0) nz = 1; 851 for ( i=0; i<mbs; i++ ) b->imax[i] = nz; 852 nz = nz*mbs; 853 } 854 else { 855 nz = 0; 856 for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 857 } 858 859 /* allocate the matrix space */ 860 len = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int); 861 b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 862 PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar)); 863 b->j = (int *) (b->a + nz*bs*bs); 864 PetscMemzero(b->j,nz*sizeof(int)); 865 b->i = b->j + nz; 866 b->singlemalloc = 1; 867 868 b->i[0] = 0; 869 for (i=1; i<mbs+1; i++) { 870 b->i[i] = b->i[i-1] + b->imax[i-1]; 871 } 872 873 /* b->ilen will count nonzeros in each block row so far. */ 874 b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); 875 PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 876 for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;} 877 878 b->bs = bs; 879 b->mbs = mbs; 880 b->nz = 0; 881 b->maxnz = nz; 882 b->sorted = 0; 883 b->roworiented = 1; 884 b->nonew = 0; 885 b->diag = 0; 886 b->solve_work = 0; 887 b->mult_work = 0; 888 b->spptr = 0; 889 890 *A = B; 891 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 892 if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 893 return 0; 894 } 895 896 int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues) 897 { 898 Mat C; 899 Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data; 900 int i,len, mbs = a->mbs, bs = a->bs,nz = a->nz; 901 902 if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix"); 903 904 *B = 0; 905 PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm); 906 PLogObjectCreate(C); 907 C->data = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c); 908 PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 909 C->destroy = MatDestroy_SeqBAIJ; 910 C->view = MatView_SeqBAIJ; 911 C->factor = A->factor; 912 c->row = 0; 913 c->col = 0; 914 C->assembled = PETSC_TRUE; 915 916 c->m = a->m; 917 c->n = a->n; 918 c->bs = a->bs; 919 c->mbs = a->mbs; 920 c->nbs = a->nbs; 921 922 c->imax = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax); 923 c->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen); 924 for ( i=0; i<mbs; i++ ) { 925 c->imax[i] = a->imax[i]; 926 c->ilen[i] = a->ilen[i]; 927 } 928 929 /* allocate the matrix space */ 930 c->singlemalloc = 1; 931 len = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int)); 932 c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 933 c->j = (int *) (c->a + nz*bs*bs); 934 c->i = c->j + nz; 935 PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int)); 936 if (mbs > 0) { 937 PetscMemcpy(c->j,a->j,nz*sizeof(int)); 938 if (cpvalues == COPY_VALUES) { 939 PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar)); 940 } 941 } 942 943 PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ)); 944 c->sorted = a->sorted; 945 c->roworiented = a->roworiented; 946 c->nonew = a->nonew; 947 948 if (a->diag) { 949 c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag); 950 PLogObjectMemory(C,(mbs+1)*sizeof(int)); 951 for ( i=0; i<mbs; i++ ) { 952 c->diag[i] = a->diag[i]; 953 } 954 } 955 else c->diag = 0; 956 c->nz = a->nz; 957 c->maxnz = a->maxnz; 958 c->solve_work = 0; 959 c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 960 c->mult_work = 0; 961 *B = C; 962 return 0; 963 } 964 965 int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A) 966 { 967 Mat_SeqBAIJ *a; 968 Mat B; 969 int i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg; 970 int *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount; 971 int kmax,jcount,block,idx,point,nzcountb,extra_rows; 972 int *masked, nmask,tmp,ishift,bs2; 973 Scalar *aa; 974 MPI_Comm comm = ((PetscObject) viewer)->comm; 975 976 ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr); 977 bs2 = bs*bs; 978 979 MPI_Comm_size(comm,&size); 980 if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor"); 981 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 982 ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 983 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object"); 984 M = header[1]; N = header[2]; nz = header[3]; 985 986 if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices"); 987 988 /* 989 This code adds extra rows to make sure the number of rows is 990 divisible by the blocksize 991 */ 992 mbs = M/bs; 993 extra_rows = bs - M + bs*(mbs); 994 if (extra_rows == bs) extra_rows = 0; 995 else mbs++; 996 if (extra_rows) { 997 PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize"); 998 } 999 1000 /* read in row lengths */ 1001 rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths); 1002 ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 1003 for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1; 1004 1005 /* read in column indices */ 1006 jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj); 1007 ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr); 1008 for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i; 1009 1010 /* loop over row lengths determining block row lengths */ 1011 browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths); 1012 PetscMemzero(browlengths,mbs*sizeof(int)); 1013 mask = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask); 1014 PetscMemzero(mask,mbs*sizeof(int)); 1015 masked = mask + mbs; 1016 rowcount = 0; nzcount = 0; 1017 for ( i=0; i<mbs; i++ ) { 1018 nmask = 0; 1019 for ( j=0; j<bs; j++ ) { 1020 kmax = rowlengths[rowcount]; 1021 for ( k=0; k<kmax; k++ ) { 1022 tmp = jj[nzcount++]/bs; 1023 if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;} 1024 } 1025 rowcount++; 1026 } 1027 browlengths[i] += nmask; 1028 /* zero out the mask elements we set */ 1029 for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1030 } 1031 1032 /* create our matrix */ 1033 ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A); 1034 CHKERRQ(ierr); 1035 B = *A; 1036 a = (Mat_SeqBAIJ *) B->data; 1037 1038 /* set matrix "i" values */ 1039 a->i[0] = 0; 1040 for ( i=1; i<= mbs; i++ ) { 1041 a->i[i] = a->i[i-1] + browlengths[i-1]; 1042 a->ilen[i-1] = browlengths[i-1]; 1043 } 1044 a->indexshift = 0; 1045 a->nz = 0; 1046 for ( i=0; i<mbs; i++ ) a->nz += browlengths[i]; 1047 1048 /* read in nonzero values */ 1049 aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa); 1050 ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr); 1051 for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0; 1052 1053 /* set "a" and "j" values into matrix */ 1054 nzcount = 0; jcount = 0; 1055 for ( i=0; i<mbs; i++ ) { 1056 nzcountb = nzcount; 1057 nmask = 0; 1058 for ( j=0; j<bs; j++ ) { 1059 kmax = rowlengths[i*bs+j]; 1060 for ( k=0; k<kmax; k++ ) { 1061 tmp = jj[nzcount++]/bs; 1062 if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;} 1063 } 1064 rowcount++; 1065 } 1066 /* sort the masked values */ 1067 PetscSortInt(nmask,masked); 1068 1069 /* set "j" values into matrix */ 1070 maskcount = 1; 1071 for ( j=0; j<nmask; j++ ) { 1072 a->j[jcount++] = masked[j]; 1073 mask[masked[j]] = maskcount++; 1074 } 1075 /* set "a" values into matrix */ 1076 ishift = bs2*a->i[i]; 1077 for ( j=0; j<bs; j++ ) { 1078 kmax = rowlengths[i*bs+j]; 1079 for ( k=0; k<kmax; k++ ) { 1080 tmp = jj[nzcountb]/bs ; 1081 block = mask[tmp] - 1; 1082 point = jj[nzcountb] - bs*tmp; 1083 idx = ishift + bs2*block + j + bs*point; 1084 a->a[idx] = aa[nzcountb++]; 1085 } 1086 } 1087 /* zero out the mask elements we set */ 1088 for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0; 1089 } 1090 if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix"); 1091 1092 PetscFree(rowlengths); 1093 PetscFree(browlengths); 1094 PetscFree(aa); 1095 PetscFree(jj); 1096 PetscFree(mask); 1097 1098 B->assembled = PETSC_TRUE; 1099 1100 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 1101 if (flg) { 1102 Viewer tviewer; 1103 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1104 ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr); 1105 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1106 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1107 } 1108 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 1109 if (flg) { 1110 Viewer tviewer; 1111 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1112 ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr); 1113 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1114 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1115 } 1116 ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr); 1117 if (flg) { 1118 Viewer tviewer; 1119 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1120 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1121 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1122 } 1123 ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr); 1124 if (flg) { 1125 Viewer tviewer; 1126 ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr); 1127 ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr); 1128 ierr = MatView(B,tviewer); CHKERRQ(ierr); 1129 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1130 } 1131 ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 1132 if (flg) { 1133 Viewer tviewer; 1134 ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr); 1135 ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr); 1136 ierr = ViewerFlush(tviewer); CHKERRQ(ierr); 1137 ierr = ViewerDestroy(tviewer); CHKERRQ(ierr); 1138 } 1139 return 0; 1140 } 1141 1142 1143 1144