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