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