1 #ifndef lint 2 static char vcid[] = "$Id: aij.c,v 1.121 1995/11/25 19:07:38 curfman 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 497 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 498 return 0; 499 } 500 501 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 502 { 503 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 504 int i,j, n,shift = a->indexshift; 505 Scalar *x, zero = 0.0; 506 507 if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix"); 508 VecSet(&zero,v); 509 VecGetArray(v,&x); VecGetLocalSize(v,&n); 510 if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 511 for ( i=0; i<a->m; i++ ) { 512 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 513 if (a->j[j]+shift == i) { 514 x[i] = a->a[j]; 515 break; 516 } 517 } 518 } 519 return 0; 520 } 521 522 /* -------------------------------------------------------*/ 523 /* Should check that shapes of vectors and matrices match */ 524 /* -------------------------------------------------------*/ 525 static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 526 { 527 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 528 Scalar *x, *y, *v, alpha; 529 int m = a->m, n, i, *idx, shift = a->indexshift; 530 531 if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix"); 532 VecGetArray(xx,&x); VecGetArray(yy,&y); 533 PetscMemzero(y,a->n*sizeof(Scalar)); 534 y = y + shift; /* shift for Fortran start by 1 indexing */ 535 for ( i=0; i<m; i++ ) { 536 idx = a->j + a->i[i] + shift; 537 v = a->a + a->i[i] + shift; 538 n = a->i[i+1] - a->i[i]; 539 alpha = x[i]; 540 while (n-->0) {y[*idx++] += alpha * *v++;} 541 } 542 PLogFlops(2*a->nz - a->n); 543 return 0; 544 } 545 546 static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 547 { 548 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 549 Scalar *x, *y, *v, alpha; 550 int m = a->m, n, i, *idx,shift = a->indexshift; 551 552 if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix"); 553 VecGetArray(xx,&x); VecGetArray(yy,&y); 554 if (zz != yy) VecCopy(zz,yy); 555 y = y + shift; /* shift for Fortran start by 1 indexing */ 556 for ( i=0; i<m; i++ ) { 557 idx = a->j + a->i[i] + shift; 558 v = a->a + a->i[i] + shift; 559 n = a->i[i+1] - a->i[i]; 560 alpha = x[i]; 561 while (n-->0) {y[*idx++] += alpha * *v++;} 562 } 563 return 0; 564 } 565 566 static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 567 { 568 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 569 Scalar *x, *y, *v, sum; 570 int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 571 572 if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix"); 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 if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix"); 597 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 598 x = x + shift; /* shift for Fortran start by 1 indexing */ 599 idx = a->j; 600 v = a->a; 601 ii = a->i; 602 for ( i=0; i<m; i++ ) { 603 n = ii[1] - ii[0]; ii++; 604 sum = y[i]; 605 /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 606 while (n--) sum += *v++ * x[*idx++]; 607 z[i] = sum; 608 } 609 PLogFlops(2*a->nz); 610 return 0; 611 } 612 613 /* 614 Adds diagonal pointers to sparse matrix structure. 615 */ 616 617 int MatMarkDiag_SeqAIJ(Mat A) 618 { 619 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 620 int i,j, *diag, m = a->m,shift = a->indexshift; 621 622 if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix"); 623 diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 624 PLogObjectMemory(A,(m+1)*sizeof(int)); 625 for ( i=0; i<a->m; i++ ) { 626 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 627 if (a->j[j]+shift == i) { 628 diag[i] = j - shift; 629 break; 630 } 631 } 632 } 633 a->diag = diag; 634 return 0; 635 } 636 637 static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 638 double fshift,int its,Vec xx) 639 { 640 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 641 Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 642 int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 643 644 VecGetArray(xx,&x); VecGetArray(bb,&b); 645 if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 646 diag = a->diag; 647 xs = x + shift; /* shifted by one for index start of a or a->j*/ 648 if (flag == SOR_APPLY_UPPER) { 649 /* apply ( U + D/omega) to the vector */ 650 bs = b + shift; 651 for ( i=0; i<m; i++ ) { 652 d = fshift + a->a[diag[i] + shift]; 653 n = a->i[i+1] - diag[i] - 1; 654 idx = a->j + diag[i] + (!shift); 655 v = a->a + diag[i] + (!shift); 656 sum = b[i]*d/omega; 657 SPARSEDENSEDOT(sum,bs,v,idx,n); 658 x[i] = sum; 659 } 660 return 0; 661 } 662 if (flag == SOR_APPLY_LOWER) { 663 SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 664 } 665 else if (flag & SOR_EISENSTAT) { 666 /* Let A = L + U + D; where L is lower trianglar, 667 U is upper triangular, E is diagonal; This routine applies 668 669 (L + E)^{-1} A (U + E)^{-1} 670 671 to a vector efficiently using Eisenstat's trick. This is for 672 the case of SSOR preconditioner, so E is D/omega where omega 673 is the relaxation factor. 674 */ 675 t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 676 scale = (2.0/omega) - 1.0; 677 678 /* x = (E + U)^{-1} b */ 679 for ( i=m-1; i>=0; i-- ) { 680 d = fshift + a->a[diag[i] + shift]; 681 n = a->i[i+1] - diag[i] - 1; 682 idx = a->j + diag[i] + (!shift); 683 v = a->a + diag[i] + (!shift); 684 sum = b[i]; 685 SPARSEDENSEMDOT(sum,xs,v,idx,n); 686 x[i] = omega*(sum/d); 687 } 688 689 /* t = b - (2*E - D)x */ 690 v = a->a; 691 for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 692 693 /* t = (E + L)^{-1}t */ 694 ts = t + shift; /* shifted by one for index start of a or a->j*/ 695 diag = a->diag; 696 for ( i=0; i<m; i++ ) { 697 d = fshift + a->a[diag[i]+shift]; 698 n = diag[i] - a->i[i]; 699 idx = a->j + a->i[i] + shift; 700 v = a->a + a->i[i] + shift; 701 sum = t[i]; 702 SPARSEDENSEMDOT(sum,ts,v,idx,n); 703 t[i] = omega*(sum/d); 704 } 705 706 /* x = x + t */ 707 for ( i=0; i<m; i++ ) { x[i] += t[i]; } 708 PetscFree(t); 709 return 0; 710 } 711 if (flag & SOR_ZERO_INITIAL_GUESS) { 712 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 713 for ( i=0; i<m; i++ ) { 714 d = fshift + a->a[diag[i]+shift]; 715 n = diag[i] - a->i[i]; 716 idx = a->j + a->i[i] + shift; 717 v = a->a + a->i[i] + shift; 718 sum = b[i]; 719 SPARSEDENSEMDOT(sum,xs,v,idx,n); 720 x[i] = omega*(sum/d); 721 } 722 xb = x; 723 } 724 else xb = b; 725 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 726 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 727 for ( i=0; i<m; i++ ) { 728 x[i] *= a->a[diag[i]+shift]; 729 } 730 } 731 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 732 for ( i=m-1; i>=0; i-- ) { 733 d = fshift + a->a[diag[i] + shift]; 734 n = a->i[i+1] - diag[i] - 1; 735 idx = a->j + diag[i] + (!shift); 736 v = a->a + diag[i] + (!shift); 737 sum = xb[i]; 738 SPARSEDENSEMDOT(sum,xs,v,idx,n); 739 x[i] = omega*(sum/d); 740 } 741 } 742 its--; 743 } 744 while (its--) { 745 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 746 for ( i=0; i<m; i++ ) { 747 d = fshift + a->a[diag[i]+shift]; 748 n = a->i[i+1] - a->i[i]; 749 idx = a->j + a->i[i] + shift; 750 v = a->a + a->i[i] + shift; 751 sum = b[i]; 752 SPARSEDENSEMDOT(sum,xs,v,idx,n); 753 x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 754 } 755 } 756 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 757 for ( i=m-1; i>=0; i-- ) { 758 d = fshift + a->a[diag[i] + shift]; 759 n = a->i[i+1] - a->i[i]; 760 idx = a->j + a->i[i] + shift; 761 v = a->a + a->i[i] + shift; 762 sum = b[i]; 763 SPARSEDENSEMDOT(sum,xs,v,idx,n); 764 x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 765 } 766 } 767 } 768 return 0; 769 } 770 771 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 772 { 773 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 774 *nz = a->nz; 775 *nzalloc = a->maxnz; 776 *mem = (int)A->mem; 777 return 0; 778 } 779 780 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 781 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 782 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 783 extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 784 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 785 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 786 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 787 788 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 789 { 790 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 791 int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 792 793 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 794 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 795 if (diag) { 796 for ( i=0; i<N; i++ ) { 797 if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 798 if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 799 a->ilen[rows[i]] = 1; 800 a->a[a->i[rows[i]]+shift] = *diag; 801 a->j[a->i[rows[i]]+shift] = rows[i]+shift; 802 } 803 else { 804 ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 805 CHKERRQ(ierr); 806 } 807 } 808 } 809 else { 810 for ( i=0; i<N; i++ ) { 811 if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 812 a->ilen[rows[i]] = 0; 813 } 814 } 815 ISRestoreIndices(is,&rows); 816 ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 817 ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 818 return 0; 819 } 820 821 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 822 { 823 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 824 *m = a->m; *n = a->n; 825 return 0; 826 } 827 828 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 829 { 830 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 831 *m = 0; *n = a->m; 832 return 0; 833 } 834 static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 835 { 836 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 837 int *itmp,i,ierr,shift = a->indexshift; 838 839 if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 840 841 if (!a->assembled) { 842 ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 843 ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 844 } 845 *nz = a->i[row+1] - a->i[row]; 846 if (v) *v = a->a + a->i[row] + shift; 847 if (idx) { 848 if (*nz) { 849 itmp = a->j + a->i[row] + shift; 850 *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 851 for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 852 } 853 else *idx = 0; 854 } 855 return 0; 856 } 857 858 static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 859 { 860 if (idx) {if (*idx) PetscFree(*idx);} 861 return 0; 862 } 863 864 static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 865 { 866 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 867 Scalar *v = a->a; 868 double sum = 0.0; 869 int i, j,shift = a->indexshift; 870 871 if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix"); 872 if (type == NORM_FROBENIUS) { 873 for (i=0; i<a->nz; i++ ) { 874 #if defined(PETSC_COMPLEX) 875 sum += real(conj(*v)*(*v)); v++; 876 #else 877 sum += (*v)*(*v); v++; 878 #endif 879 } 880 *norm = sqrt(sum); 881 } 882 else if (type == NORM_1) { 883 double *tmp; 884 int *jj = a->j; 885 tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 886 PetscMemzero(tmp,a->n*sizeof(double)); 887 *norm = 0.0; 888 for ( j=0; j<a->nz; j++ ) { 889 tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 890 } 891 for ( j=0; j<a->n; j++ ) { 892 if (tmp[j] > *norm) *norm = tmp[j]; 893 } 894 PetscFree(tmp); 895 } 896 else if (type == NORM_INFINITY) { 897 *norm = 0.0; 898 for ( j=0; j<a->m; j++ ) { 899 v = a->a + a->i[j] + shift; 900 sum = 0.0; 901 for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 902 sum += PetscAbsScalar(*v); v++; 903 } 904 if (sum > *norm) *norm = sum; 905 } 906 } 907 else { 908 SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 909 } 910 return 0; 911 } 912 913 static int MatTranspose_SeqAIJ(Mat A,Mat *B) 914 { 915 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 916 Mat C; 917 int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 918 int shift = a->indexshift; 919 Scalar *array = a->a; 920 921 if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 922 col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 923 PetscMemzero(col,(1+a->n)*sizeof(int)); 924 if (shift) { 925 for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 926 } 927 for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 928 ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 929 PetscFree(col); 930 for ( i=0; i<m; i++ ) { 931 len = ai[i+1]-ai[i]; 932 ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 933 array += len; aj += len; 934 } 935 if (shift) { 936 for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 937 } 938 939 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 940 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 941 942 if (B) { 943 *B = C; 944 } else { 945 /* This isn't really an in-place transpose */ 946 PetscFree(a->a); 947 if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 948 if (a->diag) PetscFree(a->diag); 949 if (a->ilen) PetscFree(a->ilen); 950 if (a->imax) PetscFree(a->imax); 951 if (a->solve_work) PetscFree(a->solve_work); 952 PetscFree(a); 953 PetscMemcpy(A,C,sizeof(struct _Mat)); 954 PetscHeaderDestroy(C); 955 } 956 return 0; 957 } 958 959 static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr) 960 { 961 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 962 Scalar *l,*r,x,*v; 963 int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 964 965 if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix"); 966 if (ll) { 967 VecGetArray(ll,&l); VecGetSize(ll,&m); 968 if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length"); 969 v = a->a; 970 for ( i=0; i<m; i++ ) { 971 x = l[i]; 972 M = a->i[i+1] - a->i[i]; 973 for ( j=0; j<M; j++ ) { (*v++) *= x;} 974 } 975 } 976 if (rr) { 977 VecGetArray(rr,&r); VecGetSize(rr,&n); 978 if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length"); 979 v = a->a; jj = a->j; 980 for ( i=0; i<nz; i++ ) { 981 (*v++) *= r[*jj++ + shift]; 982 } 983 } 984 return 0; 985 } 986 987 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 988 { 989 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 990 int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 991 register int sum,lensi; 992 int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 993 int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 994 Scalar *vwork,*a_new; 995 Mat C; 996 997 if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix"); 998 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 999 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 1000 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 1001 1002 if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 1003 /* special case of contiguous rows */ 1004 lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 1005 starts = lens + ncols; 1006 /* loop over new rows determining lens and starting points */ 1007 for (i=0; i<nrows; i++) { 1008 kstart = ai[irow[i]]+shift; 1009 kend = kstart + ailen[irow[i]]; 1010 for ( k=kstart; k<kend; k++ ) { 1011 if (aj[k] >= first) { 1012 starts[i] = k; 1013 break; 1014 } 1015 } 1016 sum = 0; 1017 while (k < kend) { 1018 if (aj[k++] >= first+ncols) break; 1019 sum++; 1020 } 1021 lens[i] = sum; 1022 } 1023 /* create submatrix */ 1024 if (scall == MAT_REUSE_MATRIX) { 1025 int n_cols,n_rows; 1026 ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1027 if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1028 MatZeroEntries(*B); 1029 C = *B; 1030 } 1031 else { 1032 ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 1033 } 1034 c = (Mat_SeqAIJ*) C->data; 1035 1036 /* loop over rows inserting into submatrix */ 1037 a_new = c->a; 1038 j_new = c->j; 1039 i_new = c->i; 1040 i_new[0] = -shift; 1041 for (i=0; i<nrows; i++) { 1042 ii = starts[i]; 1043 lensi = lens[i]; 1044 for ( k=0; k<lensi; k++ ) { 1045 *j_new++ = aj[ii+k] - first; 1046 } 1047 PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1048 a_new += lensi; 1049 i_new[i+1] = i_new[i] + lensi; 1050 c->ilen[i] = lensi; 1051 } 1052 PetscFree(lens); 1053 } 1054 else { 1055 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 1056 smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 1057 ssmap = smap + shift; 1058 cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 1059 lens = cwork + ncols; 1060 vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1061 PetscMemzero(smap,oldcols*sizeof(int)); 1062 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 1063 /* determine lens of each row */ 1064 for (i=0; i<nrows; i++) { 1065 kstart = a->i[irow[i]]+shift; 1066 kend = kstart + a->ilen[irow[i]]; 1067 lens[i] = 0; 1068 for ( k=kstart; k<kend; k++ ) { 1069 if (ssmap[a->j[k]]) { 1070 lens[i]++; 1071 } 1072 } 1073 } 1074 /* Create and fill new matrix */ 1075 if (scall == MAT_REUSE_MATRIX) { 1076 int n_cols,n_rows; 1077 ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1078 if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1079 MatZeroEntries(*B); 1080 C = *B; 1081 } 1082 else { 1083 ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 1084 } 1085 for (i=0; i<nrows; i++) { 1086 nznew = 0; 1087 kstart = a->i[irow[i]]+shift; 1088 kend = kstart + a->ilen[irow[i]]; 1089 for ( k=kstart; k<kend; k++ ) { 1090 if (ssmap[a->j[k]]) { 1091 cwork[nznew] = ssmap[a->j[k]] - 1; 1092 vwork[nznew++] = a->a[k]; 1093 } 1094 } 1095 ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 1096 } 1097 /* Free work space */ 1098 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 1099 PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 1100 } 1101 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1102 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1103 1104 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1105 *B = C; 1106 return 0; 1107 } 1108 1109 /* 1110 note: This can only work for identity for row and col. It would 1111 be good to check this and otherwise generate an error. 1112 */ 1113 static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1114 { 1115 Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1116 int ierr; 1117 Mat outA; 1118 1119 if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1120 1121 outA = inA; 1122 inA->factor = FACTOR_LU; 1123 a->row = row; 1124 a->col = col; 1125 1126 a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 1127 1128 if (!a->diag) { 1129 ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 1130 } 1131 ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1132 return 0; 1133 } 1134 1135 static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1136 Mat **B) 1137 { 1138 int ierr,i; 1139 1140 if (scall == MAT_INITIAL_MATRIX) { 1141 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1142 } 1143 1144 for ( i=0; i<n; i++ ) { 1145 ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1146 } 1147 return 0; 1148 } 1149 1150 static int MatIncreaseOverlap_SeqAIJ(Mat A, int n, IS *is, int ov) 1151 { 1152 int i,m,*idx,ierr; 1153 1154 for ( i=0; i<n; i++ ) { 1155 ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1156 ISGetLocalSize(is[i],&m); 1157 } 1158 SETERRQ(1,"MatIncreaseOverlap_SeqAIJ:Not implemented"); 1159 } 1160 /* -------------------------------------------------------------------*/ 1161 1162 static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 1163 MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1164 MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1165 MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 1166 MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 1167 MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 1168 MatLUFactor_SeqAIJ,0, 1169 MatRelax_SeqAIJ, 1170 MatTranspose_SeqAIJ, 1171 MatGetInfo_SeqAIJ,0, 1172 MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ, 1173 0,MatAssemblyEnd_SeqAIJ, 1174 MatCompress_SeqAIJ, 1175 MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 1176 MatGetReordering_SeqAIJ, 1177 MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 1178 MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 1179 MatILUFactorSymbolic_SeqAIJ,0, 1180 0,0,MatConvert_SeqAIJ, 1181 MatGetSubMatrix_SeqAIJ,0, 1182 MatCopyPrivate_SeqAIJ,0,0, 1183 MatILUFactor_SeqAIJ,0,0, 1184 MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1185 MatGetValues_SeqAIJ}; 1186 1187 extern int MatUseSuperLU_SeqAIJ(Mat); 1188 extern int MatUseEssl_SeqAIJ(Mat); 1189 extern int MatUseDXML_SeqAIJ(Mat); 1190 1191 /*@C 1192 MatCreateSeqAIJ - Creates a sparse matrix in AIJ format 1193 (the default uniprocessor PETSc format). 1194 1195 Input Parameters: 1196 . comm - MPI communicator, set to MPI_COMM_SELF 1197 . m - number of rows 1198 . n - number of columns 1199 . nz - number of nonzeros per row (same for all rows) 1200 . nzz - number of nonzeros per row or null (possibly different for each row) 1201 1202 Output Parameter: 1203 . A - the matrix 1204 1205 Notes: 1206 The AIJ format (also called the Yale sparse matrix format or 1207 compressed row storage), is fully compatible with standard Fortran 77 1208 storage. That is, the stored row and column indices can begin at 1209 either one (as in Fortran) or zero. See the users manual for details. 1210 1211 Specify the preallocated storage with either nz or nnz (not both). 1212 Set both nz and nnz to zero for PETSc to control dynamic memory 1213 allocation. 1214 1215 .keywords: matrix, aij, compressed row, sparse 1216 1217 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 1218 @*/ 1219 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 1220 { 1221 Mat B; 1222 Mat_SeqAIJ *b; 1223 int i,len,ierr; 1224 1225 *A = 0; 1226 PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1227 PLogObjectCreate(B); 1228 B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1229 PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1230 B->destroy = MatDestroy_SeqAIJ; 1231 B->view = MatView_SeqAIJ; 1232 B->factor = 0; 1233 B->lupivotthreshold = 1.0; 1234 OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold); 1235 b->row = 0; 1236 b->col = 0; 1237 b->indexshift = 0; 1238 if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1; 1239 1240 b->m = m; 1241 b->n = n; 1242 b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1243 if (!nnz) { 1244 if (nz <= 0) nz = 1; 1245 for ( i=0; i<m; i++ ) b->imax[i] = nz; 1246 nz = nz*m; 1247 } 1248 else { 1249 nz = 0; 1250 for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 1251 } 1252 1253 /* allocate the matrix space */ 1254 len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 1255 b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1256 b->j = (int *) (b->a + nz); 1257 PetscMemzero(b->j,nz*sizeof(int)); 1258 b->i = b->j + nz; 1259 b->singlemalloc = 1; 1260 1261 b->i[0] = -b->indexshift; 1262 for (i=1; i<m+1; i++) { 1263 b->i[i] = b->i[i-1] + b->imax[i-1]; 1264 } 1265 1266 /* b->ilen will count nonzeros in each row so far. */ 1267 b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1268 PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1269 for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 1270 1271 b->nz = 0; 1272 b->maxnz = nz; 1273 b->sorted = 0; 1274 b->roworiented = 1; 1275 b->nonew = 0; 1276 b->diag = 0; 1277 b->assembled = 0; 1278 b->solve_work = 0; 1279 b->spptr = 0; 1280 b->inode.node_count = 0; 1281 b->inode.size = 0; 1282 1283 *A = B; 1284 if (OptionsHasName(0,"-mat_aij_superlu")) { 1285 ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); 1286 } 1287 if (OptionsHasName(0,"-mat_aij_essl")) { 1288 ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); 1289 } 1290 if (OptionsHasName(0,"-mat_aij_dxml")) { 1291 if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1292 ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 1293 } 1294 1295 return 0; 1296 } 1297 1298 int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues) 1299 { 1300 Mat C; 1301 Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 1302 int i,len, m = a->m,shift = a->indexshift; 1303 1304 *B = 0; 1305 if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix"); 1306 PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1307 PLogObjectCreate(C); 1308 C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 1309 PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1310 C->destroy = MatDestroy_SeqAIJ; 1311 C->view = MatView_SeqAIJ; 1312 C->factor = A->factor; 1313 c->row = 0; 1314 c->col = 0; 1315 c->indexshift = shift; 1316 1317 c->m = a->m; 1318 c->n = a->n; 1319 1320 c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 1321 c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 1322 for ( i=0; i<m; i++ ) { 1323 c->imax[i] = a->imax[i]; 1324 c->ilen[i] = a->ilen[i]; 1325 } 1326 1327 /* allocate the matrix space */ 1328 c->singlemalloc = 1; 1329 len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 1330 c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1331 c->j = (int *) (c->a + a->i[m] + shift); 1332 c->i = c->j + a->i[m] + shift; 1333 PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 1334 if (m > 0) { 1335 PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 1336 if (cpvalues == COPY_VALUES) { 1337 PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 1338 } 1339 } 1340 1341 PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1342 c->sorted = a->sorted; 1343 c->roworiented = a->roworiented; 1344 c->nonew = a->nonew; 1345 1346 if (a->diag) { 1347 c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1348 PLogObjectMemory(C,(m+1)*sizeof(int)); 1349 for ( i=0; i<m; i++ ) { 1350 c->diag[i] = a->diag[i]; 1351 } 1352 } 1353 else c->diag = 0; 1354 if( a->inode.size){ 1355 c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1356 c->inode.node_count = a->inode.node_count; 1357 PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1358 } else { 1359 c->inode.size = 0; 1360 c->inode.node_count = 0; 1361 } 1362 c->assembled = 1; 1363 c->nz = a->nz; 1364 c->maxnz = a->maxnz; 1365 c->solve_work = 0; 1366 c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1367 1368 *B = C; 1369 return 0; 1370 } 1371 1372 int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 1373 { 1374 Mat_SeqAIJ *a; 1375 Mat B; 1376 int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1377 PetscObject vobj = (PetscObject) bview; 1378 MPI_Comm comm = vobj->comm; 1379 1380 MPI_Comm_size(comm,&size); 1381 if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 1382 ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1383 ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 1384 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 1385 M = header[1]; N = header[2]; nz = header[3]; 1386 1387 /* read in row lengths */ 1388 rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1389 ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 1390 1391 /* create our matrix */ 1392 ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1393 B = *A; 1394 a = (Mat_SeqAIJ *) B->data; 1395 shift = a->indexshift; 1396 1397 /* read in column indices and adjust for Fortran indexing*/ 1398 ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 1399 if (shift) { 1400 for ( i=0; i<nz; i++ ) { 1401 a->j[i] += 1; 1402 } 1403 } 1404 1405 /* read in nonzero values */ 1406 ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 1407 1408 /* set matrix "i" values */ 1409 a->i[0] = -shift; 1410 for ( i=1; i<= M; i++ ) { 1411 a->i[i] = a->i[i-1] + rowlengths[i-1]; 1412 a->ilen[i-1] = rowlengths[i-1]; 1413 } 1414 PetscFree(rowlengths); 1415 1416 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1417 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1418 return 0; 1419 } 1420 1421 1422 1423