1 2 #ifndef lint 3 static char vcid[] = "$Id: aij.c,v 1.170 1996/04/26 00:02:52 balay Exp bsmith $"; 4 #endif 5 6 /* 7 Defines the basic matrix operations for the AIJ (compressed row) 8 matrix storage format. 9 */ 10 #include "aij.h" 11 #include "src/vec/vecimpl.h" 12 #include "src/inline/spops.h" 13 #include "petsc.h" 14 #include "src/inline/bitarray.h" 15 16 extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 17 18 static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 19 { 20 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 21 int ierr, *ia, *ja,n,*idx,i,oshift,ishift; 22 23 /* 24 this is tacky: In the future when we have written special factorization 25 and solve routines for the identity permutation we should use a 26 stride index set instead of the general one. 27 */ 28 if (type == ORDER_NATURAL) { 29 n = a->n; 30 idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 31 for ( i=0; i<n; i++ ) idx[i] = i; 32 ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 33 ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 34 PetscFree(idx); 35 ISSetPermutation(*rperm); 36 ISSetPermutation(*cperm); 37 ISSetIdentity(*rperm); 38 ISSetIdentity(*cperm); 39 return 0; 40 } 41 42 MatReorderingRegisterAll(); 43 ishift = a->indexshift; 44 oshift = -MatReorderingIndexShift[(int)type]; 45 if (MatReorderingRequiresSymmetric[(int)type]) { 46 ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja); 47 CHKERRQ(ierr); 48 ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 49 PetscFree(ia); PetscFree(ja); 50 } else { 51 if (ishift == oshift) { 52 ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 53 } 54 else if (ishift == -1) { 55 /* temporarily subtract 1 from i and j indices */ 56 int nz = a->i[a->n] - 1; 57 for ( i=0; i<nz; i++ ) a->j[i]--; 58 for ( i=0; i<a->n+1; i++ ) a->i[i]--; 59 ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 60 for ( i=0; i<nz; i++ ) a->j[i]++; 61 for ( i=0; i<a->n+1; i++ ) a->i[i]++; 62 } else { 63 /* temporarily add 1 to i and j indices */ 64 int nz = a->i[a->n] - 1; 65 for ( i=0; i<nz; i++ ) a->j[i]++; 66 for ( i=0; i<a->n+1; i++ ) a->i[i]++; 67 ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 68 for ( i=0; i<nz; i++ ) a->j[i]--; 69 for ( i=0; i<a->n+1; i++ ) a->i[i]--; 70 } 71 } 72 return 0; 73 } 74 75 #define CHUNKSIZE 15 76 77 /* This version has row oriented v */ 78 static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 79 { 80 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 81 int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 82 int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 83 int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 84 Scalar *ap,value, *aa = a->a; 85 86 for ( k=0; k<m; k++ ) { /* loop over added rows */ 87 row = im[k]; 88 if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 89 if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 90 rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 91 rmax = imax[row]; nrow = ailen[row]; 92 low = 0; 93 for ( l=0; l<n; l++ ) { /* loop over added columns */ 94 if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 95 if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 96 col = in[l] - shift; 97 if (roworiented) { 98 value = *v++; 99 } 100 else { 101 value = v[k + l*m]; 102 } 103 if (!sorted) low = 0; high = nrow; 104 while (high-low > 5) { 105 t = (low+high)/2; 106 if (rp[t] > col) high = t; 107 else low = t; 108 } 109 for ( i=low; i<high; i++ ) { 110 if (rp[i] > col) break; 111 if (rp[i] == col) { 112 if (is == ADD_VALUES) ap[i] += value; 113 else ap[i] = value; 114 goto noinsert; 115 } 116 } 117 if (nonew) goto noinsert; 118 if (nrow >= rmax) { 119 /* there is no extra room in row, therefore enlarge */ 120 int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 121 Scalar *new_a; 122 123 /* malloc new storage space */ 124 len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 125 new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 126 new_j = (int *) (new_a + new_nz); 127 new_i = new_j + new_nz; 128 129 /* copy over old data into new slots */ 130 for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 131 for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 132 PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 133 len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 134 PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 135 len*sizeof(int)); 136 PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 137 PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 138 len*sizeof(Scalar)); 139 /* free up old matrix storage */ 140 PetscFree(a->a); 141 if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 142 aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 143 a->singlemalloc = 1; 144 145 rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 146 rmax = imax[row] = imax[row] + CHUNKSIZE; 147 PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 148 a->maxnz += CHUNKSIZE; 149 a->reallocs++; 150 } 151 N = nrow++ - 1; a->nz++; 152 /* shift up all the later entries in this row */ 153 for ( ii=N; ii>=i; ii-- ) { 154 rp[ii+1] = rp[ii]; 155 ap[ii+1] = ap[ii]; 156 } 157 rp[i] = col; 158 ap[i] = value; 159 noinsert:; 160 low = i + 1; 161 } 162 ailen[row] = nrow; 163 } 164 return 0; 165 } 166 167 static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 168 { 169 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 170 int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 171 int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 172 Scalar *ap, *aa = a->a, zero = 0.0; 173 174 for ( k=0; k<m; k++ ) { /* loop over rows */ 175 row = im[k]; 176 if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 177 if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 178 rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 179 nrow = ailen[row]; 180 for ( l=0; l<n; l++ ) { /* loop over columns */ 181 if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 182 if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 183 col = in[l] - shift; 184 high = nrow; low = 0; /* assume unsorted */ 185 while (high-low > 5) { 186 t = (low+high)/2; 187 if (rp[t] > col) high = t; 188 else low = t; 189 } 190 for ( i=low; i<high; i++ ) { 191 if (rp[i] > col) break; 192 if (rp[i] == col) { 193 *v++ = ap[i]; 194 goto finished; 195 } 196 } 197 *v++ = zero; 198 finished:; 199 } 200 } 201 return 0; 202 } 203 204 #include "draw.h" 205 #include "pinclude/pviewer.h" 206 #include "sys.h" 207 208 static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 209 { 210 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 211 int i, fd, *col_lens, ierr; 212 213 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 214 col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 215 col_lens[0] = MAT_COOKIE; 216 col_lens[1] = a->m; 217 col_lens[2] = a->n; 218 col_lens[3] = a->nz; 219 220 /* store lengths of each row and write (including header) to file */ 221 for ( i=0; i<a->m; i++ ) { 222 col_lens[4+i] = a->i[i+1] - a->i[i]; 223 } 224 ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 225 PetscFree(col_lens); 226 227 /* store column indices (zero start index) */ 228 if (a->indexshift) { 229 for ( i=0; i<a->nz; i++ ) a->j[i]--; 230 } 231 ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 232 if (a->indexshift) { 233 for ( i=0; i<a->nz; i++ ) a->j[i]++; 234 } 235 236 /* store nonzero values */ 237 ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 238 return 0; 239 } 240 241 static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 242 { 243 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 244 int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 245 FILE *fd; 246 char *outputname; 247 248 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 249 ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 250 ierr = ViewerGetFormat(viewer,&format); 251 if (format == ASCII_FORMAT_INFO) { 252 return 0; 253 } 254 else if (format == ASCII_FORMAT_INFO_DETAILED) { 255 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 256 ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 257 if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 258 else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 259 a->inode.node_count,a->inode.limit); 260 } 261 else if (format == ASCII_FORMAT_MATLAB) { 262 int nz, nzalloc, mem; 263 MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 264 fprintf(fd,"%% Size = %d %d \n",m,a->n); 265 fprintf(fd,"%% Nonzeros = %d \n",nz); 266 fprintf(fd,"zzz = zeros(%d,3);\n",nz); 267 fprintf(fd,"zzz = [\n"); 268 269 for (i=0; i<m; i++) { 270 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 271 #if defined(PETSC_COMPLEX) 272 fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]), 273 imag(a->a[j])); 274 #else 275 fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 276 #endif 277 } 278 } 279 fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 280 } 281 else if (format == ASCII_FORMAT_COMMON) { 282 for ( i=0; i<m; i++ ) { 283 fprintf(fd,"row %d:",i); 284 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 285 #if defined(PETSC_COMPLEX) 286 if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0) 287 fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 288 else if (real(a->a[j]) != 0.0) 289 fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 290 #else 291 if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 292 #endif 293 } 294 fprintf(fd,"\n"); 295 } 296 } 297 else { 298 for ( i=0; i<m; i++ ) { 299 fprintf(fd,"row %d:",i); 300 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 301 #if defined(PETSC_COMPLEX) 302 if (imag(a->a[j]) != 0.0) { 303 fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 304 } 305 else { 306 fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 307 } 308 #else 309 fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 310 #endif 311 } 312 fprintf(fd,"\n"); 313 } 314 } 315 fflush(fd); 316 return 0; 317 } 318 319 static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 320 { 321 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 322 int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 323 double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 324 Draw draw; 325 DrawButton button; 326 PetscTruth isnull; 327 328 ViewerDrawGetDraw(viewer,&draw); 329 ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 330 331 xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 332 xr += w; yr += h; xl = -w; yl = -h; 333 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 334 /* loop over matrix elements drawing boxes */ 335 color = DRAW_BLUE; 336 for ( i=0; i<m; i++ ) { 337 y_l = m - i - 1.0; y_r = y_l + 1.0; 338 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 339 x_l = a->j[j] + shift; x_r = x_l + 1.0; 340 #if defined(PETSC_COMPLEX) 341 if (real(a->a[j]) >= 0.) continue; 342 #else 343 if (a->a[j] >= 0.) continue; 344 #endif 345 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 346 } 347 } 348 color = DRAW_CYAN; 349 for ( i=0; i<m; i++ ) { 350 y_l = m - i - 1.0; y_r = y_l + 1.0; 351 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 352 x_l = a->j[j] + shift; x_r = x_l + 1.0; 353 if (a->a[j] != 0.) continue; 354 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 355 } 356 } 357 color = DRAW_RED; 358 for ( i=0; i<m; i++ ) { 359 y_l = m - i - 1.0; y_r = y_l + 1.0; 360 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 361 x_l = a->j[j] + shift; x_r = x_l + 1.0; 362 #if defined(PETSC_COMPLEX) 363 if (real(a->a[j]) <= 0.) continue; 364 #else 365 if (a->a[j] <= 0.) continue; 366 #endif 367 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 368 } 369 } 370 DrawFlush(draw); 371 DrawGetPause(draw,&pause); 372 if (pause >= 0) { PetscSleep(pause); return 0;} 373 374 /* allow the matrix to zoom or shrink */ 375 ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 376 while (button != BUTTON_RIGHT) { 377 DrawClear(draw); 378 if (button == BUTTON_LEFT) scale = .5; 379 else if (button == BUTTON_CENTER) scale = 2.; 380 xl = scale*(xl + w - xc) + xc - w*scale; 381 xr = scale*(xr - w - xc) + xc + w*scale; 382 yl = scale*(yl + h - yc) + yc - h*scale; 383 yr = scale*(yr - h - yc) + yc + h*scale; 384 w *= scale; h *= scale; 385 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 386 color = DRAW_BLUE; 387 for ( i=0; i<m; i++ ) { 388 y_l = m - i - 1.0; y_r = y_l + 1.0; 389 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 390 x_l = a->j[j] + shift; x_r = x_l + 1.0; 391 #if defined(PETSC_COMPLEX) 392 if (real(a->a[j]) >= 0.) continue; 393 #else 394 if (a->a[j] >= 0.) continue; 395 #endif 396 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 397 } 398 } 399 color = DRAW_CYAN; 400 for ( i=0; i<m; i++ ) { 401 y_l = m - i - 1.0; y_r = y_l + 1.0; 402 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 403 x_l = a->j[j] + shift; x_r = x_l + 1.0; 404 if (a->a[j] != 0.) continue; 405 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 406 } 407 } 408 color = DRAW_RED; 409 for ( i=0; i<m; i++ ) { 410 y_l = m - i - 1.0; y_r = y_l + 1.0; 411 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 412 x_l = a->j[j] + shift; x_r = x_l + 1.0; 413 #if defined(PETSC_COMPLEX) 414 if (real(a->a[j]) <= 0.) continue; 415 #else 416 if (a->a[j] <= 0.) continue; 417 #endif 418 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 419 } 420 } 421 ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 422 } 423 return 0; 424 } 425 426 static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 427 { 428 Mat A = (Mat) obj; 429 Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 430 ViewerType vtype; 431 int ierr; 432 433 if (!viewer) { 434 viewer = STDOUT_VIEWER_SELF; 435 } 436 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 437 if (vtype == MATLAB_VIEWER) { 438 return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 439 } 440 else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 441 return MatView_SeqAIJ_ASCII(A,viewer); 442 } 443 else if (vtype == BINARY_FILE_VIEWER) { 444 return MatView_SeqAIJ_Binary(A,viewer); 445 } 446 else if (vtype == DRAW_VIEWER) { 447 return MatView_SeqAIJ_Draw(A,viewer); 448 } 449 return 0; 450 } 451 452 extern int Mat_AIJ_CheckInode(Mat); 453 static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 454 { 455 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 456 int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 457 int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 458 Scalar *aa = a->a, *ap; 459 460 if (mode == FLUSH_ASSEMBLY) return 0; 461 462 for ( i=1; i<m; i++ ) { 463 /* move each row back by the amount of empty slots (fshift) before it*/ 464 fshift += imax[i-1] - ailen[i-1]; 465 if (fshift) { 466 ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 467 N = ailen[i]; 468 for ( j=0; j<N; j++ ) { 469 ip[j-fshift] = ip[j]; 470 ap[j-fshift] = ap[j]; 471 } 472 } 473 ai[i] = ai[i-1] + ailen[i-1]; 474 } 475 if (m) { 476 fshift += imax[m-1] - ailen[m-1]; 477 ai[m] = ai[m-1] + ailen[m-1]; 478 } 479 /* reset ilen and imax for each row */ 480 for ( i=0; i<m; i++ ) { 481 ailen[i] = imax[i] = ai[i+1] - ai[i]; 482 } 483 a->nz = ai[m] + shift; 484 485 /* diagonals may have moved, so kill the diagonal pointers */ 486 if (fshift && a->diag) { 487 PetscFree(a->diag); 488 PLogObjectMemory(A,-(m+1)*sizeof(int)); 489 a->diag = 0; 490 } 491 PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n", 492 fshift,a->nz,m); 493 PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n", 494 a->reallocs); 495 /* check out for identical nodes. If found, use inode functions */ 496 ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 497 return 0; 498 } 499 500 static int MatZeroEntries_SeqAIJ(Mat A) 501 { 502 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 503 PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 504 return 0; 505 } 506 507 int MatDestroy_SeqAIJ(PetscObject obj) 508 { 509 Mat A = (Mat) obj; 510 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 511 512 #if defined(PETSC_LOG) 513 PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 514 #endif 515 PetscFree(a->a); 516 if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 517 if (a->diag) PetscFree(a->diag); 518 if (a->ilen) PetscFree(a->ilen); 519 if (a->imax) PetscFree(a->imax); 520 if (a->solve_work) PetscFree(a->solve_work); 521 if (a->inode.size) PetscFree(a->inode.size); 522 PetscFree(a); 523 PLogObjectDestroy(A); 524 PetscHeaderDestroy(A); 525 return 0; 526 } 527 528 static int MatCompress_SeqAIJ(Mat A) 529 { 530 return 0; 531 } 532 533 static int MatSetOption_SeqAIJ(Mat A,MatOption op) 534 { 535 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 536 if (op == ROW_ORIENTED) a->roworiented = 1; 537 else if (op == COLUMN_ORIENTED) a->roworiented = 0; 538 else if (op == COLUMNS_SORTED) a->sorted = 1; 539 else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 540 else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 541 else if (op == ROWS_SORTED || 542 op == SYMMETRIC_MATRIX || 543 op == STRUCTURALLY_SYMMETRIC_MATRIX || 544 op == YES_NEW_DIAGONALS) 545 PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 546 else if (op == NO_NEW_DIAGONALS) 547 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 548 else if (op == INODE_LIMIT_1) a->inode.limit = 1; 549 else if (op == INODE_LIMIT_2) a->inode.limit = 2; 550 else if (op == INODE_LIMIT_3) a->inode.limit = 3; 551 else if (op == INODE_LIMIT_4) a->inode.limit = 4; 552 else if (op == INODE_LIMIT_5) a->inode.limit = 5; 553 else 554 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 555 return 0; 556 } 557 558 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 559 { 560 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 561 int i,j, n,shift = a->indexshift; 562 Scalar *x, zero = 0.0; 563 564 VecSet(&zero,v); 565 VecGetArray(v,&x); VecGetLocalSize(v,&n); 566 if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 567 for ( i=0; i<a->m; i++ ) { 568 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 569 if (a->j[j]+shift == i) { 570 x[i] = a->a[j]; 571 break; 572 } 573 } 574 } 575 return 0; 576 } 577 578 /* -------------------------------------------------------*/ 579 /* Should check that shapes of vectors and matrices match */ 580 /* -------------------------------------------------------*/ 581 int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 582 { 583 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 584 Scalar *x, *y, *v, alpha; 585 int m = a->m, n, i, *idx, shift = a->indexshift; 586 587 VecGetArray(xx,&x); VecGetArray(yy,&y); 588 PetscMemzero(y,a->n*sizeof(Scalar)); 589 y = y + shift; /* shift for Fortran start by 1 indexing */ 590 for ( i=0; i<m; i++ ) { 591 idx = a->j + a->i[i] + shift; 592 v = a->a + a->i[i] + shift; 593 n = a->i[i+1] - a->i[i]; 594 alpha = x[i]; 595 while (n-->0) {y[*idx++] += alpha * *v++;} 596 } 597 PLogFlops(2*a->nz - a->n); 598 return 0; 599 } 600 601 int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 602 { 603 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 604 Scalar *x, *y, *v, alpha; 605 int m = a->m, n, i, *idx,shift = a->indexshift; 606 607 VecGetArray(xx,&x); VecGetArray(yy,&y); 608 if (zz != yy) VecCopy(zz,yy); 609 y = y + shift; /* shift for Fortran start by 1 indexing */ 610 for ( i=0; i<m; i++ ) { 611 idx = a->j + a->i[i] + shift; 612 v = a->a + a->i[i] + shift; 613 n = a->i[i+1] - a->i[i]; 614 alpha = x[i]; 615 while (n-->0) {y[*idx++] += alpha * *v++;} 616 } 617 return 0; 618 } 619 620 int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 621 { 622 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 623 Scalar *x, *y, *v, sum; 624 int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 625 626 VecGetArray(xx,&x); VecGetArray(yy,&y); 627 x = x + shift; /* shift for Fortran start by 1 indexing */ 628 idx = a->j; 629 v = a->a; 630 ii = a->i; 631 for ( i=0; i<m; i++ ) { 632 n = ii[1] - ii[0]; ii++; 633 sum = 0.0; 634 /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 635 /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 636 while (n--) sum += *v++ * x[*idx++]; 637 y[i] = sum; 638 } 639 PLogFlops(2*a->nz - m); 640 return 0; 641 } 642 643 int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 644 { 645 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 646 Scalar *x, *y, *z, *v, sum; 647 int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 648 649 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 650 x = x + shift; /* shift for Fortran start by 1 indexing */ 651 idx = a->j; 652 v = a->a; 653 ii = a->i; 654 for ( i=0; i<m; i++ ) { 655 n = ii[1] - ii[0]; ii++; 656 sum = y[i]; 657 /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 658 while (n--) sum += *v++ * x[*idx++]; 659 z[i] = sum; 660 } 661 PLogFlops(2*a->nz); 662 return 0; 663 } 664 665 /* 666 Adds diagonal pointers to sparse matrix structure. 667 */ 668 669 int MatMarkDiag_SeqAIJ(Mat A) 670 { 671 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 672 int i,j, *diag, m = a->m,shift = a->indexshift; 673 674 diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 675 PLogObjectMemory(A,(m+1)*sizeof(int)); 676 for ( i=0; i<a->m; i++ ) { 677 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 678 if (a->j[j]+shift == i) { 679 diag[i] = j - shift; 680 break; 681 } 682 } 683 } 684 a->diag = diag; 685 return 0; 686 } 687 688 int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 689 double fshift,int its,Vec xx) 690 { 691 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 692 Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 693 int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 694 695 VecGetArray(xx,&x); VecGetArray(bb,&b); 696 if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 697 diag = a->diag; 698 xs = x + shift; /* shifted by one for index start of a or a->j*/ 699 if (flag == SOR_APPLY_UPPER) { 700 /* apply ( U + D/omega) to the vector */ 701 bs = b + shift; 702 for ( i=0; i<m; i++ ) { 703 d = fshift + a->a[diag[i] + shift]; 704 n = a->i[i+1] - diag[i] - 1; 705 idx = a->j + diag[i] + (!shift); 706 v = a->a + diag[i] + (!shift); 707 sum = b[i]*d/omega; 708 SPARSEDENSEDOT(sum,bs,v,idx,n); 709 x[i] = sum; 710 } 711 return 0; 712 } 713 if (flag == SOR_APPLY_LOWER) { 714 SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 715 } 716 else if (flag & SOR_EISENSTAT) { 717 /* Let A = L + U + D; where L is lower trianglar, 718 U is upper triangular, E is diagonal; This routine applies 719 720 (L + E)^{-1} A (U + E)^{-1} 721 722 to a vector efficiently using Eisenstat's trick. This is for 723 the case of SSOR preconditioner, so E is D/omega where omega 724 is the relaxation factor. 725 */ 726 t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 727 scale = (2.0/omega) - 1.0; 728 729 /* x = (E + U)^{-1} b */ 730 for ( i=m-1; i>=0; i-- ) { 731 d = fshift + a->a[diag[i] + shift]; 732 n = a->i[i+1] - diag[i] - 1; 733 idx = a->j + diag[i] + (!shift); 734 v = a->a + diag[i] + (!shift); 735 sum = b[i]; 736 SPARSEDENSEMDOT(sum,xs,v,idx,n); 737 x[i] = omega*(sum/d); 738 } 739 740 /* t = b - (2*E - D)x */ 741 v = a->a; 742 for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 743 744 /* t = (E + L)^{-1}t */ 745 ts = t + shift; /* shifted by one for index start of a or a->j*/ 746 diag = a->diag; 747 for ( i=0; i<m; i++ ) { 748 d = fshift + a->a[diag[i]+shift]; 749 n = diag[i] - a->i[i]; 750 idx = a->j + a->i[i] + shift; 751 v = a->a + a->i[i] + shift; 752 sum = t[i]; 753 SPARSEDENSEMDOT(sum,ts,v,idx,n); 754 t[i] = omega*(sum/d); 755 } 756 757 /* x = x + t */ 758 for ( i=0; i<m; i++ ) { x[i] += t[i]; } 759 PetscFree(t); 760 return 0; 761 } 762 if (flag & SOR_ZERO_INITIAL_GUESS) { 763 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 764 for ( i=0; i<m; i++ ) { 765 d = fshift + a->a[diag[i]+shift]; 766 n = diag[i] - a->i[i]; 767 idx = a->j + a->i[i] + shift; 768 v = a->a + a->i[i] + shift; 769 sum = b[i]; 770 SPARSEDENSEMDOT(sum,xs,v,idx,n); 771 x[i] = omega*(sum/d); 772 } 773 xb = x; 774 } 775 else xb = b; 776 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 777 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 778 for ( i=0; i<m; i++ ) { 779 x[i] *= a->a[diag[i]+shift]; 780 } 781 } 782 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 783 for ( i=m-1; i>=0; i-- ) { 784 d = fshift + a->a[diag[i] + shift]; 785 n = a->i[i+1] - diag[i] - 1; 786 idx = a->j + diag[i] + (!shift); 787 v = a->a + diag[i] + (!shift); 788 sum = xb[i]; 789 SPARSEDENSEMDOT(sum,xs,v,idx,n); 790 x[i] = omega*(sum/d); 791 } 792 } 793 its--; 794 } 795 while (its--) { 796 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 797 for ( i=0; i<m; i++ ) { 798 d = fshift + a->a[diag[i]+shift]; 799 n = a->i[i+1] - a->i[i]; 800 idx = a->j + a->i[i] + shift; 801 v = a->a + a->i[i] + shift; 802 sum = b[i]; 803 SPARSEDENSEMDOT(sum,xs,v,idx,n); 804 x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 805 } 806 } 807 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 808 for ( i=m-1; i>=0; i-- ) { 809 d = fshift + a->a[diag[i] + shift]; 810 n = a->i[i+1] - a->i[i]; 811 idx = a->j + a->i[i] + shift; 812 v = a->a + a->i[i] + shift; 813 sum = b[i]; 814 SPARSEDENSEMDOT(sum,xs,v,idx,n); 815 x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 816 } 817 } 818 } 819 return 0; 820 } 821 822 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 823 { 824 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 825 if (nz) *nz = a->nz; 826 if (nzalloc) *nzalloc = a->maxnz; 827 if (mem) *mem = (int)A->mem; 828 return 0; 829 } 830 831 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 832 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 833 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 834 extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 835 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 836 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 837 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 838 839 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 840 { 841 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 842 int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 843 844 ierr = ISGetSize(is,&N); CHKERRQ(ierr); 845 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 846 if (diag) { 847 for ( i=0; i<N; i++ ) { 848 if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 849 if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 850 a->ilen[rows[i]] = 1; 851 a->a[a->i[rows[i]]+shift] = *diag; 852 a->j[a->i[rows[i]]+shift] = rows[i]+shift; 853 } 854 else { 855 ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 856 CHKERRQ(ierr); 857 } 858 } 859 } 860 else { 861 for ( i=0; i<N; i++ ) { 862 if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 863 a->ilen[rows[i]] = 0; 864 } 865 } 866 ISRestoreIndices(is,&rows); 867 ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 868 ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 869 return 0; 870 } 871 872 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 873 { 874 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 875 *m = a->m; *n = a->n; 876 return 0; 877 } 878 879 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 880 { 881 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 882 *m = 0; *n = a->m; 883 return 0; 884 } 885 int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 886 { 887 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 888 int *itmp,i,shift = a->indexshift; 889 890 if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 891 892 *nz = a->i[row+1] - a->i[row]; 893 if (v) *v = a->a + a->i[row] + shift; 894 if (idx) { 895 itmp = a->j + a->i[row] + shift; 896 if (*nz && shift) { 897 *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 898 for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 899 } else if (*nz) { 900 *idx = itmp; 901 } 902 else *idx = 0; 903 } 904 return 0; 905 } 906 907 int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 908 { 909 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 910 if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 911 return 0; 912 } 913 914 static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 915 { 916 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 917 Scalar *v = a->a; 918 double sum = 0.0; 919 int i, j,shift = a->indexshift; 920 921 if (type == NORM_FROBENIUS) { 922 for (i=0; i<a->nz; i++ ) { 923 #if defined(PETSC_COMPLEX) 924 sum += real(conj(*v)*(*v)); v++; 925 #else 926 sum += (*v)*(*v); v++; 927 #endif 928 } 929 *norm = sqrt(sum); 930 } 931 else if (type == NORM_1) { 932 double *tmp; 933 int *jj = a->j; 934 tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 935 PetscMemzero(tmp,a->n*sizeof(double)); 936 *norm = 0.0; 937 for ( j=0; j<a->nz; j++ ) { 938 tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 939 } 940 for ( j=0; j<a->n; j++ ) { 941 if (tmp[j] > *norm) *norm = tmp[j]; 942 } 943 PetscFree(tmp); 944 } 945 else if (type == NORM_INFINITY) { 946 *norm = 0.0; 947 for ( j=0; j<a->m; j++ ) { 948 v = a->a + a->i[j] + shift; 949 sum = 0.0; 950 for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 951 sum += PetscAbsScalar(*v); v++; 952 } 953 if (sum > *norm) *norm = sum; 954 } 955 } 956 else { 957 SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 958 } 959 return 0; 960 } 961 962 static int MatTranspose_SeqAIJ(Mat A,Mat *B) 963 { 964 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 965 Mat C; 966 int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 967 int shift = a->indexshift; 968 Scalar *array = a->a; 969 970 if (B == PETSC_NULL && m != a->n) 971 SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place"); 972 col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 973 PetscMemzero(col,(1+a->n)*sizeof(int)); 974 if (shift) { 975 for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 976 } 977 for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 978 ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 979 PetscFree(col); 980 for ( i=0; i<m; i++ ) { 981 len = ai[i+1]-ai[i]; 982 ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 983 array += len; aj += len; 984 } 985 if (shift) { 986 for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 987 } 988 989 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 990 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 991 992 if (B != PETSC_NULL) { 993 *B = C; 994 } else { 995 /* This isn't really an in-place transpose */ 996 PetscFree(a->a); 997 if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 998 if (a->diag) PetscFree(a->diag); 999 if (a->ilen) PetscFree(a->ilen); 1000 if (a->imax) PetscFree(a->imax); 1001 if (a->solve_work) PetscFree(a->solve_work); 1002 if (a->inode.size) PetscFree(a->inode.size); 1003 PetscFree(a); 1004 PetscMemcpy(A,C,sizeof(struct _Mat)); 1005 PetscHeaderDestroy(C); 1006 } 1007 return 0; 1008 } 1009 1010 static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 1011 { 1012 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1013 Scalar *l,*r,x,*v; 1014 int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 1015 1016 if (ll) { 1017 VecGetArray(ll,&l); VecGetSize(ll,&m); 1018 if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 1019 v = a->a; 1020 for ( i=0; i<m; i++ ) { 1021 x = l[i]; 1022 M = a->i[i+1] - a->i[i]; 1023 for ( j=0; j<M; j++ ) { (*v++) *= x;} 1024 } 1025 PLogFlops(nz); 1026 } 1027 if (rr) { 1028 VecGetArray(rr,&r); VecGetSize(rr,&n); 1029 if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 1030 v = a->a; jj = a->j; 1031 for ( i=0; i<nz; i++ ) { 1032 (*v++) *= r[*jj++ + shift]; 1033 } 1034 PLogFlops(nz); 1035 } 1036 return 0; 1037 } 1038 1039 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 1040 { 1041 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 1042 int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 1043 int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1044 register int sum,lensi; 1045 int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 1046 int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 1047 Scalar *a_new,*mat_a; 1048 Mat C; 1049 1050 ierr = ISSorted(iscol,(PetscTruth*)&i); 1051 if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IS is not sorted"); 1052 1053 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 1054 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 1055 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 1056 1057 if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 1058 /* special case of contiguous rows */ 1059 lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 1060 starts = lens + ncols; 1061 /* loop over new rows determining lens and starting points */ 1062 for (i=0; i<nrows; i++) { 1063 kstart = ai[irow[i]]+shift; 1064 kend = kstart + ailen[irow[i]]; 1065 for ( k=kstart; k<kend; k++ ) { 1066 if (aj[k]+shift >= first) { 1067 starts[i] = k; 1068 break; 1069 } 1070 } 1071 sum = 0; 1072 while (k < kend) { 1073 if (aj[k++]+shift >= first+ncols) break; 1074 sum++; 1075 } 1076 lens[i] = sum; 1077 } 1078 /* create submatrix */ 1079 if (scall == MAT_REUSE_MATRIX) { 1080 int n_cols,n_rows; 1081 ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1082 if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1083 ierr = MatZeroEntries(*B); CHKERRQ(ierr); 1084 C = *B; 1085 } 1086 else { 1087 ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 1088 } 1089 c = (Mat_SeqAIJ*) C->data; 1090 1091 /* loop over rows inserting into submatrix */ 1092 a_new = c->a; 1093 j_new = c->j; 1094 i_new = c->i; 1095 i_new[0] = -shift; 1096 for (i=0; i<nrows; i++) { 1097 ii = starts[i]; 1098 lensi = lens[i]; 1099 for ( k=0; k<lensi; k++ ) { 1100 *j_new++ = aj[ii+k] - first; 1101 } 1102 PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1103 a_new += lensi; 1104 i_new[i+1] = i_new[i] + lensi; 1105 c->ilen[i] = lensi; 1106 } 1107 PetscFree(lens); 1108 } 1109 else { 1110 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 1111 smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 1112 ssmap = smap + shift; 1113 lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1114 PetscMemzero(smap,oldcols*sizeof(int)); 1115 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 1116 /* determine lens of each row */ 1117 for (i=0; i<nrows; i++) { 1118 kstart = ai[irow[i]]+shift; 1119 kend = kstart + a->ilen[irow[i]]; 1120 lens[i] = 0; 1121 for ( k=kstart; k<kend; k++ ) { 1122 if (ssmap[aj[k]]) { 1123 lens[i]++; 1124 } 1125 } 1126 } 1127 /* Create and fill new matrix */ 1128 if (scall == MAT_REUSE_MATRIX) { 1129 c = (Mat_SeqAIJ *)((*B)->data); 1130 1131 if (c->m != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1132 if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1133 SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros"); 1134 } 1135 PetscMemzero(c->ilen,c->m*sizeof(int)); 1136 C = *B; 1137 } 1138 else { 1139 ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 1140 } 1141 c = (Mat_SeqAIJ *)(C->data); 1142 for (i=0; i<nrows; i++) { 1143 row = irow[i]; 1144 nznew = 0; 1145 kstart = ai[row]+shift; 1146 kend = kstart + a->ilen[row]; 1147 mat_i = c->i[i]+shift; 1148 mat_j = c->j + mat_i; 1149 mat_a = c->a + mat_i; 1150 mat_ilen = c->ilen + i; 1151 for ( k=kstart; k<kend; k++ ) { 1152 if ((tcol=ssmap[a->j[k]])) { 1153 *mat_j++ = tcol - (!shift); 1154 *mat_a++ = a->a[k]; 1155 (*mat_ilen)++; 1156 1157 } 1158 } 1159 } 1160 /* Free work space */ 1161 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 1162 PetscFree(smap); PetscFree(lens); 1163 } 1164 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1165 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1166 1167 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1168 *B = C; 1169 return 0; 1170 } 1171 1172 /* 1173 note: This can only work for identity for row and col. It would 1174 be good to check this and otherwise generate an error. 1175 */ 1176 static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1177 { 1178 Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1179 int ierr; 1180 Mat outA; 1181 1182 if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1183 1184 outA = inA; 1185 inA->factor = FACTOR_LU; 1186 a->row = row; 1187 a->col = col; 1188 1189 a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 1190 1191 if (!a->diag) { 1192 ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 1193 } 1194 ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1195 return 0; 1196 } 1197 1198 #include "pinclude/plapack.h" 1199 static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1200 { 1201 Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1202 int one = 1; 1203 BLscal_( &a->nz, alpha, a->a, &one ); 1204 PLogFlops(a->nz); 1205 return 0; 1206 } 1207 1208 static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1209 Mat **B) 1210 { 1211 int ierr,i; 1212 1213 if (scall == MAT_INITIAL_MATRIX) { 1214 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1215 } 1216 1217 for ( i=0; i<n; i++ ) { 1218 ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1219 } 1220 return 0; 1221 } 1222 1223 static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 1224 { 1225 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1226 int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 1227 int start, end, *ai, *aj; 1228 char *table; 1229 shift = a->indexshift; 1230 m = a->m; 1231 ai = a->i; 1232 aj = a->j+shift; 1233 1234 if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 1235 1236 table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 1237 nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 1238 1239 for ( i=0; i<is_max; i++ ) { 1240 /* Initialise the two local arrays */ 1241 isz = 0; 1242 PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1243 1244 /* Extract the indices, assume there can be duplicate entries */ 1245 ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1246 ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1247 1248 /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1249 for ( j=0; j<n ; ++j){ 1250 if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 1251 } 1252 ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 1253 ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1254 1255 k = 0; 1256 for ( j=0; j<ov; j++){ /* for each overlap*/ 1257 n = isz; 1258 for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1259 row = nidx[k]; 1260 start = ai[row]; 1261 end = ai[row+1]; 1262 for ( l = start; l<end ; l++){ 1263 val = aj[l] + shift; 1264 if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1265 } 1266 } 1267 } 1268 ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1269 } 1270 PetscFree(table); 1271 PetscFree(nidx); 1272 return 0; 1273 } 1274 1275 int MatPrintHelp_SeqAIJ(Mat A) 1276 { 1277 static int called = 0; 1278 MPI_Comm comm = A->comm; 1279 1280 if (called) return 0; else called = 1; 1281 PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 1282 PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>\n"); 1283 PetscPrintf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 1284 PetscPrintf(comm," -mat_aij_no_inode - Do not use inodes\n"); 1285 PetscPrintf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1286 #if defined(HAVE_ESSL) 1287 PetscPrintf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1288 #endif 1289 return 0; 1290 } 1291 static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1292 /* -------------------------------------------------------------------*/ 1293 static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 1294 MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1295 MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1296 MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 1297 MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 1298 MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 1299 MatLUFactor_SeqAIJ,0, 1300 MatRelax_SeqAIJ, 1301 MatTranspose_SeqAIJ, 1302 MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1303 MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 1304 0,MatAssemblyEnd_SeqAIJ, 1305 MatCompress_SeqAIJ, 1306 MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 1307 MatGetReordering_SeqAIJ, 1308 MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 1309 MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 1310 MatILUFactorSymbolic_SeqAIJ,0, 1311 0,0,MatConvert_SeqAIJ, 1312 MatGetSubMatrix_SeqAIJ,0, 1313 MatConvertSameType_SeqAIJ,0,0, 1314 MatILUFactor_SeqAIJ,0,0, 1315 MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1316 MatGetValues_SeqAIJ,0, 1317 MatPrintHelp_SeqAIJ, 1318 MatScale_SeqAIJ}; 1319 1320 extern int MatUseSuperLU_SeqAIJ(Mat); 1321 extern int MatUseEssl_SeqAIJ(Mat); 1322 extern int MatUseDXML_SeqAIJ(Mat); 1323 1324 /*@C 1325 MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 1326 (the default parallel PETSc format). For good matrix assembly performance 1327 the user should preallocate the matrix storage by setting the parameter nz 1328 (or nzz). By setting these parameters accurately, performance can be 1329 increased by more than a factor of 50. 1330 1331 Input Parameters: 1332 . comm - MPI communicator, set to MPI_COMM_SELF 1333 . m - number of rows 1334 . n - number of columns 1335 . nz - number of nonzeros per row (same for all rows) 1336 . nzz - number of nonzeros per row or null (possibly different for each row) 1337 1338 Output Parameter: 1339 . A - the matrix 1340 1341 Notes: 1342 The AIJ format (also called the Yale sparse matrix format or 1343 compressed row storage), is fully compatible with standard Fortran 77 1344 storage. That is, the stored row and column indices can begin at 1345 either one (as in Fortran) or zero. See the users' manual for details. 1346 1347 Specify the preallocated storage with either nz or nnz (not both). 1348 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1349 allocation. For additional details, see the users' manual chapter on 1350 matrices and the file $(PETSC_DIR)/Performance. 1351 1352 By default, this format uses inodes (identical nodes) when possible, to 1353 improve numerical efficiency of Matrix vector products and solves. We 1354 search for consecutive rows with the same nonzero structure, thereby 1355 reusing matrix information to achieve increased efficiency. 1356 1357 Options Database Keys: 1358 $ -mat_aij_no_inode - Do not use inodes 1359 $ -mat_aij_inode_limit <limit> - Set inode limit. 1360 $ (max limit=5) 1361 $ -mat_aij_oneindex - Internally use indexing starting at 1 1362 $ rather than 0. Note: When calling MatSetValues(), 1363 $ the user still MUST index entries starting at 0! 1364 1365 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 1366 @*/ 1367 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 1368 { 1369 Mat B; 1370 Mat_SeqAIJ *b; 1371 int i, len, ierr, flg; 1372 1373 *A = 0; 1374 PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1375 PLogObjectCreate(B); 1376 B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1377 PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1378 PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1379 B->destroy = MatDestroy_SeqAIJ; 1380 B->view = MatView_SeqAIJ; 1381 B->factor = 0; 1382 B->lupivotthreshold = 1.0; 1383 ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 1384 &flg); CHKERRQ(ierr); 1385 b->ilu_preserve_row_sums = PETSC_FALSE; 1386 ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 1387 (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1388 b->row = 0; 1389 b->col = 0; 1390 b->indexshift = 0; 1391 b->reallocs = 0; 1392 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 1393 if (flg) b->indexshift = -1; 1394 1395 b->m = m; B->m = m; B->M = m; 1396 b->n = n; B->n = n; B->N = n; 1397 b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1398 if (nnz == PETSC_NULL) { 1399 if (nz == PETSC_DEFAULT) nz = 10; 1400 else if (nz <= 0) nz = 1; 1401 for ( i=0; i<m; i++ ) b->imax[i] = nz; 1402 nz = nz*m; 1403 } 1404 else { 1405 nz = 0; 1406 for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 1407 } 1408 1409 /* allocate the matrix space */ 1410 len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 1411 b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1412 b->j = (int *) (b->a + nz); 1413 PetscMemzero(b->j,nz*sizeof(int)); 1414 b->i = b->j + nz; 1415 b->singlemalloc = 1; 1416 1417 b->i[0] = -b->indexshift; 1418 for (i=1; i<m+1; i++) { 1419 b->i[i] = b->i[i-1] + b->imax[i-1]; 1420 } 1421 1422 /* b->ilen will count nonzeros in each row so far. */ 1423 b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1424 PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1425 for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 1426 1427 b->nz = 0; 1428 b->maxnz = nz; 1429 b->sorted = 0; 1430 b->roworiented = 1; 1431 b->nonew = 0; 1432 b->diag = 0; 1433 b->solve_work = 0; 1434 b->spptr = 0; 1435 b->inode.node_count = 0; 1436 b->inode.size = 0; 1437 b->inode.limit = 5; 1438 b->inode.max_limit = 5; 1439 1440 *A = B; 1441 /* SuperLU is not currently supported through PETSc */ 1442 #if defined(HAVE_SUPERLU) 1443 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 1444 if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 1445 #endif 1446 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 1447 if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 1448 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 1449 if (flg) { 1450 if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1451 ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 1452 } 1453 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1454 if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 1455 return 0; 1456 } 1457 1458 int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 1459 { 1460 Mat C; 1461 Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 1462 int i,len, m = a->m,shift = a->indexshift; 1463 1464 *B = 0; 1465 PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1466 PLogObjectCreate(C); 1467 C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 1468 PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1469 C->destroy = MatDestroy_SeqAIJ; 1470 C->view = MatView_SeqAIJ; 1471 C->factor = A->factor; 1472 c->row = 0; 1473 c->col = 0; 1474 c->indexshift = shift; 1475 C->assembled = PETSC_TRUE; 1476 1477 c->m = C->m = a->m; 1478 c->n = C->n = a->n; 1479 C->M = a->m; 1480 C->N = a->n; 1481 1482 c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 1483 c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 1484 for ( i=0; i<m; i++ ) { 1485 c->imax[i] = a->imax[i]; 1486 c->ilen[i] = a->ilen[i]; 1487 } 1488 1489 /* allocate the matrix space */ 1490 c->singlemalloc = 1; 1491 len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 1492 c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1493 c->j = (int *) (c->a + a->i[m] + shift); 1494 c->i = c->j + a->i[m] + shift; 1495 PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 1496 if (m > 0) { 1497 PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 1498 if (cpvalues == COPY_VALUES) { 1499 PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 1500 } 1501 } 1502 1503 PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1504 c->sorted = a->sorted; 1505 c->roworiented = a->roworiented; 1506 c->nonew = a->nonew; 1507 c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 1508 1509 if (a->diag) { 1510 c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1511 PLogObjectMemory(C,(m+1)*sizeof(int)); 1512 for ( i=0; i<m; i++ ) { 1513 c->diag[i] = a->diag[i]; 1514 } 1515 } 1516 else c->diag = 0; 1517 c->inode.limit = a->inode.limit; 1518 c->inode.max_limit = a->inode.max_limit; 1519 if (a->inode.size){ 1520 c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1521 c->inode.node_count = a->inode.node_count; 1522 PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1523 } else { 1524 c->inode.size = 0; 1525 c->inode.node_count = 0; 1526 } 1527 c->nz = a->nz; 1528 c->maxnz = a->maxnz; 1529 c->solve_work = 0; 1530 c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1531 1532 *B = C; 1533 return 0; 1534 } 1535 1536 int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 1537 { 1538 Mat_SeqAIJ *a; 1539 Mat B; 1540 int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1541 MPI_Comm comm; 1542 1543 PetscObjectGetComm((PetscObject) viewer,&comm); 1544 MPI_Comm_size(comm,&size); 1545 if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 1546 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 1547 ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 1548 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 1549 M = header[1]; N = header[2]; nz = header[3]; 1550 1551 /* read in row lengths */ 1552 rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1553 ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 1554 1555 /* create our matrix */ 1556 ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1557 B = *A; 1558 a = (Mat_SeqAIJ *) B->data; 1559 shift = a->indexshift; 1560 1561 /* read in column indices and adjust for Fortran indexing*/ 1562 ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 1563 if (shift) { 1564 for ( i=0; i<nz; i++ ) { 1565 a->j[i] += 1; 1566 } 1567 } 1568 1569 /* read in nonzero values */ 1570 ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 1571 1572 /* set matrix "i" values */ 1573 a->i[0] = -shift; 1574 for ( i=1; i<= M; i++ ) { 1575 a->i[i] = a->i[i-1] + rowlengths[i-1]; 1576 a->ilen[i-1] = rowlengths[i-1]; 1577 } 1578 PetscFree(rowlengths); 1579 1580 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1581 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1582 return 0; 1583 } 1584 1585 static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 1586 { 1587 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 1588 1589 if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type"); 1590 1591 /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 1592 if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1593 (a->indexshift != b->indexshift)) { 1594 *flg = PETSC_FALSE; return 0; 1595 } 1596 1597 /* if the a->i are the same */ 1598 if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 1599 *flg = PETSC_FALSE; return 0; 1600 } 1601 1602 /* if a->j are the same */ 1603 if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 1604 *flg = PETSC_FALSE; return 0; 1605 } 1606 1607 /* if a->a are the same */ 1608 if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 1609 *flg = PETSC_FALSE; return 0; 1610 } 1611 *flg = PETSC_TRUE; 1612 return 0; 1613 1614 } 1615