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