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