117ab2063SBarry Smith #ifndef lint 2*f31bba2fSSatish Balay static char vcid[] = "$Id: aij.c,v 1.117 1995/11/16 16:55:50 balay Exp balay $"; 317ab2063SBarry Smith #endif 417ab2063SBarry Smith 5d5d45c9bSBarry Smith /* 6d5d45c9bSBarry Smith Defines the basic matrix operations for the AIJ (compressed row) 7d5d45c9bSBarry Smith matrix storage format. 8d5d45c9bSBarry Smith */ 917ab2063SBarry Smith #include "aij.h" 1017ab2063SBarry Smith #include "vec/vecimpl.h" 1117ab2063SBarry Smith #include "inline/spops.h" 1217ab2063SBarry Smith 1317ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**); 1417ab2063SBarry Smith 15416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1617ab2063SBarry Smith { 17416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 18a2744918SBarry Smith int ierr, *ia, *ja,n,*idx,i; 19*f31bba2fSSatish Balay Viewer V1, V2; 2017ab2063SBarry Smith 21416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix"); 2217ab2063SBarry Smith 23a2744918SBarry Smith /* 24a2744918SBarry Smith this is tacky: In the future when we have written special factorization 25a2744918SBarry Smith and solve routines for the identity permutation we should use a 26a2744918SBarry Smith stride index set instead of the general one. 27a2744918SBarry Smith */ 28a2744918SBarry Smith if (type == ORDER_NATURAL) { 29a2744918SBarry Smith n = a->n; 30a2744918SBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 31a2744918SBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 32a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 33a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 34a2744918SBarry Smith PetscFree(idx); 35a2744918SBarry Smith ISSetPermutation(*rperm); 36a2744918SBarry Smith ISSetPermutation(*cperm); 37a2744918SBarry Smith ISSetIdentity(*rperm); 38a2744918SBarry Smith ISSetIdentity(*cperm); 39a2744918SBarry Smith return 0; 40a2744918SBarry Smith } 41a2744918SBarry Smith 42416022c9SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr); 43416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 4444f6d32bSSatish Balay /* ISView(*rperm, STDOUT_VIEWER_SELF);*/ 45*f31bba2fSSatish Balay 46*f31bba2fSSatish Balay ViewerFileOpenASCII(MPI_COMM_SELF,"row_is_orig", &V1); 47*f31bba2fSSatish Balay ViewerFileOpenASCII(MPI_COMM_SELF,"col_is_orig", &V2); 48*f31bba2fSSatish Balay ISView(*rperm,V1); 49*f31bba2fSSatish Balay ISView(*cperm,V2); 50*f31bba2fSSatish Balay ViewerDestroy(V1); 51*f31bba2fSSatish Balay ViewerDestroy(V2); 52*f31bba2fSSatish Balay 530452661fSBarry Smith PetscFree(ia); PetscFree(ja); 5417ab2063SBarry Smith return 0; 5517ab2063SBarry Smith } 5617ab2063SBarry Smith 57416022c9SBarry Smith #define CHUNKSIZE 10 5817ab2063SBarry Smith 5917ab2063SBarry Smith /* This version has row oriented v */ 60416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 6117ab2063SBarry Smith { 62416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 63416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 64416022c9SBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen; 65d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 66416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 6717ab2063SBarry Smith 6817ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 69416022c9SBarry Smith row = im[k]; 7017ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 71416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 7217ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 7317ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 74416022c9SBarry Smith low = 0; 7517ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 76416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 77416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 78416022c9SBarry Smith col = in[l] - shift; value = *v++; 79416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 80416022c9SBarry Smith while (high-low > 5) { 81416022c9SBarry Smith t = (low+high)/2; 82416022c9SBarry Smith if (rp[t] > col) high = t; 83416022c9SBarry Smith else low = t; 8417ab2063SBarry Smith } 85416022c9SBarry Smith for ( i=low; i<high; i++ ) { 8617ab2063SBarry Smith if (rp[i] > col) break; 8717ab2063SBarry Smith if (rp[i] == col) { 88416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 8917ab2063SBarry Smith else ap[i] = value; 9017ab2063SBarry Smith goto noinsert; 9117ab2063SBarry Smith } 9217ab2063SBarry Smith } 9317ab2063SBarry Smith if (nonew) goto noinsert; 9417ab2063SBarry Smith if (nrow >= rmax) { 9517ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 96416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 9717ab2063SBarry Smith Scalar *new_a; 9817ab2063SBarry Smith 9917ab2063SBarry Smith /* malloc new storage space */ 100416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1010452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 10217ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 10317ab2063SBarry Smith new_i = new_j + new_nz; 10417ab2063SBarry Smith 10517ab2063SBarry Smith /* copy over old data into new slots */ 10617ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 107416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 108416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 109416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 110416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 11117ab2063SBarry Smith len*sizeof(int)); 112416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 113416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 11417ab2063SBarry Smith len*sizeof(Scalar)); 11517ab2063SBarry Smith /* free up old matrix storage */ 1160452661fSBarry Smith PetscFree(a->a); 1170452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 118416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 119416022c9SBarry Smith a->singlemalloc = 1; 12017ab2063SBarry Smith 12117ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 122416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 123416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 124416022c9SBarry Smith a->maxnz += CHUNKSIZE; 12517ab2063SBarry Smith } 126416022c9SBarry Smith N = nrow++ - 1; a->nz++; 127416022c9SBarry Smith /* shift up all the later entries in this row */ 128416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 12917ab2063SBarry Smith rp[ii+1] = rp[ii]; 13017ab2063SBarry Smith ap[ii+1] = ap[ii]; 13117ab2063SBarry Smith } 13217ab2063SBarry Smith rp[i] = col; 13317ab2063SBarry Smith ap[i] = value; 13417ab2063SBarry Smith noinsert:; 135416022c9SBarry Smith low = i + 1; 13617ab2063SBarry Smith } 13717ab2063SBarry Smith ailen[row] = nrow; 13817ab2063SBarry Smith } 13917ab2063SBarry Smith return 0; 14017ab2063SBarry Smith } 14117ab2063SBarry Smith 14217ab2063SBarry Smith #include "draw.h" 14317ab2063SBarry Smith #include "pinclude/pviewer.h" 144416022c9SBarry Smith #include "sysio.h" 14517ab2063SBarry Smith 146416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 14717ab2063SBarry Smith { 148416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 149416022c9SBarry Smith int i, fd, *col_lens, ierr; 15017ab2063SBarry Smith 151416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1520452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 153416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 154416022c9SBarry Smith col_lens[1] = a->m; 155416022c9SBarry Smith col_lens[2] = a->n; 156416022c9SBarry Smith col_lens[3] = a->nz; 157416022c9SBarry Smith 158416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 159416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 160416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 16117ab2063SBarry Smith } 162416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 1630452661fSBarry Smith PetscFree(col_lens); 164416022c9SBarry Smith 165416022c9SBarry Smith /* store column indices (zero start index) */ 166416022c9SBarry Smith if (a->indexshift) { 167416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 16817ab2063SBarry Smith } 169416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 170416022c9SBarry Smith if (a->indexshift) { 171416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 17217ab2063SBarry Smith } 173416022c9SBarry Smith 174416022c9SBarry Smith /* store nonzero values */ 175416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 17617ab2063SBarry Smith return 0; 17717ab2063SBarry Smith } 178416022c9SBarry Smith 179416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 180416022c9SBarry Smith { 181416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 182416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 18317ab2063SBarry Smith FILE *fd; 18417ab2063SBarry Smith char *outputname; 18517ab2063SBarry Smith 186416022c9SBarry Smith ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 187416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 188416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 18917ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 19008480c60SBarry Smith /* no need to print additional information */ ; 19117ab2063SBarry Smith } 19217ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 19317ab2063SBarry Smith int nz, nzalloc, mem; 194416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 195416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 19617ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 19717ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 19817ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 19917ab2063SBarry Smith 20017ab2063SBarry Smith for (i=0; i<m; i++) { 201416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 20217ab2063SBarry Smith #if defined(PETSC_COMPLEX) 203416022c9SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j],real(a->a[j]), 204416022c9SBarry Smith imag(a->a[j])); 20517ab2063SBarry Smith #else 206416022c9SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j], a->a[j]); 20717ab2063SBarry Smith #endif 20817ab2063SBarry Smith } 20917ab2063SBarry Smith } 21017ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 21117ab2063SBarry Smith } 21217ab2063SBarry Smith else { 21317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 21417ab2063SBarry Smith fprintf(fd,"row %d:",i); 215416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 21617ab2063SBarry Smith #if defined(PETSC_COMPLEX) 217416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 218416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 21917ab2063SBarry Smith } 22017ab2063SBarry Smith else { 221416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 22217ab2063SBarry Smith } 22317ab2063SBarry Smith #else 224416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 22517ab2063SBarry Smith #endif 22617ab2063SBarry Smith } 22717ab2063SBarry Smith fprintf(fd,"\n"); 22817ab2063SBarry Smith } 22917ab2063SBarry Smith } 23017ab2063SBarry Smith fflush(fd); 231416022c9SBarry Smith return 0; 232416022c9SBarry Smith } 233416022c9SBarry Smith 234416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 235416022c9SBarry Smith { 236416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 237cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 238cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 239d7e8b826SBarry Smith Draw draw = (Draw) viewer; 240cddf8d76SBarry Smith DrawButton button; 241cddf8d76SBarry Smith 242416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 243416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 244416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 245416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 246cddf8d76SBarry Smith color = DRAW_BLUE; 247416022c9SBarry Smith for ( i=0; i<m; i++ ) { 248cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 249416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 250cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 251cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 252cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 253cddf8d76SBarry Smith #else 254cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 255cddf8d76SBarry Smith #endif 256cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 257cddf8d76SBarry Smith } 258cddf8d76SBarry Smith } 259cddf8d76SBarry Smith color = DRAW_CYAN; 260cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 261cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 262cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 263cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 264cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 265cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 266cddf8d76SBarry Smith } 267cddf8d76SBarry Smith } 268cddf8d76SBarry Smith color = DRAW_RED; 269cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 270cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 271cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 272cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 273cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 274cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 275cddf8d76SBarry Smith #else 276cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 277cddf8d76SBarry Smith #endif 278cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 279416022c9SBarry Smith } 280416022c9SBarry Smith } 281416022c9SBarry Smith DrawFlush(draw); 282cddf8d76SBarry Smith DrawGetPause(draw,&pause); 283cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 284cddf8d76SBarry Smith 285cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 286cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 287cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 288cddf8d76SBarry Smith DrawClear(draw); 289cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 290cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 291cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 292cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 293cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 294cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 295cddf8d76SBarry Smith w *= scale; h *= scale; 296cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 297cddf8d76SBarry Smith color = DRAW_BLUE; 298cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 299cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 300cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 301cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 302cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 303cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 304cddf8d76SBarry Smith #else 305cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 306cddf8d76SBarry Smith #endif 307cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 308cddf8d76SBarry Smith } 309cddf8d76SBarry Smith } 310cddf8d76SBarry Smith color = DRAW_CYAN; 311cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 312cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 313cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 314cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 315cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 316cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 317cddf8d76SBarry Smith } 318cddf8d76SBarry Smith } 319cddf8d76SBarry Smith color = DRAW_RED; 320cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 321cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 322cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 323cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 324cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 325cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 326cddf8d76SBarry Smith #else 327cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 328cddf8d76SBarry Smith #endif 329cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 330cddf8d76SBarry Smith } 331cddf8d76SBarry Smith } 332cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 333cddf8d76SBarry Smith } 334416022c9SBarry Smith return 0; 335416022c9SBarry Smith } 336416022c9SBarry Smith 337416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 338416022c9SBarry Smith { 339416022c9SBarry Smith Mat A = (Mat) obj; 340416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 341416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 342416022c9SBarry Smith 343416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix"); 344416022c9SBarry Smith if (!viewer) { 345416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 346416022c9SBarry Smith } 347416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 348416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 349416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 350416022c9SBarry Smith } 351416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 352416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 353416022c9SBarry Smith } 354416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 355416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 356416022c9SBarry Smith } 357416022c9SBarry Smith } 358416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 359416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 360416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 36117ab2063SBarry Smith } 36217ab2063SBarry Smith return 0; 36317ab2063SBarry Smith } 36441c01911SSatish Balay int Mat_AIJ_CheckInode(Mat); 365416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 36617ab2063SBarry Smith { 367416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 36841c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 369416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 370416022c9SBarry Smith Scalar *aa = a->a, *ap; 37117ab2063SBarry Smith 37217ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 37317ab2063SBarry Smith 37417ab2063SBarry Smith for ( i=1; i<m; i++ ) { 375416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 37617ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 37717ab2063SBarry Smith if (fshift) { 378416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 37917ab2063SBarry Smith N = ailen[i]; 38017ab2063SBarry Smith for ( j=0; j<N; j++ ) { 38117ab2063SBarry Smith ip[j-fshift] = ip[j]; 38217ab2063SBarry Smith ap[j-fshift] = ap[j]; 38317ab2063SBarry Smith } 38417ab2063SBarry Smith } 38517ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 38617ab2063SBarry Smith } 38717ab2063SBarry Smith if (m) { 38817ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 38917ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 39017ab2063SBarry Smith } 39117ab2063SBarry Smith /* reset ilen and imax for each row */ 39217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 39317ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 39417ab2063SBarry Smith } 395416022c9SBarry Smith a->nz = ai[m] + shift; 39617ab2063SBarry Smith 39717ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 398416022c9SBarry Smith if (fshift && a->diag) { 3990452661fSBarry Smith PetscFree(a->diag); 400416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 401416022c9SBarry Smith a->diag = 0; 40217ab2063SBarry Smith } 40376dd722bSSatish Balay /* check out for identical nodes. If found,use inode functions*/ 40441c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 405416022c9SBarry Smith a->assembled = 1; 40617ab2063SBarry Smith return 0; 40717ab2063SBarry Smith } 40817ab2063SBarry Smith 409416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 41017ab2063SBarry Smith { 411416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 412cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 41317ab2063SBarry Smith return 0; 41417ab2063SBarry Smith } 415416022c9SBarry Smith 41617ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 41717ab2063SBarry Smith { 418416022c9SBarry Smith Mat A = (Mat) obj; 419416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 420d5d45c9bSBarry Smith 42117ab2063SBarry Smith #if defined(PETSC_LOG) 422416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 42317ab2063SBarry Smith #endif 4240452661fSBarry Smith PetscFree(a->a); 4250452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4260452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4270452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4280452661fSBarry Smith if (a->imax) PetscFree(a->imax); 4290452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 43076dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 4310452661fSBarry Smith PetscFree(a); 432416022c9SBarry Smith PLogObjectDestroy(A); 4330452661fSBarry Smith PetscHeaderDestroy(A); 43417ab2063SBarry Smith return 0; 43517ab2063SBarry Smith } 43617ab2063SBarry Smith 437416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 43817ab2063SBarry Smith { 43917ab2063SBarry Smith return 0; 44017ab2063SBarry Smith } 44117ab2063SBarry Smith 442416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 44317ab2063SBarry Smith { 444416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 445416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 446416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 447416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 448416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 449e2f28af5SBarry Smith else if (op == ROWS_SORTED || 450e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 451e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 452e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 453df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 454df719601SLois Curfman McInnes else if (op == COLUMN_ORIENTED) 4554d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");} 456df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4574d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 458e2f28af5SBarry Smith else 4594d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 46017ab2063SBarry Smith return 0; 46117ab2063SBarry Smith } 46217ab2063SBarry Smith 463416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 46417ab2063SBarry Smith { 465416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 466416022c9SBarry Smith int i,j, n,shift = a->indexshift; 46717ab2063SBarry Smith Scalar *x, zero = 0.0; 46817ab2063SBarry Smith 469416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix"); 47017ab2063SBarry Smith VecSet(&zero,v); 47117ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 472416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 473416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 474416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 475416022c9SBarry Smith if (a->j[j]+shift == i) { 476416022c9SBarry Smith x[i] = a->a[j]; 47717ab2063SBarry Smith break; 47817ab2063SBarry Smith } 47917ab2063SBarry Smith } 48017ab2063SBarry Smith } 48117ab2063SBarry Smith return 0; 48217ab2063SBarry Smith } 48317ab2063SBarry Smith 48417ab2063SBarry Smith /* -------------------------------------------------------*/ 48517ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 48617ab2063SBarry Smith /* -------------------------------------------------------*/ 487416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 48817ab2063SBarry Smith { 489416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 49017ab2063SBarry Smith Scalar *x, *y, *v, alpha; 491416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 49217ab2063SBarry Smith 493416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix"); 49417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 495cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 49617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 49717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 498416022c9SBarry Smith idx = a->j + a->i[i] + shift; 499416022c9SBarry Smith v = a->a + a->i[i] + shift; 500416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 50117ab2063SBarry Smith alpha = x[i]; 50217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 50317ab2063SBarry Smith } 504416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 50517ab2063SBarry Smith return 0; 50617ab2063SBarry Smith } 507d5d45c9bSBarry Smith 508416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 50917ab2063SBarry Smith { 510416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 51117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 512416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 51317ab2063SBarry Smith 514416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix"); 51517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 51617ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 51717ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 51817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 519416022c9SBarry Smith idx = a->j + a->i[i] + shift; 520416022c9SBarry Smith v = a->a + a->i[i] + shift; 521416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 52217ab2063SBarry Smith alpha = x[i]; 52317ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 52417ab2063SBarry Smith } 52517ab2063SBarry Smith return 0; 52617ab2063SBarry Smith } 52717ab2063SBarry Smith 528416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 52917ab2063SBarry Smith { 530416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 53117ab2063SBarry Smith Scalar *x, *y, *v, sum; 532416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 53317ab2063SBarry Smith 534416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix"); 53517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 53617ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 537416022c9SBarry Smith idx = a->j; 538416022c9SBarry Smith v = a->a; 539416022c9SBarry Smith ii = a->i; 54017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 541416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 54217ab2063SBarry Smith sum = 0.0; 54317ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 54417ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 545416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 54617ab2063SBarry Smith y[i] = sum; 54717ab2063SBarry Smith } 548416022c9SBarry Smith PLogFlops(2*a->nz - m); 54917ab2063SBarry Smith return 0; 55017ab2063SBarry Smith } 55117ab2063SBarry Smith 552416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 55317ab2063SBarry Smith { 554416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 55517ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 556cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 55717ab2063SBarry Smith 55848d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix"); 55917ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 56017ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 561cddf8d76SBarry Smith idx = a->j; 562cddf8d76SBarry Smith v = a->a; 563cddf8d76SBarry Smith ii = a->i; 56417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 565cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 56617ab2063SBarry Smith sum = y[i]; 567cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 568cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 56917ab2063SBarry Smith z[i] = sum; 57017ab2063SBarry Smith } 571416022c9SBarry Smith PLogFlops(2*a->nz); 57217ab2063SBarry Smith return 0; 57317ab2063SBarry Smith } 57417ab2063SBarry Smith 57517ab2063SBarry Smith /* 57617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 57717ab2063SBarry Smith */ 57817ab2063SBarry Smith 579416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 58017ab2063SBarry Smith { 581416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 582416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 58317ab2063SBarry Smith 584416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix"); 5850452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 586416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 587416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 588416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 589416022c9SBarry Smith if (a->j[j]+shift == i) { 59017ab2063SBarry Smith diag[i] = j - shift; 59117ab2063SBarry Smith break; 59217ab2063SBarry Smith } 59317ab2063SBarry Smith } 59417ab2063SBarry Smith } 595416022c9SBarry Smith a->diag = diag; 59617ab2063SBarry Smith return 0; 59717ab2063SBarry Smith } 59817ab2063SBarry Smith 599416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 60017ab2063SBarry Smith double fshift,int its,Vec xx) 60117ab2063SBarry Smith { 602416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 603416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 604d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 60517ab2063SBarry Smith 60617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 607416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 608416022c9SBarry Smith diag = a->diag; 609416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 61017ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 61117ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 61217ab2063SBarry Smith bs = b + shift; 61317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 614416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 615416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 616416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 617416022c9SBarry Smith v = a->a + diag[i] + (!shift); 61817ab2063SBarry Smith sum = b[i]*d/omega; 61917ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 62017ab2063SBarry Smith x[i] = sum; 62117ab2063SBarry Smith } 62217ab2063SBarry Smith return 0; 62317ab2063SBarry Smith } 62417ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 62517ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 62617ab2063SBarry Smith } 627416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 62817ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 62917ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 63017ab2063SBarry Smith 63117ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 63217ab2063SBarry Smith 63317ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 63417ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 63517ab2063SBarry Smith is the relaxation factor. 63617ab2063SBarry Smith */ 6370452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 63817ab2063SBarry Smith scale = (2.0/omega) - 1.0; 63917ab2063SBarry Smith 64017ab2063SBarry Smith /* x = (E + U)^{-1} b */ 64117ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 642416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 643416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 644416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 645416022c9SBarry Smith v = a->a + diag[i] + (!shift); 64617ab2063SBarry Smith sum = b[i]; 64717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 64817ab2063SBarry Smith x[i] = omega*(sum/d); 64917ab2063SBarry Smith } 65017ab2063SBarry Smith 65117ab2063SBarry Smith /* t = b - (2*E - D)x */ 652416022c9SBarry Smith v = a->a; 65317ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 65417ab2063SBarry Smith 65517ab2063SBarry Smith /* t = (E + L)^{-1}t */ 656416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 657416022c9SBarry Smith diag = a->diag; 65817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 659416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 660416022c9SBarry Smith n = diag[i] - a->i[i]; 661416022c9SBarry Smith idx = a->j + a->i[i] + shift; 662416022c9SBarry Smith v = a->a + a->i[i] + shift; 66317ab2063SBarry Smith sum = t[i]; 66417ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 66517ab2063SBarry Smith t[i] = omega*(sum/d); 66617ab2063SBarry Smith } 66717ab2063SBarry Smith 66817ab2063SBarry Smith /* x = x + t */ 66917ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 6700452661fSBarry Smith PetscFree(t); 67117ab2063SBarry Smith return 0; 67217ab2063SBarry Smith } 67317ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 67417ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 67517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 676416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 677416022c9SBarry Smith n = diag[i] - a->i[i]; 678416022c9SBarry Smith idx = a->j + a->i[i] + shift; 679416022c9SBarry Smith v = a->a + a->i[i] + shift; 68017ab2063SBarry Smith sum = b[i]; 68117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 68217ab2063SBarry Smith x[i] = omega*(sum/d); 68317ab2063SBarry Smith } 68417ab2063SBarry Smith xb = x; 68517ab2063SBarry Smith } 68617ab2063SBarry Smith else xb = b; 68717ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 68817ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 68917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 690416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 69117ab2063SBarry Smith } 69217ab2063SBarry Smith } 69317ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 69417ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 695416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 696416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 697416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 698416022c9SBarry Smith v = a->a + diag[i] + (!shift); 69917ab2063SBarry Smith sum = xb[i]; 70017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 70117ab2063SBarry Smith x[i] = omega*(sum/d); 70217ab2063SBarry Smith } 70317ab2063SBarry Smith } 70417ab2063SBarry Smith its--; 70517ab2063SBarry Smith } 70617ab2063SBarry Smith while (its--) { 70717ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 70817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 709416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 710416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 711416022c9SBarry Smith idx = a->j + a->i[i] + shift; 712416022c9SBarry Smith v = a->a + a->i[i] + shift; 71317ab2063SBarry Smith sum = b[i]; 71417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 71517ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 71617ab2063SBarry Smith } 71717ab2063SBarry Smith } 71817ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 71917ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 720416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 721416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 722416022c9SBarry Smith idx = a->j + a->i[i] + shift; 723416022c9SBarry Smith v = a->a + a->i[i] + shift; 72417ab2063SBarry Smith sum = b[i]; 72517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 72617ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 72717ab2063SBarry Smith } 72817ab2063SBarry Smith } 72917ab2063SBarry Smith } 73017ab2063SBarry Smith return 0; 73117ab2063SBarry Smith } 73217ab2063SBarry Smith 733d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 73417ab2063SBarry Smith { 735416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 736416022c9SBarry Smith *nz = a->nz; 737416022c9SBarry Smith *nzalloc = a->maxnz; 738416022c9SBarry Smith *mem = (int)A->mem; 73917ab2063SBarry Smith return 0; 74017ab2063SBarry Smith } 74117ab2063SBarry Smith 74217ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 74317ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 74417ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 74517ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 74617ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 74717ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 74817ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 74917ab2063SBarry Smith 75017ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 75117ab2063SBarry Smith { 752416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 753416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 75417ab2063SBarry Smith 75517ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 75617ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 75717ab2063SBarry Smith if (diag) { 75817ab2063SBarry Smith for ( i=0; i<N; i++ ) { 759416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 760416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 761416022c9SBarry Smith a->ilen[rows[i]] = 1; 762416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 763416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 76417ab2063SBarry Smith } 76517ab2063SBarry Smith else { 76617ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 76717ab2063SBarry Smith CHKERRQ(ierr); 76817ab2063SBarry Smith } 76917ab2063SBarry Smith } 77017ab2063SBarry Smith } 77117ab2063SBarry Smith else { 77217ab2063SBarry Smith for ( i=0; i<N; i++ ) { 773416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 774416022c9SBarry Smith a->ilen[rows[i]] = 0; 77517ab2063SBarry Smith } 77617ab2063SBarry Smith } 77717ab2063SBarry Smith ISRestoreIndices(is,&rows); 77817ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 77917ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 78017ab2063SBarry Smith return 0; 78117ab2063SBarry Smith } 78217ab2063SBarry Smith 783416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 78417ab2063SBarry Smith { 785416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 786416022c9SBarry Smith *m = a->m; *n = a->n; 78717ab2063SBarry Smith return 0; 78817ab2063SBarry Smith } 78917ab2063SBarry Smith 790416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 79117ab2063SBarry Smith { 792416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 793416022c9SBarry Smith *m = 0; *n = a->m; 79417ab2063SBarry Smith return 0; 79517ab2063SBarry Smith } 796416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 79717ab2063SBarry Smith { 798416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 799416022c9SBarry Smith int *itmp,i,ierr,shift = a->indexshift; 80017ab2063SBarry Smith 801416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 80217ab2063SBarry Smith 803416022c9SBarry Smith if (!a->assembled) { 804416022c9SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 805416022c9SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 80617ab2063SBarry Smith } 807416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 808416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 80917ab2063SBarry Smith if (idx) { 81017ab2063SBarry Smith if (*nz) { 811416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8120452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 81317ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 81417ab2063SBarry Smith } 81517ab2063SBarry Smith else *idx = 0; 81617ab2063SBarry Smith } 81717ab2063SBarry Smith return 0; 81817ab2063SBarry Smith } 81917ab2063SBarry Smith 820416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 82117ab2063SBarry Smith { 8220452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 82317ab2063SBarry Smith return 0; 82417ab2063SBarry Smith } 82517ab2063SBarry Smith 826cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 82717ab2063SBarry Smith { 828416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 829416022c9SBarry Smith Scalar *v = a->a; 83017ab2063SBarry Smith double sum = 0.0; 831416022c9SBarry Smith int i, j,shift = a->indexshift; 83217ab2063SBarry Smith 833416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix"); 83417ab2063SBarry Smith if (type == NORM_FROBENIUS) { 835416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 83617ab2063SBarry Smith #if defined(PETSC_COMPLEX) 83717ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 83817ab2063SBarry Smith #else 83917ab2063SBarry Smith sum += (*v)*(*v); v++; 84017ab2063SBarry Smith #endif 84117ab2063SBarry Smith } 84217ab2063SBarry Smith *norm = sqrt(sum); 84317ab2063SBarry Smith } 84417ab2063SBarry Smith else if (type == NORM_1) { 84517ab2063SBarry Smith double *tmp; 846416022c9SBarry Smith int *jj = a->j; 8470452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 848cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 84917ab2063SBarry Smith *norm = 0.0; 850416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 851a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 85217ab2063SBarry Smith } 853416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 85417ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 85517ab2063SBarry Smith } 8560452661fSBarry Smith PetscFree(tmp); 85717ab2063SBarry Smith } 85817ab2063SBarry Smith else if (type == NORM_INFINITY) { 85917ab2063SBarry Smith *norm = 0.0; 860416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 861416022c9SBarry Smith v = a->a + a->i[j] + shift; 86217ab2063SBarry Smith sum = 0.0; 863416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 864cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 86517ab2063SBarry Smith } 86617ab2063SBarry Smith if (sum > *norm) *norm = sum; 86717ab2063SBarry Smith } 86817ab2063SBarry Smith } 86917ab2063SBarry Smith else { 87048d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 87117ab2063SBarry Smith } 87217ab2063SBarry Smith return 0; 87317ab2063SBarry Smith } 87417ab2063SBarry Smith 875416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 87617ab2063SBarry Smith { 877416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 878416022c9SBarry Smith Mat C; 879416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 880416022c9SBarry Smith int shift = a->indexshift; 881d5d45c9bSBarry Smith Scalar *array = a->a; 88217ab2063SBarry Smith 883416022c9SBarry Smith if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 8840452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 885cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 88617ab2063SBarry Smith if (shift) { 88717ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 88817ab2063SBarry Smith } 88917ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 890416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 8910452661fSBarry Smith PetscFree(col); 89217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 89317ab2063SBarry Smith len = ai[i+1]-ai[i]; 894416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 89517ab2063SBarry Smith array += len; aj += len; 89617ab2063SBarry Smith } 89717ab2063SBarry Smith if (shift) { 89817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 89917ab2063SBarry Smith } 90017ab2063SBarry Smith 901416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 902416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 90317ab2063SBarry Smith 904416022c9SBarry Smith if (B) { 905416022c9SBarry Smith *B = C; 90617ab2063SBarry Smith } else { 907416022c9SBarry Smith /* This isn't really an in-place transpose */ 9080452661fSBarry Smith PetscFree(a->a); 9090452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9100452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9110452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9120452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9130452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9140452661fSBarry Smith PetscFree(a); 915416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9160452661fSBarry Smith PetscHeaderDestroy(C); 91717ab2063SBarry Smith } 91817ab2063SBarry Smith return 0; 91917ab2063SBarry Smith } 92017ab2063SBarry Smith 921416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr) 92217ab2063SBarry Smith { 923416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 92417ab2063SBarry Smith Scalar *l,*r,x,*v; 925d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 92617ab2063SBarry Smith 92748d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix"); 92817ab2063SBarry Smith if (ll) { 92917ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 930416022c9SBarry Smith if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length"); 931416022c9SBarry Smith v = a->a; 93217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 93317ab2063SBarry Smith x = l[i]; 934416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 93517ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 93617ab2063SBarry Smith } 93717ab2063SBarry Smith } 93817ab2063SBarry Smith if (rr) { 93917ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 940416022c9SBarry Smith if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length"); 941416022c9SBarry Smith v = a->a; jj = a->j; 94217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 94317ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 94417ab2063SBarry Smith } 94517ab2063SBarry Smith } 94617ab2063SBarry Smith return 0; 94717ab2063SBarry Smith } 94817ab2063SBarry Smith 949cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 95017ab2063SBarry Smith { 951db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 95202834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 953a2744918SBarry Smith register int sum,lensi; 95402834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 955a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 956db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 957416022c9SBarry Smith Mat C; 95817ab2063SBarry Smith 959416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix"); 96017ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 96117ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 96217ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 96317ab2063SBarry Smith 96402834360SBarry Smith if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 96502834360SBarry Smith /* special case of contiguous rows */ 9660452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 96702834360SBarry Smith starts = lens + ncols; 96802834360SBarry Smith /* loop over new rows determining lens and starting points */ 96902834360SBarry Smith for (i=0; i<nrows; i++) { 970a2744918SBarry Smith kstart = ai[irow[i]]+shift; 971a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 97202834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 973a2744918SBarry Smith if (aj[k] >= first) { 97402834360SBarry Smith starts[i] = k; 97502834360SBarry Smith break; 97602834360SBarry Smith } 97702834360SBarry Smith } 978a2744918SBarry Smith sum = 0; 97902834360SBarry Smith while (k < kend) { 980a2744918SBarry Smith if (aj[k++] >= first+ncols) break; 981a2744918SBarry Smith sum++; 98202834360SBarry Smith } 983a2744918SBarry Smith lens[i] = sum; 98402834360SBarry Smith } 98502834360SBarry Smith /* create submatrix */ 986cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 98708480c60SBarry Smith int n_cols,n_rows; 98808480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 98908480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 99008480c60SBarry Smith MatZeroEntries(*B); 99108480c60SBarry Smith C = *B; 99208480c60SBarry Smith } 99308480c60SBarry Smith else { 99402834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 99508480c60SBarry Smith } 996db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 997db02288aSLois Curfman McInnes 99802834360SBarry Smith /* loop over rows inserting into submatrix */ 999db02288aSLois Curfman McInnes a_new = c->a; 1000db02288aSLois Curfman McInnes j_new = c->j; 1001db02288aSLois Curfman McInnes i_new = c->i; 1002db02288aSLois Curfman McInnes i_new[0] = -shift; 100302834360SBarry Smith for (i=0; i<nrows; i++) { 1004a2744918SBarry Smith ii = starts[i]; 1005a2744918SBarry Smith lensi = lens[i]; 1006a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1007a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 100802834360SBarry Smith } 1009a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1010a2744918SBarry Smith a_new += lensi; 1011a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1012a2744918SBarry Smith c->ilen[i] = lensi; 101302834360SBarry Smith } 10140452661fSBarry Smith PetscFree(lens); 101502834360SBarry Smith } 101602834360SBarry Smith else { 101702834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 10180452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 101902834360SBarry Smith ssmap = smap + shift; 10200452661fSBarry Smith cwork = (int *) PetscMalloc((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork); 102102834360SBarry Smith lens = cwork + ncols; 10220452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1023cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 102417ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 102502834360SBarry Smith /* determine lens of each row */ 102602834360SBarry Smith for (i=0; i<nrows; i++) { 102702834360SBarry Smith kstart = a->i[irow[i]]+shift; 102802834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 102902834360SBarry Smith lens[i] = 0; 103002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 103102834360SBarry Smith if (ssmap[a->j[k]]) { 103202834360SBarry Smith lens[i]++; 103302834360SBarry Smith } 103402834360SBarry Smith } 103502834360SBarry Smith } 103617ab2063SBarry Smith /* Create and fill new matrix */ 1037a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 103808480c60SBarry Smith int n_cols,n_rows; 103908480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 104008480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 104108480c60SBarry Smith MatZeroEntries(*B); 104208480c60SBarry Smith C = *B; 104308480c60SBarry Smith } 104408480c60SBarry Smith else { 104502834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 104608480c60SBarry Smith } 104717ab2063SBarry Smith for (i=0; i<nrows; i++) { 104817ab2063SBarry Smith nznew = 0; 1049416022c9SBarry Smith kstart = a->i[irow[i]]+shift; 1050416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 105117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 105202834360SBarry Smith if (ssmap[a->j[k]]) { 105302834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1054416022c9SBarry Smith vwork[nznew++] = a->a[k]; 105517ab2063SBarry Smith } 105617ab2063SBarry Smith } 1057416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 105817ab2063SBarry Smith } 105902834360SBarry Smith /* Free work space */ 106002834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 10610452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 106202834360SBarry Smith } 1063416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1064416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 106517ab2063SBarry Smith 106617ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1067416022c9SBarry Smith *B = C; 106817ab2063SBarry Smith return 0; 106917ab2063SBarry Smith } 107017ab2063SBarry Smith 1071a871dcd8SBarry Smith /* 107263b91edcSBarry Smith note: This can only work for identity for row and col. It would 107363b91edcSBarry Smith be good to check this and otherwise generate an error. 1074a871dcd8SBarry Smith */ 107563b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1076a871dcd8SBarry Smith { 107763b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 107808480c60SBarry Smith int ierr; 107963b91edcSBarry Smith Mat outA; 108063b91edcSBarry Smith 1081a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1082a871dcd8SBarry Smith 108363b91edcSBarry Smith outA = inA; 108463b91edcSBarry Smith inA->factor = FACTOR_LU; 108563b91edcSBarry Smith a->row = row; 108663b91edcSBarry Smith a->col = col; 108763b91edcSBarry Smith 10880452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 108963b91edcSBarry Smith 109008480c60SBarry Smith if (!a->diag) { 109108480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 109263b91edcSBarry Smith } 109363b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1094a871dcd8SBarry Smith return 0; 1095a871dcd8SBarry Smith } 1096a871dcd8SBarry Smith 1097cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1098cddf8d76SBarry Smith Mat **B) 1099cddf8d76SBarry Smith { 1100cddf8d76SBarry Smith int ierr,i; 1101cddf8d76SBarry Smith 1102cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 11030452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1104cddf8d76SBarry Smith } 1105cddf8d76SBarry Smith 1106cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1107cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1108cddf8d76SBarry Smith } 1109cddf8d76SBarry Smith return 0; 1110cddf8d76SBarry Smith } 1111cddf8d76SBarry Smith 111217ab2063SBarry Smith /* -------------------------------------------------------------------*/ 111317ab2063SBarry Smith 111417ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 111517ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1116416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1117416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 111817ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 111917ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 112017ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 112117ab2063SBarry Smith MatRelax_SeqAIJ, 112217ab2063SBarry Smith MatTranspose_SeqAIJ, 112317ab2063SBarry Smith MatGetInfo_SeqAIJ,0, 112417ab2063SBarry Smith MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ, 112517ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 112617ab2063SBarry Smith MatCompress_SeqAIJ, 112717ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 112817ab2063SBarry Smith MatGetReordering_SeqAIJ, 112917ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 113017ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 113117ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 113217ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1133416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 1134a871dcd8SBarry Smith MatCopyPrivate_SeqAIJ,0,0, 1135cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 1136cddf8d76SBarry Smith MatGetSubMatrices_SeqAIJ}; 113717ab2063SBarry Smith 113817ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 113917ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 114017ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 114117ab2063SBarry Smith 114217ab2063SBarry Smith /*@C 114317ab2063SBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ format 114417ab2063SBarry Smith (the default uniprocessor PETSc format). 114517ab2063SBarry Smith 114617ab2063SBarry Smith Input Parameters: 114717ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 114817ab2063SBarry Smith . m - number of rows 114917ab2063SBarry Smith . n - number of columns 115017ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 115117ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 115217ab2063SBarry Smith 115317ab2063SBarry Smith Output Parameter: 1154416022c9SBarry Smith . A - the matrix 115517ab2063SBarry Smith 115617ab2063SBarry Smith Notes: 115717ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 115817ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 11590002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 11600002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 116117ab2063SBarry Smith 116217ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 116317ab2063SBarry Smith Set both nz and nnz to zero for PETSc to control dynamic memory 116417ab2063SBarry Smith allocation. 116517ab2063SBarry Smith 116617ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse 116717ab2063SBarry Smith 116817ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 116917ab2063SBarry Smith @*/ 1170416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 117117ab2063SBarry Smith { 1172416022c9SBarry Smith Mat B; 1173416022c9SBarry Smith Mat_SeqAIJ *b; 117417ab2063SBarry Smith int i,len,ierr; 1175d5d45c9bSBarry Smith 1176416022c9SBarry Smith *A = 0; 11770452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1178416022c9SBarry Smith PLogObjectCreate(B); 11790452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1180416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1181416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1182416022c9SBarry Smith B->view = MatView_SeqAIJ; 1183416022c9SBarry Smith B->factor = 0; 1184416022c9SBarry Smith B->lupivotthreshold = 1.0; 1185416022c9SBarry Smith OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold); 1186416022c9SBarry Smith b->row = 0; 1187416022c9SBarry Smith b->col = 0; 1188416022c9SBarry Smith b->indexshift = 0; 1189416022c9SBarry Smith if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1; 119017ab2063SBarry Smith 1191416022c9SBarry Smith b->m = m; 1192416022c9SBarry Smith b->n = n; 11930452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 119417ab2063SBarry Smith if (!nnz) { 119517ab2063SBarry Smith if (nz <= 0) nz = 1; 1196416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 119717ab2063SBarry Smith nz = nz*m; 119817ab2063SBarry Smith } 119917ab2063SBarry Smith else { 120017ab2063SBarry Smith nz = 0; 1201416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 120217ab2063SBarry Smith } 120317ab2063SBarry Smith 120417ab2063SBarry Smith /* allocate the matrix space */ 1205416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 12060452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1207416022c9SBarry Smith b->j = (int *) (b->a + nz); 1208cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1209416022c9SBarry Smith b->i = b->j + nz; 1210416022c9SBarry Smith b->singlemalloc = 1; 121117ab2063SBarry Smith 1212416022c9SBarry Smith b->i[0] = -b->indexshift; 121317ab2063SBarry Smith for (i=1; i<m+1; i++) { 1214416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 121517ab2063SBarry Smith } 121617ab2063SBarry Smith 1217416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 12180452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1219416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1220416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 122117ab2063SBarry Smith 1222416022c9SBarry Smith b->nz = 0; 1223416022c9SBarry Smith b->maxnz = nz; 1224416022c9SBarry Smith b->sorted = 0; 1225416022c9SBarry Smith b->roworiented = 1; 1226416022c9SBarry Smith b->nonew = 0; 1227416022c9SBarry Smith b->diag = 0; 1228416022c9SBarry Smith b->assembled = 0; 1229416022c9SBarry Smith b->solve_work = 0; 123071bd300dSLois Curfman McInnes b->spptr = 0; 1231754ec7b1SSatish Balay b->inode.node_count = 0; 1232754ec7b1SSatish Balay b->inode.size = 0; 123317ab2063SBarry Smith 1234416022c9SBarry Smith *A = B; 123517ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_superlu")) { 1236416022c9SBarry Smith ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); 123717ab2063SBarry Smith } 123817ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_essl")) { 1239416022c9SBarry Smith ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); 124017ab2063SBarry Smith } 124117ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_dxml")) { 1242416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1243416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 124417ab2063SBarry Smith } 124517ab2063SBarry Smith 124617ab2063SBarry Smith return 0; 124717ab2063SBarry Smith } 124817ab2063SBarry Smith 124908480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues) 125017ab2063SBarry Smith { 1251416022c9SBarry Smith Mat C; 1252416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 125308480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 125417ab2063SBarry Smith 12554043dd9cSLois Curfman McInnes *B = 0; 1256416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix"); 12570452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1258416022c9SBarry Smith PLogObjectCreate(C); 12590452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 126041c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1261416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1262416022c9SBarry Smith C->view = MatView_SeqAIJ; 1263416022c9SBarry Smith C->factor = A->factor; 1264416022c9SBarry Smith c->row = 0; 1265416022c9SBarry Smith c->col = 0; 1266416022c9SBarry Smith c->indexshift = shift; 126717ab2063SBarry Smith 1268416022c9SBarry Smith c->m = a->m; 1269416022c9SBarry Smith c->n = a->n; 127017ab2063SBarry Smith 12710452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 12720452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 127317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1274416022c9SBarry Smith c->imax[i] = a->imax[i]; 1275416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 127617ab2063SBarry Smith } 127717ab2063SBarry Smith 127817ab2063SBarry Smith /* allocate the matrix space */ 1279416022c9SBarry Smith c->singlemalloc = 1; 1280416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 12810452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1282416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1283416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1284416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 128517ab2063SBarry Smith if (m > 0) { 1286416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 128708480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1288416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 128917ab2063SBarry Smith } 129008480c60SBarry Smith } 129117ab2063SBarry Smith 1292416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1293416022c9SBarry Smith c->sorted = a->sorted; 1294416022c9SBarry Smith c->roworiented = a->roworiented; 1295416022c9SBarry Smith c->nonew = a->nonew; 129617ab2063SBarry Smith 1297416022c9SBarry Smith if (a->diag) { 12980452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1299416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 130017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1301416022c9SBarry Smith c->diag[i] = a->diag[i]; 130217ab2063SBarry Smith } 130317ab2063SBarry Smith } 1304416022c9SBarry Smith else c->diag = 0; 1305754ec7b1SSatish Balay if( a->inode.size){ 1306754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1307754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1308754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1309754ec7b1SSatish Balay } else { 1310754ec7b1SSatish Balay c->inode.size = 0; 1311754ec7b1SSatish Balay c->inode.node_count = 0; 1312754ec7b1SSatish Balay } 1313416022c9SBarry Smith c->assembled = 1; 1314416022c9SBarry Smith c->nz = a->nz; 1315416022c9SBarry Smith c->maxnz = a->maxnz; 1316416022c9SBarry Smith c->solve_work = 0; 131776dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1318754ec7b1SSatish Balay 1319416022c9SBarry Smith *B = C; 132017ab2063SBarry Smith return 0; 132117ab2063SBarry Smith } 132217ab2063SBarry Smith 1323416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 132417ab2063SBarry Smith { 1325416022c9SBarry Smith Mat_SeqAIJ *a; 1326416022c9SBarry Smith Mat B; 132717699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 132817ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 132917ab2063SBarry Smith MPI_Comm comm = vobj->comm; 133017ab2063SBarry Smith 133117699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 133217699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 133317ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1334416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 133548d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 133617ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 133717ab2063SBarry Smith 133817ab2063SBarry Smith /* read in row lengths */ 13390452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1340416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 134117ab2063SBarry Smith 134217ab2063SBarry Smith /* create our matrix */ 1343416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1344416022c9SBarry Smith B = *A; 1345416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1346416022c9SBarry Smith shift = a->indexshift; 134717ab2063SBarry Smith 134817ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1349416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 135017ab2063SBarry Smith if (shift) { 135117ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1352416022c9SBarry Smith a->j[i] += 1; 135317ab2063SBarry Smith } 135417ab2063SBarry Smith } 135517ab2063SBarry Smith 135617ab2063SBarry Smith /* read in nonzero values */ 1357416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 135817ab2063SBarry Smith 135917ab2063SBarry Smith /* set matrix "i" values */ 1360416022c9SBarry Smith a->i[0] = -shift; 136117ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1362416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1363416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 136417ab2063SBarry Smith } 13650452661fSBarry Smith PetscFree(rowlengths); 136617ab2063SBarry Smith 1367416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1368416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 136917ab2063SBarry Smith return 0; 137017ab2063SBarry Smith } 137117ab2063SBarry Smith 137217ab2063SBarry Smith 137317ab2063SBarry Smith 1374