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