117ab2063SBarry Smith #ifndef lint 2*44f6d32bSSatish Balay static char vcid[] = "$Id: aij.c,v 1.114 1995/11/09 22:28:49 bsmith 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; 1917ab2063SBarry Smith 20416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix"); 2117ab2063SBarry Smith 22a2744918SBarry Smith /* 23a2744918SBarry Smith this is tacky: In the future when we have written special factorization 24a2744918SBarry Smith and solve routines for the identity permutation we should use a 25a2744918SBarry Smith stride index set instead of the general one. 26a2744918SBarry Smith */ 27a2744918SBarry Smith if (type == ORDER_NATURAL) { 28a2744918SBarry Smith n = a->n; 29a2744918SBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 30a2744918SBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 31a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 32a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 33a2744918SBarry Smith PetscFree(idx); 34a2744918SBarry Smith ISSetPermutation(*rperm); 35a2744918SBarry Smith ISSetPermutation(*cperm); 36a2744918SBarry Smith ISSetIdentity(*rperm); 37a2744918SBarry Smith ISSetIdentity(*cperm); 38a2744918SBarry Smith return 0; 39a2744918SBarry Smith } 40a2744918SBarry Smith 41416022c9SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr); 42416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 43*44f6d32bSSatish Balay /* ISView(*rperm, STDOUT_VIEWER_SELF);*/ 440452661fSBarry Smith PetscFree(ia); PetscFree(ja); 4517ab2063SBarry Smith return 0; 4617ab2063SBarry Smith } 4717ab2063SBarry Smith 48416022c9SBarry Smith #define CHUNKSIZE 10 4917ab2063SBarry Smith 5017ab2063SBarry Smith /* This version has row oriented v */ 51416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 5217ab2063SBarry Smith { 53416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 54416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 55416022c9SBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen; 56d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 57416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 5817ab2063SBarry Smith 5917ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 60416022c9SBarry Smith row = im[k]; 6117ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 62416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 6317ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 6417ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 65416022c9SBarry Smith low = 0; 6617ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 67416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 68416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 69416022c9SBarry Smith col = in[l] - shift; value = *v++; 70416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 71416022c9SBarry Smith while (high-low > 5) { 72416022c9SBarry Smith t = (low+high)/2; 73416022c9SBarry Smith if (rp[t] > col) high = t; 74416022c9SBarry Smith else low = t; 7517ab2063SBarry Smith } 76416022c9SBarry Smith for ( i=low; i<high; i++ ) { 7717ab2063SBarry Smith if (rp[i] > col) break; 7817ab2063SBarry Smith if (rp[i] == col) { 79416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 8017ab2063SBarry Smith else ap[i] = value; 8117ab2063SBarry Smith goto noinsert; 8217ab2063SBarry Smith } 8317ab2063SBarry Smith } 8417ab2063SBarry Smith if (nonew) goto noinsert; 8517ab2063SBarry Smith if (nrow >= rmax) { 8617ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 87416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 8817ab2063SBarry Smith Scalar *new_a; 8917ab2063SBarry Smith 9017ab2063SBarry Smith /* malloc new storage space */ 91416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 920452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 9317ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 9417ab2063SBarry Smith new_i = new_j + new_nz; 9517ab2063SBarry Smith 9617ab2063SBarry Smith /* copy over old data into new slots */ 9717ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 98416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 99416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 100416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 101416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 10217ab2063SBarry Smith len*sizeof(int)); 103416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 104416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 10517ab2063SBarry Smith len*sizeof(Scalar)); 10617ab2063SBarry Smith /* free up old matrix storage */ 1070452661fSBarry Smith PetscFree(a->a); 1080452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 109416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 110416022c9SBarry Smith a->singlemalloc = 1; 11117ab2063SBarry Smith 11217ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 113416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 114416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 115416022c9SBarry Smith a->maxnz += CHUNKSIZE; 11617ab2063SBarry Smith } 117416022c9SBarry Smith N = nrow++ - 1; a->nz++; 118416022c9SBarry Smith /* shift up all the later entries in this row */ 119416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 12017ab2063SBarry Smith rp[ii+1] = rp[ii]; 12117ab2063SBarry Smith ap[ii+1] = ap[ii]; 12217ab2063SBarry Smith } 12317ab2063SBarry Smith rp[i] = col; 12417ab2063SBarry Smith ap[i] = value; 12517ab2063SBarry Smith noinsert:; 126416022c9SBarry Smith low = i + 1; 12717ab2063SBarry Smith } 12817ab2063SBarry Smith ailen[row] = nrow; 12917ab2063SBarry Smith } 13017ab2063SBarry Smith return 0; 13117ab2063SBarry Smith } 13217ab2063SBarry Smith 13317ab2063SBarry Smith #include "draw.h" 13417ab2063SBarry Smith #include "pinclude/pviewer.h" 135416022c9SBarry Smith #include "sysio.h" 13617ab2063SBarry Smith 137416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 13817ab2063SBarry Smith { 139416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 140416022c9SBarry Smith int i, fd, *col_lens, ierr; 14117ab2063SBarry Smith 142416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1430452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 144416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 145416022c9SBarry Smith col_lens[1] = a->m; 146416022c9SBarry Smith col_lens[2] = a->n; 147416022c9SBarry Smith col_lens[3] = a->nz; 148416022c9SBarry Smith 149416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 150416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 151416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 15217ab2063SBarry Smith } 153416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 1540452661fSBarry Smith PetscFree(col_lens); 155416022c9SBarry Smith 156416022c9SBarry Smith /* store column indices (zero start index) */ 157416022c9SBarry Smith if (a->indexshift) { 158416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 15917ab2063SBarry Smith } 160416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 161416022c9SBarry Smith if (a->indexshift) { 162416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 16317ab2063SBarry Smith } 164416022c9SBarry Smith 165416022c9SBarry Smith /* store nonzero values */ 166416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 16717ab2063SBarry Smith return 0; 16817ab2063SBarry Smith } 169416022c9SBarry Smith 170416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 171416022c9SBarry Smith { 172416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 173416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 17417ab2063SBarry Smith FILE *fd; 17517ab2063SBarry Smith char *outputname; 17617ab2063SBarry Smith 177416022c9SBarry Smith ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 178416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 179416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 18017ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 18108480c60SBarry Smith /* no need to print additional information */ ; 18217ab2063SBarry Smith } 18317ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 18417ab2063SBarry Smith int nz, nzalloc, mem; 185416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 186416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 18717ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 18817ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 18917ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 19017ab2063SBarry Smith 19117ab2063SBarry Smith for (i=0; i<m; i++) { 192416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 19317ab2063SBarry Smith #if defined(PETSC_COMPLEX) 194416022c9SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j],real(a->a[j]), 195416022c9SBarry Smith imag(a->a[j])); 19617ab2063SBarry Smith #else 197416022c9SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j], a->a[j]); 19817ab2063SBarry Smith #endif 19917ab2063SBarry Smith } 20017ab2063SBarry Smith } 20117ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 20217ab2063SBarry Smith } 20317ab2063SBarry Smith else { 20417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 20517ab2063SBarry Smith fprintf(fd,"row %d:",i); 206416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 20717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 208416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 209416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 21017ab2063SBarry Smith } 21117ab2063SBarry Smith else { 212416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 21317ab2063SBarry Smith } 21417ab2063SBarry Smith #else 215416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 21617ab2063SBarry Smith #endif 21717ab2063SBarry Smith } 21817ab2063SBarry Smith fprintf(fd,"\n"); 21917ab2063SBarry Smith } 22017ab2063SBarry Smith } 22117ab2063SBarry Smith fflush(fd); 222416022c9SBarry Smith return 0; 223416022c9SBarry Smith } 224416022c9SBarry Smith 225416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 226416022c9SBarry Smith { 227416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 228cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 229cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 230d7e8b826SBarry Smith Draw draw = (Draw) viewer; 231cddf8d76SBarry Smith DrawButton button; 232cddf8d76SBarry Smith 233416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 234416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 235416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 236416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 237cddf8d76SBarry Smith color = DRAW_BLUE; 238416022c9SBarry Smith for ( i=0; i<m; i++ ) { 239cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 240416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 241cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 242cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 243cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 244cddf8d76SBarry Smith #else 245cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 246cddf8d76SBarry Smith #endif 247cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 248cddf8d76SBarry Smith } 249cddf8d76SBarry Smith } 250cddf8d76SBarry Smith color = DRAW_CYAN; 251cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 252cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 253cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 254cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 255cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 256cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 257cddf8d76SBarry Smith } 258cddf8d76SBarry Smith } 259cddf8d76SBarry Smith color = DRAW_RED; 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 defined(PETSC_COMPLEX) 265cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 266cddf8d76SBarry Smith #else 267cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 268cddf8d76SBarry Smith #endif 269cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 270416022c9SBarry Smith } 271416022c9SBarry Smith } 272416022c9SBarry Smith DrawFlush(draw); 273cddf8d76SBarry Smith DrawGetPause(draw,&pause); 274cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 275cddf8d76SBarry Smith 276cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 277cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 278cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 279cddf8d76SBarry Smith DrawClear(draw); 280cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 281cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 282cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 283cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 284cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 285cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 286cddf8d76SBarry Smith w *= scale; h *= scale; 287cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 288cddf8d76SBarry Smith color = DRAW_BLUE; 289cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 290cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 291cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 292cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 293cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 294cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 295cddf8d76SBarry Smith #else 296cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 297cddf8d76SBarry Smith #endif 298cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 299cddf8d76SBarry Smith } 300cddf8d76SBarry Smith } 301cddf8d76SBarry Smith color = DRAW_CYAN; 302cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 303cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 304cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 305cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 306cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 307cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 308cddf8d76SBarry Smith } 309cddf8d76SBarry Smith } 310cddf8d76SBarry Smith color = DRAW_RED; 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 defined(PETSC_COMPLEX) 316cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 317cddf8d76SBarry Smith #else 318cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 319cddf8d76SBarry Smith #endif 320cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 321cddf8d76SBarry Smith } 322cddf8d76SBarry Smith } 323cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 324cddf8d76SBarry Smith } 325416022c9SBarry Smith return 0; 326416022c9SBarry Smith } 327416022c9SBarry Smith 328416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 329416022c9SBarry Smith { 330416022c9SBarry Smith Mat A = (Mat) obj; 331416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 332416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 333416022c9SBarry Smith 334416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix"); 335416022c9SBarry Smith if (!viewer) { 336416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 337416022c9SBarry Smith } 338416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 339416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 340416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 341416022c9SBarry Smith } 342416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 343416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 344416022c9SBarry Smith } 345416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 346416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 347416022c9SBarry Smith } 348416022c9SBarry Smith } 349416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 350416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 351416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 35217ab2063SBarry Smith } 35317ab2063SBarry Smith return 0; 35417ab2063SBarry Smith } 35541c01911SSatish Balay int Mat_AIJ_CheckInode(Mat); 356416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 35717ab2063SBarry Smith { 358416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 35941c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 360416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 361416022c9SBarry Smith Scalar *aa = a->a, *ap; 36217ab2063SBarry Smith 36317ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 36417ab2063SBarry Smith 36517ab2063SBarry Smith for ( i=1; i<m; i++ ) { 366416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 36717ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 36817ab2063SBarry Smith if (fshift) { 369416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 37017ab2063SBarry Smith N = ailen[i]; 37117ab2063SBarry Smith for ( j=0; j<N; j++ ) { 37217ab2063SBarry Smith ip[j-fshift] = ip[j]; 37317ab2063SBarry Smith ap[j-fshift] = ap[j]; 37417ab2063SBarry Smith } 37517ab2063SBarry Smith } 37617ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 37717ab2063SBarry Smith } 37817ab2063SBarry Smith if (m) { 37917ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 38017ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 38117ab2063SBarry Smith } 38217ab2063SBarry Smith /* reset ilen and imax for each row */ 38317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 38417ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 38517ab2063SBarry Smith } 386416022c9SBarry Smith a->nz = ai[m] + shift; 38717ab2063SBarry Smith 38817ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 389416022c9SBarry Smith if (fshift && a->diag) { 3900452661fSBarry Smith PetscFree(a->diag); 391416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 392416022c9SBarry Smith a->diag = 0; 39317ab2063SBarry Smith } 39476dd722bSSatish Balay /* check out for identical nodes. If found,use inode functions*/ 39541c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 396416022c9SBarry Smith a->assembled = 1; 39717ab2063SBarry Smith return 0; 39817ab2063SBarry Smith } 39917ab2063SBarry Smith 400416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 40117ab2063SBarry Smith { 402416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 403cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 40417ab2063SBarry Smith return 0; 40517ab2063SBarry Smith } 406416022c9SBarry Smith 40717ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 40817ab2063SBarry Smith { 409416022c9SBarry Smith Mat A = (Mat) obj; 410416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 411d5d45c9bSBarry Smith 41217ab2063SBarry Smith #if defined(PETSC_LOG) 413416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 41417ab2063SBarry Smith #endif 4150452661fSBarry Smith PetscFree(a->a); 4160452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4170452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4180452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4190452661fSBarry Smith if (a->imax) PetscFree(a->imax); 4200452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 42176dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 4220452661fSBarry Smith PetscFree(a); 423416022c9SBarry Smith PLogObjectDestroy(A); 4240452661fSBarry Smith PetscHeaderDestroy(A); 42517ab2063SBarry Smith return 0; 42617ab2063SBarry Smith } 42717ab2063SBarry Smith 428416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 42917ab2063SBarry Smith { 43017ab2063SBarry Smith return 0; 43117ab2063SBarry Smith } 43217ab2063SBarry Smith 433416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 43417ab2063SBarry Smith { 435416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 436416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 437416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 438416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 439416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 440e2f28af5SBarry Smith else if (op == ROWS_SORTED || 441e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 442e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 443e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 444df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 445df719601SLois Curfman McInnes else if (op == COLUMN_ORIENTED) 4464d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");} 447df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4484d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 449e2f28af5SBarry Smith else 4504d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 45117ab2063SBarry Smith return 0; 45217ab2063SBarry Smith } 45317ab2063SBarry Smith 454416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 45517ab2063SBarry Smith { 456416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 457416022c9SBarry Smith int i,j, n,shift = a->indexshift; 45817ab2063SBarry Smith Scalar *x, zero = 0.0; 45917ab2063SBarry Smith 460416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix"); 46117ab2063SBarry Smith VecSet(&zero,v); 46217ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 463416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 464416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 465416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 466416022c9SBarry Smith if (a->j[j]+shift == i) { 467416022c9SBarry Smith x[i] = a->a[j]; 46817ab2063SBarry Smith break; 46917ab2063SBarry Smith } 47017ab2063SBarry Smith } 47117ab2063SBarry Smith } 47217ab2063SBarry Smith return 0; 47317ab2063SBarry Smith } 47417ab2063SBarry Smith 47517ab2063SBarry Smith /* -------------------------------------------------------*/ 47617ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 47717ab2063SBarry Smith /* -------------------------------------------------------*/ 478416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 47917ab2063SBarry Smith { 480416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 48117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 482416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 48317ab2063SBarry Smith 484416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix"); 48517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 486cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 48717ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 48817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 489416022c9SBarry Smith idx = a->j + a->i[i] + shift; 490416022c9SBarry Smith v = a->a + a->i[i] + shift; 491416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 49217ab2063SBarry Smith alpha = x[i]; 49317ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 49417ab2063SBarry Smith } 495416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 49617ab2063SBarry Smith return 0; 49717ab2063SBarry Smith } 498d5d45c9bSBarry Smith 499416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 50017ab2063SBarry Smith { 501416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 50217ab2063SBarry Smith Scalar *x, *y, *v, alpha; 503416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 50417ab2063SBarry Smith 505416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix"); 50617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 50717ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 50817ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 50917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 510416022c9SBarry Smith idx = a->j + a->i[i] + shift; 511416022c9SBarry Smith v = a->a + a->i[i] + shift; 512416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 51317ab2063SBarry Smith alpha = x[i]; 51417ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 51517ab2063SBarry Smith } 51617ab2063SBarry Smith return 0; 51717ab2063SBarry Smith } 51817ab2063SBarry Smith 519416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 52017ab2063SBarry Smith { 521416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 52217ab2063SBarry Smith Scalar *x, *y, *v, sum; 523416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 52417ab2063SBarry Smith 525416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix"); 52617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 52717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 528416022c9SBarry Smith idx = a->j; 529416022c9SBarry Smith v = a->a; 530416022c9SBarry Smith ii = a->i; 53117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 532416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 53317ab2063SBarry Smith sum = 0.0; 53417ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 53517ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 536416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 53717ab2063SBarry Smith y[i] = sum; 53817ab2063SBarry Smith } 539416022c9SBarry Smith PLogFlops(2*a->nz - m); 54017ab2063SBarry Smith return 0; 54117ab2063SBarry Smith } 54217ab2063SBarry Smith 543416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 54417ab2063SBarry Smith { 545416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 54617ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 547cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 54817ab2063SBarry Smith 54948d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix"); 55017ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 55117ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 552cddf8d76SBarry Smith idx = a->j; 553cddf8d76SBarry Smith v = a->a; 554cddf8d76SBarry Smith ii = a->i; 55517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 556cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 55717ab2063SBarry Smith sum = y[i]; 558cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 559cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 56017ab2063SBarry Smith z[i] = sum; 56117ab2063SBarry Smith } 562416022c9SBarry Smith PLogFlops(2*a->nz); 56317ab2063SBarry Smith return 0; 56417ab2063SBarry Smith } 56517ab2063SBarry Smith 56617ab2063SBarry Smith /* 56717ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 56817ab2063SBarry Smith */ 56917ab2063SBarry Smith 570416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 57117ab2063SBarry Smith { 572416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 573416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 57417ab2063SBarry Smith 575416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix"); 5760452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 577416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 578416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 579416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 580416022c9SBarry Smith if (a->j[j]+shift == i) { 58117ab2063SBarry Smith diag[i] = j - shift; 58217ab2063SBarry Smith break; 58317ab2063SBarry Smith } 58417ab2063SBarry Smith } 58517ab2063SBarry Smith } 586416022c9SBarry Smith a->diag = diag; 58717ab2063SBarry Smith return 0; 58817ab2063SBarry Smith } 58917ab2063SBarry Smith 590416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 59117ab2063SBarry Smith double fshift,int its,Vec xx) 59217ab2063SBarry Smith { 593416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 594416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 595d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 59617ab2063SBarry Smith 59717ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 598416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 599416022c9SBarry Smith diag = a->diag; 600416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 60117ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 60217ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 60317ab2063SBarry Smith bs = b + shift; 60417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 605416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 606416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 607416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 608416022c9SBarry Smith v = a->a + diag[i] + (!shift); 60917ab2063SBarry Smith sum = b[i]*d/omega; 61017ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 61117ab2063SBarry Smith x[i] = sum; 61217ab2063SBarry Smith } 61317ab2063SBarry Smith return 0; 61417ab2063SBarry Smith } 61517ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 61617ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 61717ab2063SBarry Smith } 618416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 61917ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 62017ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 62117ab2063SBarry Smith 62217ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 62317ab2063SBarry Smith 62417ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 62517ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 62617ab2063SBarry Smith is the relaxation factor. 62717ab2063SBarry Smith */ 6280452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 62917ab2063SBarry Smith scale = (2.0/omega) - 1.0; 63017ab2063SBarry Smith 63117ab2063SBarry Smith /* x = (E + U)^{-1} b */ 63217ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 633416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 634416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 635416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 636416022c9SBarry Smith v = a->a + diag[i] + (!shift); 63717ab2063SBarry Smith sum = b[i]; 63817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 63917ab2063SBarry Smith x[i] = omega*(sum/d); 64017ab2063SBarry Smith } 64117ab2063SBarry Smith 64217ab2063SBarry Smith /* t = b - (2*E - D)x */ 643416022c9SBarry Smith v = a->a; 64417ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 64517ab2063SBarry Smith 64617ab2063SBarry Smith /* t = (E + L)^{-1}t */ 647416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 648416022c9SBarry Smith diag = a->diag; 64917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 650416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 651416022c9SBarry Smith n = diag[i] - a->i[i]; 652416022c9SBarry Smith idx = a->j + a->i[i] + shift; 653416022c9SBarry Smith v = a->a + a->i[i] + shift; 65417ab2063SBarry Smith sum = t[i]; 65517ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 65617ab2063SBarry Smith t[i] = omega*(sum/d); 65717ab2063SBarry Smith } 65817ab2063SBarry Smith 65917ab2063SBarry Smith /* x = x + t */ 66017ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 6610452661fSBarry Smith PetscFree(t); 66217ab2063SBarry Smith return 0; 66317ab2063SBarry Smith } 66417ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 66517ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 66617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 667416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 668416022c9SBarry Smith n = diag[i] - a->i[i]; 669416022c9SBarry Smith idx = a->j + a->i[i] + shift; 670416022c9SBarry Smith v = a->a + a->i[i] + shift; 67117ab2063SBarry Smith sum = b[i]; 67217ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 67317ab2063SBarry Smith x[i] = omega*(sum/d); 67417ab2063SBarry Smith } 67517ab2063SBarry Smith xb = x; 67617ab2063SBarry Smith } 67717ab2063SBarry Smith else xb = b; 67817ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 67917ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 68017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 681416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 68217ab2063SBarry Smith } 68317ab2063SBarry Smith } 68417ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 68517ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 686416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 687416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 688416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 689416022c9SBarry Smith v = a->a + diag[i] + (!shift); 69017ab2063SBarry Smith sum = xb[i]; 69117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 69217ab2063SBarry Smith x[i] = omega*(sum/d); 69317ab2063SBarry Smith } 69417ab2063SBarry Smith } 69517ab2063SBarry Smith its--; 69617ab2063SBarry Smith } 69717ab2063SBarry Smith while (its--) { 69817ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 69917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 700416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 701416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 702416022c9SBarry Smith idx = a->j + a->i[i] + shift; 703416022c9SBarry Smith v = a->a + a->i[i] + shift; 70417ab2063SBarry Smith sum = b[i]; 70517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 70617ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 70717ab2063SBarry Smith } 70817ab2063SBarry Smith } 70917ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 71017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 711416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 712416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 713416022c9SBarry Smith idx = a->j + a->i[i] + shift; 714416022c9SBarry Smith v = a->a + a->i[i] + shift; 71517ab2063SBarry Smith sum = b[i]; 71617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 71717ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 71817ab2063SBarry Smith } 71917ab2063SBarry Smith } 72017ab2063SBarry Smith } 72117ab2063SBarry Smith return 0; 72217ab2063SBarry Smith } 72317ab2063SBarry Smith 724d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 72517ab2063SBarry Smith { 726416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 727416022c9SBarry Smith *nz = a->nz; 728416022c9SBarry Smith *nzalloc = a->maxnz; 729416022c9SBarry Smith *mem = (int)A->mem; 73017ab2063SBarry Smith return 0; 73117ab2063SBarry Smith } 73217ab2063SBarry Smith 73317ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 73417ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 73517ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 73617ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 73717ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 73817ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 73917ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 74017ab2063SBarry Smith 74117ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 74217ab2063SBarry Smith { 743416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 744416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 74517ab2063SBarry Smith 74617ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 74717ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 74817ab2063SBarry Smith if (diag) { 74917ab2063SBarry Smith for ( i=0; i<N; i++ ) { 750416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 751416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 752416022c9SBarry Smith a->ilen[rows[i]] = 1; 753416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 754416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 75517ab2063SBarry Smith } 75617ab2063SBarry Smith else { 75717ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 75817ab2063SBarry Smith CHKERRQ(ierr); 75917ab2063SBarry Smith } 76017ab2063SBarry Smith } 76117ab2063SBarry Smith } 76217ab2063SBarry Smith else { 76317ab2063SBarry Smith for ( i=0; i<N; i++ ) { 764416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 765416022c9SBarry Smith a->ilen[rows[i]] = 0; 76617ab2063SBarry Smith } 76717ab2063SBarry Smith } 76817ab2063SBarry Smith ISRestoreIndices(is,&rows); 76917ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 77017ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 77117ab2063SBarry Smith return 0; 77217ab2063SBarry Smith } 77317ab2063SBarry Smith 774416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 77517ab2063SBarry Smith { 776416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 777416022c9SBarry Smith *m = a->m; *n = a->n; 77817ab2063SBarry Smith return 0; 77917ab2063SBarry Smith } 78017ab2063SBarry Smith 781416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 78217ab2063SBarry Smith { 783416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 784416022c9SBarry Smith *m = 0; *n = a->m; 78517ab2063SBarry Smith return 0; 78617ab2063SBarry Smith } 787416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 78817ab2063SBarry Smith { 789416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 790416022c9SBarry Smith int *itmp,i,ierr,shift = a->indexshift; 79117ab2063SBarry Smith 792416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 79317ab2063SBarry Smith 794416022c9SBarry Smith if (!a->assembled) { 795416022c9SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 796416022c9SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 79717ab2063SBarry Smith } 798416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 799416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 80017ab2063SBarry Smith if (idx) { 80117ab2063SBarry Smith if (*nz) { 802416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8030452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 80417ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 80517ab2063SBarry Smith } 80617ab2063SBarry Smith else *idx = 0; 80717ab2063SBarry Smith } 80817ab2063SBarry Smith return 0; 80917ab2063SBarry Smith } 81017ab2063SBarry Smith 811416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 81217ab2063SBarry Smith { 8130452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 81417ab2063SBarry Smith return 0; 81517ab2063SBarry Smith } 81617ab2063SBarry Smith 817cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 81817ab2063SBarry Smith { 819416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 820416022c9SBarry Smith Scalar *v = a->a; 82117ab2063SBarry Smith double sum = 0.0; 822416022c9SBarry Smith int i, j,shift = a->indexshift; 82317ab2063SBarry Smith 824416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix"); 82517ab2063SBarry Smith if (type == NORM_FROBENIUS) { 826416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 82717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 82817ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 82917ab2063SBarry Smith #else 83017ab2063SBarry Smith sum += (*v)*(*v); v++; 83117ab2063SBarry Smith #endif 83217ab2063SBarry Smith } 83317ab2063SBarry Smith *norm = sqrt(sum); 83417ab2063SBarry Smith } 83517ab2063SBarry Smith else if (type == NORM_1) { 83617ab2063SBarry Smith double *tmp; 837416022c9SBarry Smith int *jj = a->j; 8380452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 839cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 84017ab2063SBarry Smith *norm = 0.0; 841416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 842a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 84317ab2063SBarry Smith } 844416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 84517ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 84617ab2063SBarry Smith } 8470452661fSBarry Smith PetscFree(tmp); 84817ab2063SBarry Smith } 84917ab2063SBarry Smith else if (type == NORM_INFINITY) { 85017ab2063SBarry Smith *norm = 0.0; 851416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 852416022c9SBarry Smith v = a->a + a->i[j] + shift; 85317ab2063SBarry Smith sum = 0.0; 854416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 855cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 85617ab2063SBarry Smith } 85717ab2063SBarry Smith if (sum > *norm) *norm = sum; 85817ab2063SBarry Smith } 85917ab2063SBarry Smith } 86017ab2063SBarry Smith else { 86148d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 86217ab2063SBarry Smith } 86317ab2063SBarry Smith return 0; 86417ab2063SBarry Smith } 86517ab2063SBarry Smith 866416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 86717ab2063SBarry Smith { 868416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 869416022c9SBarry Smith Mat C; 870416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 871416022c9SBarry Smith int shift = a->indexshift; 872d5d45c9bSBarry Smith Scalar *array = a->a; 87317ab2063SBarry Smith 874416022c9SBarry Smith if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 8750452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 876cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 87717ab2063SBarry Smith if (shift) { 87817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 87917ab2063SBarry Smith } 88017ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 881416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 8820452661fSBarry Smith PetscFree(col); 88317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 88417ab2063SBarry Smith len = ai[i+1]-ai[i]; 885416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 88617ab2063SBarry Smith array += len; aj += len; 88717ab2063SBarry Smith } 88817ab2063SBarry Smith if (shift) { 88917ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 89017ab2063SBarry Smith } 89117ab2063SBarry Smith 892416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 893416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 89417ab2063SBarry Smith 895416022c9SBarry Smith if (B) { 896416022c9SBarry Smith *B = C; 89717ab2063SBarry Smith } else { 898416022c9SBarry Smith /* This isn't really an in-place transpose */ 8990452661fSBarry Smith PetscFree(a->a); 9000452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9010452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9020452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9030452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9040452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9050452661fSBarry Smith PetscFree(a); 906416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9070452661fSBarry Smith PetscHeaderDestroy(C); 90817ab2063SBarry Smith } 90917ab2063SBarry Smith return 0; 91017ab2063SBarry Smith } 91117ab2063SBarry Smith 912416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr) 91317ab2063SBarry Smith { 914416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 91517ab2063SBarry Smith Scalar *l,*r,x,*v; 916d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 91717ab2063SBarry Smith 91848d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix"); 91917ab2063SBarry Smith if (ll) { 92017ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 921416022c9SBarry Smith if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length"); 922416022c9SBarry Smith v = a->a; 92317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 92417ab2063SBarry Smith x = l[i]; 925416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 92617ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 92717ab2063SBarry Smith } 92817ab2063SBarry Smith } 92917ab2063SBarry Smith if (rr) { 93017ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 931416022c9SBarry Smith if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length"); 932416022c9SBarry Smith v = a->a; jj = a->j; 93317ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 93417ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 93517ab2063SBarry Smith } 93617ab2063SBarry Smith } 93717ab2063SBarry Smith return 0; 93817ab2063SBarry Smith } 93917ab2063SBarry Smith 940cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 94117ab2063SBarry Smith { 942db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 94302834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 944a2744918SBarry Smith register int sum,lensi; 94502834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 946a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 947db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 948416022c9SBarry Smith Mat C; 94917ab2063SBarry Smith 950416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix"); 95117ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 95217ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 95317ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 95417ab2063SBarry Smith 95502834360SBarry Smith if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 95602834360SBarry Smith /* special case of contiguous rows */ 9570452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 95802834360SBarry Smith starts = lens + ncols; 95902834360SBarry Smith /* loop over new rows determining lens and starting points */ 96002834360SBarry Smith for (i=0; i<nrows; i++) { 961a2744918SBarry Smith kstart = ai[irow[i]]+shift; 962a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 96302834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 964a2744918SBarry Smith if (aj[k] >= first) { 96502834360SBarry Smith starts[i] = k; 96602834360SBarry Smith break; 96702834360SBarry Smith } 96802834360SBarry Smith } 969a2744918SBarry Smith sum = 0; 97002834360SBarry Smith while (k < kend) { 971a2744918SBarry Smith if (aj[k++] >= first+ncols) break; 972a2744918SBarry Smith sum++; 97302834360SBarry Smith } 974a2744918SBarry Smith lens[i] = sum; 97502834360SBarry Smith } 97602834360SBarry Smith /* create submatrix */ 977cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 97808480c60SBarry Smith int n_cols,n_rows; 97908480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 98008480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 98108480c60SBarry Smith MatZeroEntries(*B); 98208480c60SBarry Smith C = *B; 98308480c60SBarry Smith } 98408480c60SBarry Smith else { 98502834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 98608480c60SBarry Smith } 987db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 988db02288aSLois Curfman McInnes 98902834360SBarry Smith /* loop over rows inserting into submatrix */ 990db02288aSLois Curfman McInnes a_new = c->a; 991db02288aSLois Curfman McInnes j_new = c->j; 992db02288aSLois Curfman McInnes i_new = c->i; 993db02288aSLois Curfman McInnes i_new[0] = -shift; 99402834360SBarry Smith for (i=0; i<nrows; i++) { 995a2744918SBarry Smith ii = starts[i]; 996a2744918SBarry Smith lensi = lens[i]; 997a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 998a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 99902834360SBarry Smith } 1000a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1001a2744918SBarry Smith a_new += lensi; 1002a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1003a2744918SBarry Smith c->ilen[i] = lensi; 100402834360SBarry Smith } 10050452661fSBarry Smith PetscFree(lens); 100602834360SBarry Smith } 100702834360SBarry Smith else { 100802834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 10090452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 101002834360SBarry Smith ssmap = smap + shift; 10110452661fSBarry Smith cwork = (int *) PetscMalloc((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork); 101202834360SBarry Smith lens = cwork + ncols; 10130452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1014cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 101517ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 101602834360SBarry Smith /* determine lens of each row */ 101702834360SBarry Smith for (i=0; i<nrows; i++) { 101802834360SBarry Smith kstart = a->i[irow[i]]+shift; 101902834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 102002834360SBarry Smith lens[i] = 0; 102102834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 102202834360SBarry Smith if (ssmap[a->j[k]]) { 102302834360SBarry Smith lens[i]++; 102402834360SBarry Smith } 102502834360SBarry Smith } 102602834360SBarry Smith } 102717ab2063SBarry Smith /* Create and fill new matrix */ 1028a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 102908480c60SBarry Smith int n_cols,n_rows; 103008480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 103108480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 103208480c60SBarry Smith MatZeroEntries(*B); 103308480c60SBarry Smith C = *B; 103408480c60SBarry Smith } 103508480c60SBarry Smith else { 103602834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 103708480c60SBarry Smith } 103817ab2063SBarry Smith for (i=0; i<nrows; i++) { 103917ab2063SBarry Smith nznew = 0; 1040416022c9SBarry Smith kstart = a->i[irow[i]]+shift; 1041416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 104217ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 104302834360SBarry Smith if (ssmap[a->j[k]]) { 104402834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1045416022c9SBarry Smith vwork[nznew++] = a->a[k]; 104617ab2063SBarry Smith } 104717ab2063SBarry Smith } 1048416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 104917ab2063SBarry Smith } 105002834360SBarry Smith /* Free work space */ 105102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 10520452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 105302834360SBarry Smith } 1054416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1055416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 105617ab2063SBarry Smith 105717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1058416022c9SBarry Smith *B = C; 105917ab2063SBarry Smith return 0; 106017ab2063SBarry Smith } 106117ab2063SBarry Smith 1062a871dcd8SBarry Smith /* 106363b91edcSBarry Smith note: This can only work for identity for row and col. It would 106463b91edcSBarry Smith be good to check this and otherwise generate an error. 1065a871dcd8SBarry Smith */ 106663b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1067a871dcd8SBarry Smith { 106863b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 106908480c60SBarry Smith int ierr; 107063b91edcSBarry Smith Mat outA; 107163b91edcSBarry Smith 1072a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1073a871dcd8SBarry Smith 107463b91edcSBarry Smith outA = inA; 107563b91edcSBarry Smith inA->factor = FACTOR_LU; 107663b91edcSBarry Smith a->row = row; 107763b91edcSBarry Smith a->col = col; 107863b91edcSBarry Smith 10790452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 108063b91edcSBarry Smith 108108480c60SBarry Smith if (!a->diag) { 108208480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 108363b91edcSBarry Smith } 108463b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1085a871dcd8SBarry Smith return 0; 1086a871dcd8SBarry Smith } 1087a871dcd8SBarry Smith 1088cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1089cddf8d76SBarry Smith Mat **B) 1090cddf8d76SBarry Smith { 1091cddf8d76SBarry Smith int ierr,i; 1092cddf8d76SBarry Smith 1093cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 10940452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1095cddf8d76SBarry Smith } 1096cddf8d76SBarry Smith 1097cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1098cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1099cddf8d76SBarry Smith } 1100cddf8d76SBarry Smith return 0; 1101cddf8d76SBarry Smith } 1102cddf8d76SBarry Smith 110317ab2063SBarry Smith /* -------------------------------------------------------------------*/ 110417ab2063SBarry Smith 110517ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 110617ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1107416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1108416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 110917ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 111017ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 111117ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 111217ab2063SBarry Smith MatRelax_SeqAIJ, 111317ab2063SBarry Smith MatTranspose_SeqAIJ, 111417ab2063SBarry Smith MatGetInfo_SeqAIJ,0, 111517ab2063SBarry Smith MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ, 111617ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 111717ab2063SBarry Smith MatCompress_SeqAIJ, 111817ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 111917ab2063SBarry Smith MatGetReordering_SeqAIJ, 112017ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 112117ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 112217ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 112317ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1124416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 1125a871dcd8SBarry Smith MatCopyPrivate_SeqAIJ,0,0, 1126cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 1127cddf8d76SBarry Smith MatGetSubMatrices_SeqAIJ}; 112817ab2063SBarry Smith 112917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 113017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 113117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 113217ab2063SBarry Smith 113317ab2063SBarry Smith /*@C 113417ab2063SBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ format 113517ab2063SBarry Smith (the default uniprocessor PETSc format). 113617ab2063SBarry Smith 113717ab2063SBarry Smith Input Parameters: 113817ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 113917ab2063SBarry Smith . m - number of rows 114017ab2063SBarry Smith . n - number of columns 114117ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 114217ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 114317ab2063SBarry Smith 114417ab2063SBarry Smith Output Parameter: 1145416022c9SBarry Smith . A - the matrix 114617ab2063SBarry Smith 114717ab2063SBarry Smith Notes: 114817ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 114917ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 11500002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 11510002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 115217ab2063SBarry Smith 115317ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 115417ab2063SBarry Smith Set both nz and nnz to zero for PETSc to control dynamic memory 115517ab2063SBarry Smith allocation. 115617ab2063SBarry Smith 115717ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse 115817ab2063SBarry Smith 115917ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 116017ab2063SBarry Smith @*/ 1161416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 116217ab2063SBarry Smith { 1163416022c9SBarry Smith Mat B; 1164416022c9SBarry Smith Mat_SeqAIJ *b; 116517ab2063SBarry Smith int i,len,ierr; 1166d5d45c9bSBarry Smith 1167416022c9SBarry Smith *A = 0; 11680452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1169416022c9SBarry Smith PLogObjectCreate(B); 11700452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1171416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1172416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1173416022c9SBarry Smith B->view = MatView_SeqAIJ; 1174416022c9SBarry Smith B->factor = 0; 1175416022c9SBarry Smith B->lupivotthreshold = 1.0; 1176416022c9SBarry Smith OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold); 1177416022c9SBarry Smith b->row = 0; 1178416022c9SBarry Smith b->col = 0; 1179416022c9SBarry Smith b->indexshift = 0; 1180416022c9SBarry Smith if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1; 118117ab2063SBarry Smith 1182416022c9SBarry Smith b->m = m; 1183416022c9SBarry Smith b->n = n; 11840452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 118517ab2063SBarry Smith if (!nnz) { 118617ab2063SBarry Smith if (nz <= 0) nz = 1; 1187416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 118817ab2063SBarry Smith nz = nz*m; 118917ab2063SBarry Smith } 119017ab2063SBarry Smith else { 119117ab2063SBarry Smith nz = 0; 1192416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 119317ab2063SBarry Smith } 119417ab2063SBarry Smith 119517ab2063SBarry Smith /* allocate the matrix space */ 1196416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 11970452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1198416022c9SBarry Smith b->j = (int *) (b->a + nz); 1199cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1200416022c9SBarry Smith b->i = b->j + nz; 1201416022c9SBarry Smith b->singlemalloc = 1; 120217ab2063SBarry Smith 1203416022c9SBarry Smith b->i[0] = -b->indexshift; 120417ab2063SBarry Smith for (i=1; i<m+1; i++) { 1205416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 120617ab2063SBarry Smith } 120717ab2063SBarry Smith 1208416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 12090452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1210416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1211416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 121217ab2063SBarry Smith 1213416022c9SBarry Smith b->nz = 0; 1214416022c9SBarry Smith b->maxnz = nz; 1215416022c9SBarry Smith b->sorted = 0; 1216416022c9SBarry Smith b->roworiented = 1; 1217416022c9SBarry Smith b->nonew = 0; 1218416022c9SBarry Smith b->diag = 0; 1219416022c9SBarry Smith b->assembled = 0; 1220416022c9SBarry Smith b->solve_work = 0; 122171bd300dSLois Curfman McInnes b->spptr = 0; 1222754ec7b1SSatish Balay b->inode.node_count = 0; 1223754ec7b1SSatish Balay b->inode.size = 0; 122417ab2063SBarry Smith 1225416022c9SBarry Smith *A = B; 122617ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_superlu")) { 1227416022c9SBarry Smith ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); 122817ab2063SBarry Smith } 122917ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_essl")) { 1230416022c9SBarry Smith ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); 123117ab2063SBarry Smith } 123217ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_dxml")) { 1233416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1234416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 123517ab2063SBarry Smith } 123617ab2063SBarry Smith 123717ab2063SBarry Smith return 0; 123817ab2063SBarry Smith } 123917ab2063SBarry Smith 124008480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues) 124117ab2063SBarry Smith { 1242416022c9SBarry Smith Mat C; 1243416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 124408480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 124517ab2063SBarry Smith 12464043dd9cSLois Curfman McInnes *B = 0; 1247416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix"); 12480452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1249416022c9SBarry Smith PLogObjectCreate(C); 12500452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 125141c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1252416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1253416022c9SBarry Smith C->view = MatView_SeqAIJ; 1254416022c9SBarry Smith C->factor = A->factor; 1255416022c9SBarry Smith c->row = 0; 1256416022c9SBarry Smith c->col = 0; 1257416022c9SBarry Smith c->indexshift = shift; 125817ab2063SBarry Smith 1259416022c9SBarry Smith c->m = a->m; 1260416022c9SBarry Smith c->n = a->n; 126117ab2063SBarry Smith 12620452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 12630452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 126417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1265416022c9SBarry Smith c->imax[i] = a->imax[i]; 1266416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 126717ab2063SBarry Smith } 126817ab2063SBarry Smith 126917ab2063SBarry Smith /* allocate the matrix space */ 1270416022c9SBarry Smith c->singlemalloc = 1; 1271416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 12720452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1273416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1274416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1275416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 127617ab2063SBarry Smith if (m > 0) { 1277416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 127808480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1279416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 128017ab2063SBarry Smith } 128108480c60SBarry Smith } 128217ab2063SBarry Smith 1283416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1284416022c9SBarry Smith c->sorted = a->sorted; 1285416022c9SBarry Smith c->roworiented = a->roworiented; 1286416022c9SBarry Smith c->nonew = a->nonew; 128717ab2063SBarry Smith 1288416022c9SBarry Smith if (a->diag) { 12890452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1290416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 129117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1292416022c9SBarry Smith c->diag[i] = a->diag[i]; 129317ab2063SBarry Smith } 129417ab2063SBarry Smith } 1295416022c9SBarry Smith else c->diag = 0; 1296754ec7b1SSatish Balay if( a->inode.size){ 1297754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1298754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1299754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1300754ec7b1SSatish Balay } else { 1301754ec7b1SSatish Balay c->inode.size = 0; 1302754ec7b1SSatish Balay c->inode.node_count = 0; 1303754ec7b1SSatish Balay } 1304416022c9SBarry Smith c->assembled = 1; 1305416022c9SBarry Smith c->nz = a->nz; 1306416022c9SBarry Smith c->maxnz = a->maxnz; 1307416022c9SBarry Smith c->solve_work = 0; 130876dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1309754ec7b1SSatish Balay 1310416022c9SBarry Smith *B = C; 131117ab2063SBarry Smith return 0; 131217ab2063SBarry Smith } 131317ab2063SBarry Smith 1314416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 131517ab2063SBarry Smith { 1316416022c9SBarry Smith Mat_SeqAIJ *a; 1317416022c9SBarry Smith Mat B; 131817699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 131917ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 132017ab2063SBarry Smith MPI_Comm comm = vobj->comm; 132117ab2063SBarry Smith 132217699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 132317699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 132417ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1325416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 132648d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 132717ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 132817ab2063SBarry Smith 132917ab2063SBarry Smith /* read in row lengths */ 13300452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1331416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 133217ab2063SBarry Smith 133317ab2063SBarry Smith /* create our matrix */ 1334416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1335416022c9SBarry Smith B = *A; 1336416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1337416022c9SBarry Smith shift = a->indexshift; 133817ab2063SBarry Smith 133917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1340416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 134117ab2063SBarry Smith if (shift) { 134217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1343416022c9SBarry Smith a->j[i] += 1; 134417ab2063SBarry Smith } 134517ab2063SBarry Smith } 134617ab2063SBarry Smith 134717ab2063SBarry Smith /* read in nonzero values */ 1348416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 134917ab2063SBarry Smith 135017ab2063SBarry Smith /* set matrix "i" values */ 1351416022c9SBarry Smith a->i[0] = -shift; 135217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1353416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1354416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 135517ab2063SBarry Smith } 13560452661fSBarry Smith PetscFree(rowlengths); 135717ab2063SBarry Smith 1358416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1359416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 136017ab2063SBarry Smith return 0; 136117ab2063SBarry Smith } 136217ab2063SBarry Smith 136317ab2063SBarry Smith 136417ab2063SBarry Smith 1365