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