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