117ab2063SBarry Smith #ifndef lint 2*0452661fSBarry Smith static char vcid[] = "$Id: aij.c,v 1.107 1995/11/01 19:10:07 bsmith Exp bsmith $"; 317ab2063SBarry Smith #endif 417ab2063SBarry Smith 517ab2063SBarry Smith #include "aij.h" 617ab2063SBarry Smith #include "vec/vecimpl.h" 717ab2063SBarry Smith #include "inline/spops.h" 817ab2063SBarry Smith 917ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**); 1017ab2063SBarry Smith 11416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1217ab2063SBarry Smith { 13416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1417ab2063SBarry Smith int ierr, *ia, *ja; 1517ab2063SBarry Smith 16416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix"); 1717ab2063SBarry Smith 18416022c9SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr); 19416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 20*0452661fSBarry Smith PetscFree(ia); PetscFree(ja); 2117ab2063SBarry Smith return 0; 2217ab2063SBarry Smith } 2317ab2063SBarry Smith 24416022c9SBarry Smith #define CHUNKSIZE 10 2517ab2063SBarry Smith 2617ab2063SBarry Smith /* This version has row oriented v */ 27416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 2817ab2063SBarry Smith { 29416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 30416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 31416022c9SBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen; 32416022c9SBarry Smith int *aj = a->j, nonew = a->nonew; 33416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 34416022c9SBarry Smith int shift = a->indexshift; 3517ab2063SBarry Smith 3617ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 37416022c9SBarry Smith row = im[k]; 3817ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 39416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 4017ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 4117ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 42416022c9SBarry Smith low = 0; 4317ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 44416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 45416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 46416022c9SBarry Smith col = in[l] - shift; value = *v++; 47416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 48416022c9SBarry Smith while (high-low > 5) { 49416022c9SBarry Smith t = (low+high)/2; 50416022c9SBarry Smith if (rp[t] > col) high = t; 51416022c9SBarry Smith else low = t; 5217ab2063SBarry Smith } 53416022c9SBarry Smith for ( i=low; i<high; i++ ) { 5417ab2063SBarry Smith if (rp[i] > col) break; 5517ab2063SBarry Smith if (rp[i] == col) { 56416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 5717ab2063SBarry Smith else ap[i] = value; 5817ab2063SBarry Smith goto noinsert; 5917ab2063SBarry Smith } 6017ab2063SBarry Smith } 6117ab2063SBarry Smith if (nonew) goto noinsert; 6217ab2063SBarry Smith if (nrow >= rmax) { 6317ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 64416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 6517ab2063SBarry Smith Scalar *new_a; 6617ab2063SBarry Smith 6717ab2063SBarry Smith /* malloc new storage space */ 68416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 69*0452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 7017ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 7117ab2063SBarry Smith new_i = new_j + new_nz; 7217ab2063SBarry Smith 7317ab2063SBarry Smith /* copy over old data into new slots */ 7417ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 75416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 76416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 77416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 78416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 7917ab2063SBarry Smith len*sizeof(int)); 80416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 81416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 8217ab2063SBarry Smith len*sizeof(Scalar)); 8317ab2063SBarry Smith /* free up old matrix storage */ 84*0452661fSBarry Smith PetscFree(a->a); 85*0452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 86416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 87416022c9SBarry Smith a->singlemalloc = 1; 8817ab2063SBarry Smith 8917ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 90416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 91416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 92416022c9SBarry Smith a->maxnz += CHUNKSIZE; 9317ab2063SBarry Smith } 94416022c9SBarry Smith N = nrow++ - 1; a->nz++; 95416022c9SBarry Smith /* shift up all the later entries in this row */ 96416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 9717ab2063SBarry Smith rp[ii+1] = rp[ii]; 9817ab2063SBarry Smith ap[ii+1] = ap[ii]; 9917ab2063SBarry Smith } 10017ab2063SBarry Smith rp[i] = col; 10117ab2063SBarry Smith ap[i] = value; 10217ab2063SBarry Smith noinsert:; 103416022c9SBarry Smith low = i + 1; 10417ab2063SBarry Smith } 10517ab2063SBarry Smith ailen[row] = nrow; 10617ab2063SBarry Smith } 10717ab2063SBarry Smith return 0; 10817ab2063SBarry Smith } 10917ab2063SBarry Smith 11017ab2063SBarry Smith #include "draw.h" 11117ab2063SBarry Smith #include "pinclude/pviewer.h" 112416022c9SBarry Smith #include "sysio.h" 11317ab2063SBarry Smith 114416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 11517ab2063SBarry Smith { 116416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 117416022c9SBarry Smith int i, fd, *col_lens, ierr; 11817ab2063SBarry Smith 119416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 120*0452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 121416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 122416022c9SBarry Smith col_lens[1] = a->m; 123416022c9SBarry Smith col_lens[2] = a->n; 124416022c9SBarry Smith col_lens[3] = a->nz; 125416022c9SBarry Smith 126416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 127416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 128416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 12917ab2063SBarry Smith } 130416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 131*0452661fSBarry Smith PetscFree(col_lens); 132416022c9SBarry Smith 133416022c9SBarry Smith /* store column indices (zero start index) */ 134416022c9SBarry Smith if (a->indexshift) { 135416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 13617ab2063SBarry Smith } 137416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 138416022c9SBarry Smith if (a->indexshift) { 139416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 14017ab2063SBarry Smith } 141416022c9SBarry Smith 142416022c9SBarry Smith /* store nonzero values */ 143416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 14417ab2063SBarry Smith return 0; 14517ab2063SBarry Smith } 146416022c9SBarry Smith 147416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 148416022c9SBarry Smith { 149416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 150416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 15117ab2063SBarry Smith FILE *fd; 15217ab2063SBarry Smith char *outputname; 15317ab2063SBarry Smith 154416022c9SBarry Smith ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 155416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 156416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 15717ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 15808480c60SBarry Smith /* no need to print additional information */ ; 15917ab2063SBarry Smith } 16017ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 16117ab2063SBarry Smith int nz, nzalloc, mem; 162416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 163416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 16417ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 16517ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 16617ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 16717ab2063SBarry Smith 16817ab2063SBarry Smith for (i=0; i<m; i++) { 169416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 17017ab2063SBarry Smith #if defined(PETSC_COMPLEX) 171416022c9SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j],real(a->a[j]), 172416022c9SBarry Smith imag(a->a[j])); 17317ab2063SBarry Smith #else 174416022c9SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j], a->a[j]); 17517ab2063SBarry Smith #endif 17617ab2063SBarry Smith } 17717ab2063SBarry Smith } 17817ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 17917ab2063SBarry Smith } 18017ab2063SBarry Smith else { 18117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 18217ab2063SBarry Smith fprintf(fd,"row %d:",i); 183416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 18417ab2063SBarry Smith #if defined(PETSC_COMPLEX) 185416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 186416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 18717ab2063SBarry Smith } 18817ab2063SBarry Smith else { 189416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 19017ab2063SBarry Smith } 19117ab2063SBarry Smith #else 192416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 19317ab2063SBarry Smith #endif 19417ab2063SBarry Smith } 19517ab2063SBarry Smith fprintf(fd,"\n"); 19617ab2063SBarry Smith } 19717ab2063SBarry Smith } 19817ab2063SBarry Smith fflush(fd); 199416022c9SBarry Smith return 0; 200416022c9SBarry Smith } 201416022c9SBarry Smith 202416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 203416022c9SBarry Smith { 204416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 205cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 206cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 207416022c9SBarry Smith DrawCtx draw = (DrawCtx) viewer; 208cddf8d76SBarry Smith DrawButton button; 209cddf8d76SBarry Smith 210416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 211416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 212416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 213416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 214cddf8d76SBarry Smith color = DRAW_BLUE; 215416022c9SBarry Smith for ( i=0; i<m; i++ ) { 216cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 217416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 218cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 219cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 220cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 221cddf8d76SBarry Smith #else 222cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 223cddf8d76SBarry Smith #endif 224cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 225cddf8d76SBarry Smith } 226cddf8d76SBarry Smith } 227cddf8d76SBarry Smith color = DRAW_CYAN; 228cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 229cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 230cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 231cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 232cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 233cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 234cddf8d76SBarry Smith } 235cddf8d76SBarry Smith } 236cddf8d76SBarry Smith color = DRAW_RED; 237cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 238cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 239cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 240cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 241cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 242cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 243cddf8d76SBarry Smith #else 244cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 245cddf8d76SBarry Smith #endif 246cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 247416022c9SBarry Smith } 248416022c9SBarry Smith } 249416022c9SBarry Smith DrawFlush(draw); 250cddf8d76SBarry Smith DrawGetPause(draw,&pause); 251cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 252cddf8d76SBarry Smith 253cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 254cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 255cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 256cddf8d76SBarry Smith DrawClear(draw); 257cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 258cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 259cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 260cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 261cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 262cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 263cddf8d76SBarry Smith w *= scale; h *= scale; 264cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 265cddf8d76SBarry Smith color = DRAW_BLUE; 266cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 267cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 268cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 269cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 270cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 271cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 272cddf8d76SBarry Smith #else 273cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 274cddf8d76SBarry Smith #endif 275cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 276cddf8d76SBarry Smith } 277cddf8d76SBarry Smith } 278cddf8d76SBarry Smith color = DRAW_CYAN; 279cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 280cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 281cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 282cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 283cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 284cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 285cddf8d76SBarry Smith } 286cddf8d76SBarry Smith } 287cddf8d76SBarry Smith color = DRAW_RED; 288cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 289cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 290cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 291cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 292cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 293cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 294cddf8d76SBarry Smith #else 295cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 296cddf8d76SBarry Smith #endif 297cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 298cddf8d76SBarry Smith } 299cddf8d76SBarry Smith } 300cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 301cddf8d76SBarry Smith } 302416022c9SBarry Smith return 0; 303416022c9SBarry Smith } 304416022c9SBarry Smith 305416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 306416022c9SBarry Smith { 307416022c9SBarry Smith Mat A = (Mat) obj; 308416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 309416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 310416022c9SBarry Smith 311416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix"); 312416022c9SBarry Smith if (!viewer) { 313416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 314416022c9SBarry Smith } 315416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 316416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 317416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 318416022c9SBarry Smith } 319416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 320416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 321416022c9SBarry Smith } 322416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 323416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 324416022c9SBarry Smith } 325416022c9SBarry Smith } 326416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 327416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 328416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 32917ab2063SBarry Smith } 33017ab2063SBarry Smith return 0; 33117ab2063SBarry Smith } 33217ab2063SBarry Smith 333416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 33417ab2063SBarry Smith { 335416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 336416022c9SBarry Smith int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax; 337416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 338416022c9SBarry Smith Scalar *aa = a->a, *ap; 33917ab2063SBarry Smith 34017ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 34117ab2063SBarry Smith 34217ab2063SBarry Smith for ( i=1; i<m; i++ ) { 343416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 34417ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 34517ab2063SBarry Smith if (fshift) { 346416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 34717ab2063SBarry Smith N = ailen[i]; 34817ab2063SBarry Smith for ( j=0; j<N; j++ ) { 34917ab2063SBarry Smith ip[j-fshift] = ip[j]; 35017ab2063SBarry Smith ap[j-fshift] = ap[j]; 35117ab2063SBarry Smith } 35217ab2063SBarry Smith } 35317ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 35417ab2063SBarry Smith } 35517ab2063SBarry Smith if (m) { 35617ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 35717ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 35817ab2063SBarry Smith } 35917ab2063SBarry Smith /* reset ilen and imax for each row */ 36017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 36117ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 36217ab2063SBarry Smith } 363416022c9SBarry Smith a->nz = ai[m] + shift; 36417ab2063SBarry Smith 36517ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 366416022c9SBarry Smith if (fshift && a->diag) { 367*0452661fSBarry Smith PetscFree(a->diag); 368416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 369416022c9SBarry Smith a->diag = 0; 37017ab2063SBarry Smith } 371416022c9SBarry Smith a->assembled = 1; 37217ab2063SBarry Smith return 0; 37317ab2063SBarry Smith } 37417ab2063SBarry Smith 375416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 37617ab2063SBarry Smith { 377416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 378cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 37917ab2063SBarry Smith return 0; 38017ab2063SBarry Smith } 381416022c9SBarry Smith 38217ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 38317ab2063SBarry Smith { 384416022c9SBarry Smith Mat A = (Mat) obj; 385416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 38617ab2063SBarry Smith #if defined(PETSC_LOG) 387416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 38817ab2063SBarry Smith #endif 389*0452661fSBarry Smith PetscFree(a->a); 390*0452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 391*0452661fSBarry Smith if (a->diag) PetscFree(a->diag); 392*0452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 393*0452661fSBarry Smith if (a->imax) PetscFree(a->imax); 394*0452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 395*0452661fSBarry Smith PetscFree(a); 396416022c9SBarry Smith PLogObjectDestroy(A); 397*0452661fSBarry Smith PetscHeaderDestroy(A); 39817ab2063SBarry Smith return 0; 39917ab2063SBarry Smith } 40017ab2063SBarry Smith 401416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 40217ab2063SBarry Smith { 40317ab2063SBarry Smith return 0; 40417ab2063SBarry Smith } 40517ab2063SBarry Smith 406416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 40717ab2063SBarry Smith { 408416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 409416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 410416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 411416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 412416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 413e2f28af5SBarry Smith else if (op == ROWS_SORTED || 414e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 415e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 416e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 417df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 418df719601SLois Curfman McInnes else if (op == COLUMN_ORIENTED) 4194d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");} 420df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4214d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 422e2f28af5SBarry Smith else 4234d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 42417ab2063SBarry Smith return 0; 42517ab2063SBarry Smith } 42617ab2063SBarry Smith 427416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 42817ab2063SBarry Smith { 429416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 430416022c9SBarry Smith int i,j, n,shift = a->indexshift; 43117ab2063SBarry Smith Scalar *x, zero = 0.0; 43217ab2063SBarry Smith 433416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix"); 43417ab2063SBarry Smith VecSet(&zero,v); 43517ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 436416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 437416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 438416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 439416022c9SBarry Smith if (a->j[j]+shift == i) { 440416022c9SBarry Smith x[i] = a->a[j]; 44117ab2063SBarry Smith break; 44217ab2063SBarry Smith } 44317ab2063SBarry Smith } 44417ab2063SBarry Smith } 44517ab2063SBarry Smith return 0; 44617ab2063SBarry Smith } 44717ab2063SBarry Smith 44817ab2063SBarry Smith /* -------------------------------------------------------*/ 44917ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 45017ab2063SBarry Smith /* -------------------------------------------------------*/ 451416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 45217ab2063SBarry Smith { 453416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 45417ab2063SBarry Smith Scalar *x, *y, *v, alpha; 455416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 45617ab2063SBarry Smith 457416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix"); 45817ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 459cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 46017ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 46117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 462416022c9SBarry Smith idx = a->j + a->i[i] + shift; 463416022c9SBarry Smith v = a->a + a->i[i] + shift; 464416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 46517ab2063SBarry Smith alpha = x[i]; 46617ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 46717ab2063SBarry Smith } 468416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 46917ab2063SBarry Smith return 0; 47017ab2063SBarry Smith } 471416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 47217ab2063SBarry Smith { 473416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 47417ab2063SBarry Smith Scalar *x, *y, *v, alpha; 475416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 47617ab2063SBarry Smith 477416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix"); 47817ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 47917ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 48017ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 48117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 482416022c9SBarry Smith idx = a->j + a->i[i] + shift; 483416022c9SBarry Smith v = a->a + a->i[i] + shift; 484416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 48517ab2063SBarry Smith alpha = x[i]; 48617ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 48717ab2063SBarry Smith } 48817ab2063SBarry Smith return 0; 48917ab2063SBarry Smith } 49017ab2063SBarry Smith 491416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 49217ab2063SBarry Smith { 493416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 49417ab2063SBarry Smith Scalar *x, *y, *v, sum; 495416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 49617ab2063SBarry Smith 497416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix"); 49817ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 49917ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 500416022c9SBarry Smith idx = a->j; 501416022c9SBarry Smith v = a->a; 502416022c9SBarry Smith ii = a->i; 50317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 504416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 50517ab2063SBarry Smith sum = 0.0; 50617ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 50717ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 508416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 50917ab2063SBarry Smith y[i] = sum; 51017ab2063SBarry Smith } 511416022c9SBarry Smith PLogFlops(2*a->nz - m); 51217ab2063SBarry Smith return 0; 51317ab2063SBarry Smith } 51417ab2063SBarry Smith 515416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 51617ab2063SBarry Smith { 517416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 51817ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 519cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 52017ab2063SBarry Smith 52148d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix"); 52217ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 52317ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 524cddf8d76SBarry Smith idx = a->j; 525cddf8d76SBarry Smith v = a->a; 526cddf8d76SBarry Smith ii = a->i; 52717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 528cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 52917ab2063SBarry Smith sum = y[i]; 530cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 531cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 53217ab2063SBarry Smith z[i] = sum; 53317ab2063SBarry Smith } 534416022c9SBarry Smith PLogFlops(2*a->nz); 53517ab2063SBarry Smith return 0; 53617ab2063SBarry Smith } 53717ab2063SBarry Smith 53817ab2063SBarry Smith /* 53917ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 54017ab2063SBarry Smith */ 54117ab2063SBarry Smith 542416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 54317ab2063SBarry Smith { 544416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 545416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 54617ab2063SBarry Smith 547416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix"); 548*0452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 549416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 550416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 551416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 552416022c9SBarry Smith if (a->j[j]+shift == i) { 55317ab2063SBarry Smith diag[i] = j - shift; 55417ab2063SBarry Smith break; 55517ab2063SBarry Smith } 55617ab2063SBarry Smith } 55717ab2063SBarry Smith } 558416022c9SBarry Smith a->diag = diag; 55917ab2063SBarry Smith return 0; 56017ab2063SBarry Smith } 56117ab2063SBarry Smith 562416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 56317ab2063SBarry Smith double fshift,int its,Vec xx) 56417ab2063SBarry Smith { 565416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 566416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 567416022c9SBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i; 568416022c9SBarry Smith int shift = a->indexshift; 56917ab2063SBarry Smith 57017ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 571416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 572416022c9SBarry Smith diag = a->diag; 573416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 57417ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 57517ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 57617ab2063SBarry Smith bs = b + shift; 57717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 578416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 579416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 580416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 581416022c9SBarry Smith v = a->a + diag[i] + (!shift); 58217ab2063SBarry Smith sum = b[i]*d/omega; 58317ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 58417ab2063SBarry Smith x[i] = sum; 58517ab2063SBarry Smith } 58617ab2063SBarry Smith return 0; 58717ab2063SBarry Smith } 58817ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 58917ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 59017ab2063SBarry Smith } 591416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 59217ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 59317ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 59417ab2063SBarry Smith 59517ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 59617ab2063SBarry Smith 59717ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 59817ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 59917ab2063SBarry Smith is the relaxation factor. 60017ab2063SBarry Smith */ 601*0452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 60217ab2063SBarry Smith scale = (2.0/omega) - 1.0; 60317ab2063SBarry Smith 60417ab2063SBarry Smith /* x = (E + U)^{-1} b */ 60517ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 606416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 607416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 608416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 609416022c9SBarry Smith v = a->a + diag[i] + (!shift); 61017ab2063SBarry Smith sum = b[i]; 61117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 61217ab2063SBarry Smith x[i] = omega*(sum/d); 61317ab2063SBarry Smith } 61417ab2063SBarry Smith 61517ab2063SBarry Smith /* t = b - (2*E - D)x */ 616416022c9SBarry Smith v = a->a; 61717ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 61817ab2063SBarry Smith 61917ab2063SBarry Smith /* t = (E + L)^{-1}t */ 620416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 621416022c9SBarry Smith diag = a->diag; 62217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 623416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 624416022c9SBarry Smith n = diag[i] - a->i[i]; 625416022c9SBarry Smith idx = a->j + a->i[i] + shift; 626416022c9SBarry Smith v = a->a + a->i[i] + shift; 62717ab2063SBarry Smith sum = t[i]; 62817ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 62917ab2063SBarry Smith t[i] = omega*(sum/d); 63017ab2063SBarry Smith } 63117ab2063SBarry Smith 63217ab2063SBarry Smith /* x = x + t */ 63317ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 634*0452661fSBarry Smith PetscFree(t); 63517ab2063SBarry Smith return 0; 63617ab2063SBarry Smith } 63717ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 63817ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 63917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 640416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 641416022c9SBarry Smith n = diag[i] - a->i[i]; 642416022c9SBarry Smith idx = a->j + a->i[i] + shift; 643416022c9SBarry Smith v = a->a + a->i[i] + shift; 64417ab2063SBarry Smith sum = b[i]; 64517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 64617ab2063SBarry Smith x[i] = omega*(sum/d); 64717ab2063SBarry Smith } 64817ab2063SBarry Smith xb = x; 64917ab2063SBarry Smith } 65017ab2063SBarry Smith else xb = b; 65117ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 65217ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 65317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 654416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 65517ab2063SBarry Smith } 65617ab2063SBarry Smith } 65717ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 65817ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 659416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 660416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 661416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 662416022c9SBarry Smith v = a->a + diag[i] + (!shift); 66317ab2063SBarry Smith sum = xb[i]; 66417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 66517ab2063SBarry Smith x[i] = omega*(sum/d); 66617ab2063SBarry Smith } 66717ab2063SBarry Smith } 66817ab2063SBarry Smith its--; 66917ab2063SBarry Smith } 67017ab2063SBarry Smith while (its--) { 67117ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 67217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 673416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 674416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 675416022c9SBarry Smith idx = a->j + a->i[i] + shift; 676416022c9SBarry Smith v = a->a + a->i[i] + shift; 67717ab2063SBarry Smith sum = b[i]; 67817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 67917ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 68017ab2063SBarry Smith } 68117ab2063SBarry Smith } 68217ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 68317ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 684416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 685416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 686416022c9SBarry Smith idx = a->j + a->i[i] + shift; 687416022c9SBarry Smith v = a->a + a->i[i] + shift; 68817ab2063SBarry Smith sum = b[i]; 68917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 69017ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 69117ab2063SBarry Smith } 69217ab2063SBarry Smith } 69317ab2063SBarry Smith } 69417ab2063SBarry Smith return 0; 69517ab2063SBarry Smith } 69617ab2063SBarry Smith 697416022c9SBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz, 69817ab2063SBarry Smith int *nzalloc,int *mem) 69917ab2063SBarry Smith { 700416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 701416022c9SBarry Smith *nz = a->nz; 702416022c9SBarry Smith *nzalloc = a->maxnz; 703416022c9SBarry Smith *mem = (int)A->mem; 70417ab2063SBarry Smith return 0; 70517ab2063SBarry Smith } 70617ab2063SBarry Smith 70717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 70817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 70917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 71017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 71117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 71217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 71317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 71417ab2063SBarry Smith 71517ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 71617ab2063SBarry Smith { 717416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 718416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 71917ab2063SBarry Smith 72017ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 72117ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 72217ab2063SBarry Smith if (diag) { 72317ab2063SBarry Smith for ( i=0; i<N; i++ ) { 724416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 725416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 726416022c9SBarry Smith a->ilen[rows[i]] = 1; 727416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 728416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 72917ab2063SBarry Smith } 73017ab2063SBarry Smith else { 73117ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 73217ab2063SBarry Smith CHKERRQ(ierr); 73317ab2063SBarry Smith } 73417ab2063SBarry Smith } 73517ab2063SBarry Smith } 73617ab2063SBarry Smith else { 73717ab2063SBarry Smith for ( i=0; i<N; i++ ) { 738416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 739416022c9SBarry Smith a->ilen[rows[i]] = 0; 74017ab2063SBarry Smith } 74117ab2063SBarry Smith } 74217ab2063SBarry Smith ISRestoreIndices(is,&rows); 74317ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 74417ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 74517ab2063SBarry Smith return 0; 74617ab2063SBarry Smith } 74717ab2063SBarry Smith 748416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 74917ab2063SBarry Smith { 750416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 751416022c9SBarry Smith *m = a->m; *n = a->n; 75217ab2063SBarry Smith return 0; 75317ab2063SBarry Smith } 75417ab2063SBarry Smith 755416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 75617ab2063SBarry Smith { 757416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 758416022c9SBarry Smith *m = 0; *n = a->m; 75917ab2063SBarry Smith return 0; 76017ab2063SBarry Smith } 761416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 76217ab2063SBarry Smith { 763416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 764416022c9SBarry Smith int *itmp,i,ierr,shift = a->indexshift; 76517ab2063SBarry Smith 766416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 76717ab2063SBarry Smith 768416022c9SBarry Smith if (!a->assembled) { 769416022c9SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 770416022c9SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 77117ab2063SBarry Smith } 772416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 773416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 77417ab2063SBarry Smith if (idx) { 77517ab2063SBarry Smith if (*nz) { 776416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 777*0452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 77817ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 77917ab2063SBarry Smith } 78017ab2063SBarry Smith else *idx = 0; 78117ab2063SBarry Smith } 78217ab2063SBarry Smith return 0; 78317ab2063SBarry Smith } 78417ab2063SBarry Smith 785416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 78617ab2063SBarry Smith { 787*0452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 78817ab2063SBarry Smith return 0; 78917ab2063SBarry Smith } 79017ab2063SBarry Smith 791cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 79217ab2063SBarry Smith { 793416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 794416022c9SBarry Smith Scalar *v = a->a; 79517ab2063SBarry Smith double sum = 0.0; 796416022c9SBarry Smith int i, j,shift = a->indexshift; 79717ab2063SBarry Smith 798416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix"); 79917ab2063SBarry Smith if (type == NORM_FROBENIUS) { 800416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 80117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 80217ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 80317ab2063SBarry Smith #else 80417ab2063SBarry Smith sum += (*v)*(*v); v++; 80517ab2063SBarry Smith #endif 80617ab2063SBarry Smith } 80717ab2063SBarry Smith *norm = sqrt(sum); 80817ab2063SBarry Smith } 80917ab2063SBarry Smith else if (type == NORM_1) { 81017ab2063SBarry Smith double *tmp; 811416022c9SBarry Smith int *jj = a->j; 812*0452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 813cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 81417ab2063SBarry Smith *norm = 0.0; 815416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 816cddf8d76SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v++); 81717ab2063SBarry Smith } 818416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 81917ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 82017ab2063SBarry Smith } 821*0452661fSBarry Smith PetscFree(tmp); 82217ab2063SBarry Smith } 82317ab2063SBarry Smith else if (type == NORM_INFINITY) { 82417ab2063SBarry Smith *norm = 0.0; 825416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 826416022c9SBarry Smith v = a->a + a->i[j] + shift; 82717ab2063SBarry Smith sum = 0.0; 828416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 829cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 83017ab2063SBarry Smith } 83117ab2063SBarry Smith if (sum > *norm) *norm = sum; 83217ab2063SBarry Smith } 83317ab2063SBarry Smith } 83417ab2063SBarry Smith else { 83548d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 83617ab2063SBarry Smith } 83717ab2063SBarry Smith return 0; 83817ab2063SBarry Smith } 83917ab2063SBarry Smith 840416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 84117ab2063SBarry Smith { 842416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 843416022c9SBarry Smith Mat C; 844416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 845416022c9SBarry Smith Scalar *array = a->a; 846416022c9SBarry Smith int shift = a->indexshift; 84717ab2063SBarry Smith 848416022c9SBarry Smith if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 849*0452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 850cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 85117ab2063SBarry Smith if (shift) { 85217ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 85317ab2063SBarry Smith } 85417ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 855416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 856*0452661fSBarry Smith PetscFree(col); 85717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 85817ab2063SBarry Smith len = ai[i+1]-ai[i]; 859416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 86017ab2063SBarry Smith array += len; aj += len; 86117ab2063SBarry Smith } 86217ab2063SBarry Smith if (shift) { 86317ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 86417ab2063SBarry Smith } 86517ab2063SBarry Smith 866416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 867416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 86817ab2063SBarry Smith 869416022c9SBarry Smith if (B) { 870416022c9SBarry Smith *B = C; 87117ab2063SBarry Smith } else { 872416022c9SBarry Smith /* This isn't really an in-place transpose */ 873*0452661fSBarry Smith PetscFree(a->a); 874*0452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 875*0452661fSBarry Smith if (a->diag) PetscFree(a->diag); 876*0452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 877*0452661fSBarry Smith if (a->imax) PetscFree(a->imax); 878*0452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 879*0452661fSBarry Smith PetscFree(a); 880416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 881*0452661fSBarry Smith PetscHeaderDestroy(C); 88217ab2063SBarry Smith } 88317ab2063SBarry Smith return 0; 88417ab2063SBarry Smith } 88517ab2063SBarry Smith 886416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr) 88717ab2063SBarry Smith { 888416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 88917ab2063SBarry Smith Scalar *l,*r,x,*v; 890416022c9SBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj; 891416022c9SBarry Smith int shift = a->indexshift; 89217ab2063SBarry Smith 89348d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix"); 89417ab2063SBarry Smith if (ll) { 89517ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 896416022c9SBarry Smith if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length"); 897416022c9SBarry Smith v = a->a; 89817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 89917ab2063SBarry Smith x = l[i]; 900416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 90117ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 90217ab2063SBarry Smith } 90317ab2063SBarry Smith } 90417ab2063SBarry Smith if (rr) { 90517ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 906416022c9SBarry Smith if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length"); 907416022c9SBarry Smith v = a->a; jj = a->j; 90817ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 90917ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 91017ab2063SBarry Smith } 91117ab2063SBarry Smith } 91217ab2063SBarry Smith return 0; 91317ab2063SBarry Smith } 91417ab2063SBarry Smith 915cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 91617ab2063SBarry Smith { 917db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 91802834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 91902834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 920db02288aSLois Curfman McInnes int first,step,*starts,*j_new,*i_new; 921db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 922416022c9SBarry Smith Mat C; 92317ab2063SBarry Smith 924416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix"); 92517ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 92617ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 92717ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 92817ab2063SBarry Smith 92902834360SBarry Smith if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 93002834360SBarry Smith /* special case of contiguous rows */ 931*0452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 93202834360SBarry Smith starts = lens + ncols; 93302834360SBarry Smith /* loop over new rows determining lens and starting points */ 93402834360SBarry Smith for (i=0; i<nrows; i++) { 93502834360SBarry Smith kstart = a->i[irow[i]]+shift; 93602834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 93702834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 93802834360SBarry Smith if (a->j[k] >= first) { 93902834360SBarry Smith starts[i] = k; 94002834360SBarry Smith break; 94102834360SBarry Smith } 94202834360SBarry Smith } 94302834360SBarry Smith lens[i] = 0; 94402834360SBarry Smith while (k < kend) { 94502834360SBarry Smith if (a->j[k++] >= first+ncols) break; 94602834360SBarry Smith lens[i]++; 94702834360SBarry Smith } 94802834360SBarry Smith } 94902834360SBarry Smith /* create submatrix */ 950cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 95108480c60SBarry Smith int n_cols,n_rows; 95208480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 95308480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 95408480c60SBarry Smith MatZeroEntries(*B); 95508480c60SBarry Smith C = *B; 95608480c60SBarry Smith } 95708480c60SBarry Smith else { 95802834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 95908480c60SBarry Smith } 960db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 961db02288aSLois Curfman McInnes 96202834360SBarry Smith /* loop over rows inserting into submatrix */ 963db02288aSLois Curfman McInnes a_new = c->a; 964db02288aSLois Curfman McInnes j_new = c->j; 965db02288aSLois Curfman McInnes i_new = c->i; 966db02288aSLois Curfman McInnes i_new[0] = -shift; 96702834360SBarry Smith for (i=0; i<nrows; i++) { 96802834360SBarry Smith for ( k=0; k<lens[i]; k++ ) { 969db02288aSLois Curfman McInnes *j_new++ = a->j[k+starts[i]] - first; 97002834360SBarry Smith } 971db02288aSLois Curfman McInnes PetscMemcpy(a_new,a->a + starts[i],lens[i]*sizeof(Scalar)); 972db02288aSLois Curfman McInnes a_new += lens[i]; 973db02288aSLois Curfman McInnes i_new[i+1] = i_new[i] + lens[i]; 9741987afe7SBarry Smith c->ilen[i] = lens[i]; 97502834360SBarry Smith } 976*0452661fSBarry Smith PetscFree(lens); 97702834360SBarry Smith } 97802834360SBarry Smith else { 97902834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 980*0452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 98102834360SBarry Smith ssmap = smap + shift; 982*0452661fSBarry Smith cwork = (int *) PetscMalloc((1+2*ncols)*sizeof(int)); CHKPTRQ(cwork); 98302834360SBarry Smith lens = cwork + ncols; 984*0452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 985cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 98617ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 98702834360SBarry Smith /* determine lens of each row */ 98802834360SBarry Smith for (i=0; i<nrows; i++) { 98902834360SBarry Smith kstart = a->i[irow[i]]+shift; 99002834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 99102834360SBarry Smith lens[i] = 0; 99202834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 99302834360SBarry Smith if (ssmap[a->j[k]]) { 99402834360SBarry Smith lens[i]++; 99502834360SBarry Smith } 99602834360SBarry Smith } 99702834360SBarry Smith } 99817ab2063SBarry Smith /* Create and fill new matrix */ 99908480c60SBarry Smith if (MatValidMatrix(*B)) { 100008480c60SBarry Smith int n_cols,n_rows; 100108480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 100208480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 100308480c60SBarry Smith MatZeroEntries(*B); 100408480c60SBarry Smith C = *B; 100508480c60SBarry Smith } 100608480c60SBarry Smith else { 100702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 100808480c60SBarry Smith } 100917ab2063SBarry Smith for (i=0; i<nrows; i++) { 101017ab2063SBarry Smith nznew = 0; 1011416022c9SBarry Smith kstart = a->i[irow[i]]+shift; 1012416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 101317ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 101402834360SBarry Smith if (ssmap[a->j[k]]) { 101502834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1016416022c9SBarry Smith vwork[nznew++] = a->a[k]; 101717ab2063SBarry Smith } 101817ab2063SBarry Smith } 1019416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 102017ab2063SBarry Smith } 102102834360SBarry Smith /* Free work space */ 102202834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 1023*0452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 102402834360SBarry Smith } 1025416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1026416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 102717ab2063SBarry Smith 102817ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1029416022c9SBarry Smith *B = C; 103017ab2063SBarry Smith return 0; 103117ab2063SBarry Smith } 103217ab2063SBarry Smith 1033a871dcd8SBarry Smith /* 103463b91edcSBarry Smith note: This can only work for identity for row and col. It would 103563b91edcSBarry Smith be good to check this and otherwise generate an error. 1036a871dcd8SBarry Smith */ 103763b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1038a871dcd8SBarry Smith { 103963b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 104008480c60SBarry Smith int ierr; 104163b91edcSBarry Smith Mat outA; 104263b91edcSBarry Smith 1043a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1044a871dcd8SBarry Smith 104563b91edcSBarry Smith outA = inA; 104663b91edcSBarry Smith inA->factor = FACTOR_LU; 104763b91edcSBarry Smith a->row = row; 104863b91edcSBarry Smith a->col = col; 104963b91edcSBarry Smith 1050*0452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 105163b91edcSBarry Smith 105208480c60SBarry Smith if (!a->diag) { 105308480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 105463b91edcSBarry Smith } 105563b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1056a871dcd8SBarry Smith return 0; 1057a871dcd8SBarry Smith } 1058a871dcd8SBarry Smith 1059cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1060cddf8d76SBarry Smith Mat **B) 1061cddf8d76SBarry Smith { 1062cddf8d76SBarry Smith int ierr,i; 1063cddf8d76SBarry Smith 1064cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1065*0452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1066cddf8d76SBarry Smith } 1067cddf8d76SBarry Smith 1068cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1069cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1070cddf8d76SBarry Smith } 1071cddf8d76SBarry Smith return 0; 1072cddf8d76SBarry Smith } 1073cddf8d76SBarry Smith 107417ab2063SBarry Smith /* -------------------------------------------------------------------*/ 107517ab2063SBarry Smith 107617ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 107717ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1078416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1079416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 108017ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 108117ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 108217ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 108317ab2063SBarry Smith MatRelax_SeqAIJ, 108417ab2063SBarry Smith MatTranspose_SeqAIJ, 108517ab2063SBarry Smith MatGetInfo_SeqAIJ,0, 108617ab2063SBarry Smith MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ, 108717ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 108817ab2063SBarry Smith MatCompress_SeqAIJ, 108917ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 109017ab2063SBarry Smith MatGetReordering_SeqAIJ, 109117ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 109217ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 109317ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 109417ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1095416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 1096a871dcd8SBarry Smith MatCopyPrivate_SeqAIJ,0,0, 1097cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 1098cddf8d76SBarry Smith MatGetSubMatrices_SeqAIJ}; 109917ab2063SBarry Smith 110017ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 110117ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 110217ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 110317ab2063SBarry Smith 110417ab2063SBarry Smith /*@C 110517ab2063SBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ format 110617ab2063SBarry Smith (the default uniprocessor PETSc format). 110717ab2063SBarry Smith 110817ab2063SBarry Smith Input Parameters: 110917ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 111017ab2063SBarry Smith . m - number of rows 111117ab2063SBarry Smith . n - number of columns 111217ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 111317ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 111417ab2063SBarry Smith 111517ab2063SBarry Smith Output Parameter: 1116416022c9SBarry Smith . A - the matrix 111717ab2063SBarry Smith 111817ab2063SBarry Smith Notes: 111917ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 112017ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 11210002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 11220002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 112317ab2063SBarry Smith 112417ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 112517ab2063SBarry Smith Set both nz and nnz to zero for PETSc to control dynamic memory 112617ab2063SBarry Smith allocation. 112717ab2063SBarry Smith 112817ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse 112917ab2063SBarry Smith 113017ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 113117ab2063SBarry Smith @*/ 1132416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 113317ab2063SBarry Smith { 1134416022c9SBarry Smith Mat B; 1135416022c9SBarry Smith Mat_SeqAIJ *b; 113617ab2063SBarry Smith int i,len,ierr; 1137416022c9SBarry Smith *A = 0; 1138*0452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1139416022c9SBarry Smith PLogObjectCreate(B); 1140*0452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1141416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1142416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1143416022c9SBarry Smith B->view = MatView_SeqAIJ; 1144416022c9SBarry Smith B->factor = 0; 1145416022c9SBarry Smith B->lupivotthreshold = 1.0; 1146416022c9SBarry Smith OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold); 1147416022c9SBarry Smith b->row = 0; 1148416022c9SBarry Smith b->col = 0; 1149416022c9SBarry Smith b->indexshift = 0; 1150416022c9SBarry Smith if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1; 115117ab2063SBarry Smith 1152416022c9SBarry Smith b->m = m; 1153416022c9SBarry Smith b->n = n; 1154*0452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 115517ab2063SBarry Smith if (!nnz) { 115617ab2063SBarry Smith if (nz <= 0) nz = 1; 1157416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 115817ab2063SBarry Smith nz = nz*m; 115917ab2063SBarry Smith } 116017ab2063SBarry Smith else { 116117ab2063SBarry Smith nz = 0; 1162416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 116317ab2063SBarry Smith } 116417ab2063SBarry Smith 116517ab2063SBarry Smith /* allocate the matrix space */ 1166416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 1167*0452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1168416022c9SBarry Smith b->j = (int *) (b->a + nz); 1169cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1170416022c9SBarry Smith b->i = b->j + nz; 1171416022c9SBarry Smith b->singlemalloc = 1; 117217ab2063SBarry Smith 1173416022c9SBarry Smith b->i[0] = -b->indexshift; 117417ab2063SBarry Smith for (i=1; i<m+1; i++) { 1175416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 117617ab2063SBarry Smith } 117717ab2063SBarry Smith 1178416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 1179*0452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1180416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1181416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 118217ab2063SBarry Smith 1183416022c9SBarry Smith b->nz = 0; 1184416022c9SBarry Smith b->maxnz = nz; 1185416022c9SBarry Smith b->sorted = 0; 1186416022c9SBarry Smith b->roworiented = 1; 1187416022c9SBarry Smith b->nonew = 0; 1188416022c9SBarry Smith b->diag = 0; 1189416022c9SBarry Smith b->assembled = 0; 1190416022c9SBarry Smith b->solve_work = 0; 119171bd300dSLois Curfman McInnes b->spptr = 0; 119217ab2063SBarry Smith 1193416022c9SBarry Smith *A = B; 119417ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_superlu")) { 1195416022c9SBarry Smith ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); 119617ab2063SBarry Smith } 119717ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_essl")) { 1198416022c9SBarry Smith ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); 119917ab2063SBarry Smith } 120017ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_dxml")) { 1201416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1202416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 120317ab2063SBarry Smith } 120417ab2063SBarry Smith 120517ab2063SBarry Smith return 0; 120617ab2063SBarry Smith } 120717ab2063SBarry Smith 120808480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues) 120917ab2063SBarry Smith { 1210416022c9SBarry Smith Mat C; 1211416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 121208480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 121317ab2063SBarry Smith 12144043dd9cSLois Curfman McInnes *B = 0; 1215416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix"); 1216*0452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1217416022c9SBarry Smith PLogObjectCreate(C); 1218*0452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 1219416022c9SBarry Smith PetscMemcpy(&C->ops,&MatOps,sizeof(struct _MatOps)); 1220416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1221416022c9SBarry Smith C->view = MatView_SeqAIJ; 1222416022c9SBarry Smith C->factor = A->factor; 1223416022c9SBarry Smith c->row = 0; 1224416022c9SBarry Smith c->col = 0; 1225416022c9SBarry Smith c->indexshift = shift; 122617ab2063SBarry Smith 1227416022c9SBarry Smith c->m = a->m; 1228416022c9SBarry Smith c->n = a->n; 122917ab2063SBarry Smith 1230*0452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 1231*0452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 123217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1233416022c9SBarry Smith c->imax[i] = a->imax[i]; 1234416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 123517ab2063SBarry Smith } 123617ab2063SBarry Smith 123717ab2063SBarry Smith /* allocate the matrix space */ 1238416022c9SBarry Smith c->singlemalloc = 1; 1239416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 1240*0452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1241416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1242416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1243416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 124417ab2063SBarry Smith if (m > 0) { 1245416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 124608480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1247416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 124817ab2063SBarry Smith } 124908480c60SBarry Smith } 125017ab2063SBarry Smith 1251416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1252416022c9SBarry Smith c->sorted = a->sorted; 1253416022c9SBarry Smith c->roworiented = a->roworiented; 1254416022c9SBarry Smith c->nonew = a->nonew; 125517ab2063SBarry Smith 1256416022c9SBarry Smith if (a->diag) { 1257*0452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1258416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 125917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1260416022c9SBarry Smith c->diag[i] = a->diag[i]; 126117ab2063SBarry Smith } 126217ab2063SBarry Smith } 1263416022c9SBarry Smith else c->diag = 0; 1264416022c9SBarry Smith c->assembled = 1; 1265416022c9SBarry Smith c->nz = a->nz; 1266416022c9SBarry Smith c->maxnz = a->maxnz; 1267416022c9SBarry Smith c->solve_work = 0; 12684043dd9cSLois Curfman McInnes c->spptr = 0; 1269416022c9SBarry Smith *B = C; 127017ab2063SBarry Smith return 0; 127117ab2063SBarry Smith } 127217ab2063SBarry Smith 1273416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 127417ab2063SBarry Smith { 1275416022c9SBarry Smith Mat_SeqAIJ *a; 1276416022c9SBarry Smith Mat B; 127717699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 127817ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 127917ab2063SBarry Smith MPI_Comm comm = vobj->comm; 128017ab2063SBarry Smith 128117699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 128217699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 128317ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1284416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 128548d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 128617ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 128717ab2063SBarry Smith 128817ab2063SBarry Smith /* read in row lengths */ 1289*0452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1290416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 129117ab2063SBarry Smith 129217ab2063SBarry Smith /* create our matrix */ 1293416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1294416022c9SBarry Smith B = *A; 1295416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1296416022c9SBarry Smith shift = a->indexshift; 129717ab2063SBarry Smith 129817ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1299416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 130017ab2063SBarry Smith if (shift) { 130117ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1302416022c9SBarry Smith a->j[i] += 1; 130317ab2063SBarry Smith } 130417ab2063SBarry Smith } 130517ab2063SBarry Smith 130617ab2063SBarry Smith /* read in nonzero values */ 1307416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 130817ab2063SBarry Smith 130917ab2063SBarry Smith /* set matrix "i" values */ 1310416022c9SBarry Smith a->i[0] = -shift; 131117ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1312416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1313416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 131417ab2063SBarry Smith } 1315*0452661fSBarry Smith PetscFree(rowlengths); 131617ab2063SBarry Smith 1317416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1318416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 131917ab2063SBarry Smith return 0; 132017ab2063SBarry Smith } 132117ab2063SBarry Smith 132217ab2063SBarry Smith 132317ab2063SBarry Smith 1324