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