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