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