117ab2063SBarry Smith #ifndef lint 2*b810aeb4SBarry Smith static char vcid[] = "$Id: aij.c,v 1.141 1996/01/26 19:33:31 bsmith Exp bsmith $"; 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" 12e4d965acSSatish Balay #include "petsc.h" 1317ab2063SBarry Smith 1417ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**); 1517ab2063SBarry Smith 16416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1717ab2063SBarry Smith { 18416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 19a2744918SBarry Smith int ierr, *ia, *ja,n,*idx,i; 203d54f168SSatish Balay /*Viewer V1, V2;*/ 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); 43f31bba2fSSatish Balay 440452661fSBarry Smith PetscFree(ia); PetscFree(ja); 4517ab2063SBarry Smith return 0; 4617ab2063SBarry Smith } 4717ab2063SBarry Smith 48227d817aSBarry Smith #define CHUNKSIZE 15 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; 554b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 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"); 694b0e389bSBarry Smith col = in[l] - shift; 704b0e389bSBarry Smith if (roworiented) { 714b0e389bSBarry Smith value = *v++; 724b0e389bSBarry Smith } 734b0e389bSBarry Smith else { 744b0e389bSBarry Smith value = v[k + l*m]; 754b0e389bSBarry Smith } 76416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 77416022c9SBarry Smith while (high-low > 5) { 78416022c9SBarry Smith t = (low+high)/2; 79416022c9SBarry Smith if (rp[t] > col) high = t; 80416022c9SBarry Smith else low = t; 8117ab2063SBarry Smith } 82416022c9SBarry Smith for ( i=low; i<high; i++ ) { 8317ab2063SBarry Smith if (rp[i] > col) break; 8417ab2063SBarry Smith if (rp[i] == col) { 85416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 8617ab2063SBarry Smith else ap[i] = value; 8717ab2063SBarry Smith goto noinsert; 8817ab2063SBarry Smith } 8917ab2063SBarry Smith } 9017ab2063SBarry Smith if (nonew) goto noinsert; 9117ab2063SBarry Smith if (nrow >= rmax) { 9217ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 93416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 9417ab2063SBarry Smith Scalar *new_a; 9517ab2063SBarry Smith 9617ab2063SBarry Smith /* malloc new storage space */ 97416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 980452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 9917ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 10017ab2063SBarry Smith new_i = new_j + new_nz; 10117ab2063SBarry Smith 10217ab2063SBarry Smith /* copy over old data into new slots */ 10317ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 104416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 105416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 106416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 107416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 10817ab2063SBarry Smith len*sizeof(int)); 109416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 110416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 11117ab2063SBarry Smith len*sizeof(Scalar)); 11217ab2063SBarry Smith /* free up old matrix storage */ 1130452661fSBarry Smith PetscFree(a->a); 1140452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 115416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 116416022c9SBarry Smith a->singlemalloc = 1; 11717ab2063SBarry Smith 11817ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 119416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 120416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 121416022c9SBarry Smith a->maxnz += CHUNKSIZE; 122*b810aeb4SBarry Smith a->reallocs++; 12317ab2063SBarry Smith } 124416022c9SBarry Smith N = nrow++ - 1; a->nz++; 125416022c9SBarry Smith /* shift up all the later entries in this row */ 126416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 12717ab2063SBarry Smith rp[ii+1] = rp[ii]; 12817ab2063SBarry Smith ap[ii+1] = ap[ii]; 12917ab2063SBarry Smith } 13017ab2063SBarry Smith rp[i] = col; 13117ab2063SBarry Smith ap[i] = value; 13217ab2063SBarry Smith noinsert:; 133416022c9SBarry Smith low = i + 1; 13417ab2063SBarry Smith } 13517ab2063SBarry Smith ailen[row] = nrow; 13617ab2063SBarry Smith } 13717ab2063SBarry Smith return 0; 13817ab2063SBarry Smith } 13917ab2063SBarry Smith 1407eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 1417eb43aa7SLois Curfman McInnes { 1427eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 143b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 1447eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 1457eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 1467eb43aa7SLois Curfman McInnes 1477eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 1487eb43aa7SLois Curfman McInnes row = im[k]; 1497eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 1507eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 1517eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 1527eb43aa7SLois Curfman McInnes nrow = ailen[row]; 1537eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 1547eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 1557eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 1567eb43aa7SLois Curfman McInnes col = in[l] - shift; 1577eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 1587eb43aa7SLois Curfman McInnes while (high-low > 5) { 1597eb43aa7SLois Curfman McInnes t = (low+high)/2; 1607eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 1617eb43aa7SLois Curfman McInnes else low = t; 1627eb43aa7SLois Curfman McInnes } 1637eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 1647eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 1657eb43aa7SLois Curfman McInnes if (rp[i] == col) { 166b49de8d1SLois Curfman McInnes *v++ = ap[i]; 1677eb43aa7SLois Curfman McInnes goto finished; 1687eb43aa7SLois Curfman McInnes } 1697eb43aa7SLois Curfman McInnes } 170b49de8d1SLois Curfman McInnes *v++ = zero; 1717eb43aa7SLois Curfman McInnes finished:; 1727eb43aa7SLois Curfman McInnes } 1737eb43aa7SLois Curfman McInnes } 1747eb43aa7SLois Curfman McInnes return 0; 1757eb43aa7SLois Curfman McInnes } 1767eb43aa7SLois Curfman McInnes 17717ab2063SBarry Smith #include "draw.h" 17817ab2063SBarry Smith #include "pinclude/pviewer.h" 179416022c9SBarry Smith #include "sysio.h" 18017ab2063SBarry Smith 181416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 18217ab2063SBarry Smith { 183416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 184416022c9SBarry Smith int i, fd, *col_lens, ierr; 18517ab2063SBarry Smith 186416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1870452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 188416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 189416022c9SBarry Smith col_lens[1] = a->m; 190416022c9SBarry Smith col_lens[2] = a->n; 191416022c9SBarry Smith col_lens[3] = a->nz; 192416022c9SBarry Smith 193416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 194416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 195416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 19617ab2063SBarry Smith } 197416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 1980452661fSBarry Smith PetscFree(col_lens); 199416022c9SBarry Smith 200416022c9SBarry Smith /* store column indices (zero start index) */ 201416022c9SBarry Smith if (a->indexshift) { 202416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 20317ab2063SBarry Smith } 204416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 205416022c9SBarry Smith if (a->indexshift) { 206416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 20717ab2063SBarry Smith } 208416022c9SBarry Smith 209416022c9SBarry Smith /* store nonzero values */ 210416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 21117ab2063SBarry Smith return 0; 21217ab2063SBarry Smith } 213416022c9SBarry Smith 214416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 215416022c9SBarry Smith { 216416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 217416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 21817ab2063SBarry Smith FILE *fd; 21917ab2063SBarry Smith char *outputname; 22017ab2063SBarry Smith 221227d817aSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 222416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 223416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 22417ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 22508480c60SBarry Smith /* no need to print additional information */ ; 22617ab2063SBarry Smith } 22717ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 22817ab2063SBarry Smith int nz, nzalloc, mem; 229416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 230416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 23117ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 23217ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 23317ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 23417ab2063SBarry Smith 23517ab2063SBarry Smith for (i=0; i<m; i++) { 236416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 23717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 238416022c9SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j],real(a->a[j]), 239416022c9SBarry Smith imag(a->a[j])); 24017ab2063SBarry Smith #else 241416022c9SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j], a->a[j]); 24217ab2063SBarry Smith #endif 24317ab2063SBarry Smith } 24417ab2063SBarry Smith } 24517ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 24617ab2063SBarry Smith } 24717ab2063SBarry Smith else { 24817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 24917ab2063SBarry Smith fprintf(fd,"row %d:",i); 250416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 25117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 252416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 253416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 25417ab2063SBarry Smith } 25517ab2063SBarry Smith else { 256416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 25717ab2063SBarry Smith } 25817ab2063SBarry Smith #else 259416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 26017ab2063SBarry Smith #endif 26117ab2063SBarry Smith } 26217ab2063SBarry Smith fprintf(fd,"\n"); 26317ab2063SBarry Smith } 26417ab2063SBarry Smith } 26517ab2063SBarry Smith fflush(fd); 266416022c9SBarry Smith return 0; 267416022c9SBarry Smith } 268416022c9SBarry Smith 269416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 270416022c9SBarry Smith { 271416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 272cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 273cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 274d7e8b826SBarry Smith Draw draw = (Draw) viewer; 275cddf8d76SBarry Smith DrawButton button; 276cddf8d76SBarry Smith 277416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 278416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 279416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 280416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 281cddf8d76SBarry Smith color = DRAW_BLUE; 282416022c9SBarry Smith for ( i=0; i<m; i++ ) { 283cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 284416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 285cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 286cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 287cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 288cddf8d76SBarry Smith #else 289cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 290cddf8d76SBarry Smith #endif 291cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 292cddf8d76SBarry Smith } 293cddf8d76SBarry Smith } 294cddf8d76SBarry Smith color = DRAW_CYAN; 295cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 296cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 297cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 298cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 299cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 300cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 301cddf8d76SBarry Smith } 302cddf8d76SBarry Smith } 303cddf8d76SBarry Smith color = DRAW_RED; 304cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 305cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 306cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 307cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 308cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 309cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 310cddf8d76SBarry Smith #else 311cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 312cddf8d76SBarry Smith #endif 313cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 314416022c9SBarry Smith } 315416022c9SBarry Smith } 316416022c9SBarry Smith DrawFlush(draw); 317cddf8d76SBarry Smith DrawGetPause(draw,&pause); 318cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 319cddf8d76SBarry Smith 320cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 321cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 322cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 323cddf8d76SBarry Smith DrawClear(draw); 324cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 325cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 326cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 327cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 328cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 329cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 330cddf8d76SBarry Smith w *= scale; h *= scale; 331cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 332cddf8d76SBarry Smith color = DRAW_BLUE; 333cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 334cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 335cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 336cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 337cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 338cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 339cddf8d76SBarry Smith #else 340cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 341cddf8d76SBarry Smith #endif 342cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 343cddf8d76SBarry Smith } 344cddf8d76SBarry Smith } 345cddf8d76SBarry Smith color = DRAW_CYAN; 346cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 347cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 348cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 349cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 350cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 351cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 352cddf8d76SBarry Smith } 353cddf8d76SBarry Smith } 354cddf8d76SBarry Smith color = DRAW_RED; 355cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 356cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 357cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 358cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 359cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 360cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 361cddf8d76SBarry Smith #else 362cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 363cddf8d76SBarry Smith #endif 364cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 365cddf8d76SBarry Smith } 366cddf8d76SBarry Smith } 367cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 368cddf8d76SBarry Smith } 369416022c9SBarry Smith return 0; 370416022c9SBarry Smith } 371416022c9SBarry Smith 372416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 373416022c9SBarry Smith { 374416022c9SBarry Smith Mat A = (Mat) obj; 375416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 376416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 377416022c9SBarry Smith 378416022c9SBarry Smith if (!viewer) { 379416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 380416022c9SBarry Smith } 381416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 382416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 383416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 384416022c9SBarry Smith } 385416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 386416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 387416022c9SBarry Smith } 388416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 389416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 390416022c9SBarry Smith } 391416022c9SBarry Smith } 392416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 393416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 394416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 39517ab2063SBarry Smith } 39617ab2063SBarry Smith return 0; 39717ab2063SBarry Smith } 398c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 399416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 40017ab2063SBarry Smith { 401416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 40241c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 403416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 404416022c9SBarry Smith Scalar *aa = a->a, *ap; 40517ab2063SBarry Smith 40617ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 40717ab2063SBarry Smith 40817ab2063SBarry Smith for ( i=1; i<m; i++ ) { 409416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 41017ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 41117ab2063SBarry Smith if (fshift) { 412416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 41317ab2063SBarry Smith N = ailen[i]; 41417ab2063SBarry Smith for ( j=0; j<N; j++ ) { 41517ab2063SBarry Smith ip[j-fshift] = ip[j]; 41617ab2063SBarry Smith ap[j-fshift] = ap[j]; 41717ab2063SBarry Smith } 41817ab2063SBarry Smith } 41917ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 42017ab2063SBarry Smith } 42117ab2063SBarry Smith if (m) { 42217ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 42317ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 42417ab2063SBarry Smith } 42517ab2063SBarry Smith /* reset ilen and imax for each row */ 42617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 42717ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 42817ab2063SBarry Smith } 429416022c9SBarry Smith a->nz = ai[m] + shift; 43017ab2063SBarry Smith 43117ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 432416022c9SBarry Smith if (fshift && a->diag) { 4330452661fSBarry Smith PetscFree(a->diag); 434416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 435416022c9SBarry Smith a->diag = 0; 43617ab2063SBarry Smith } 437*b810aeb4SBarry Smith PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n", 438*b810aeb4SBarry Smith fshift,a->nz,m); 439*b810aeb4SBarry Smith PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n", 440*b810aeb4SBarry Smith a->reallocs); 44176dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 44241c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 44317ab2063SBarry Smith return 0; 44417ab2063SBarry Smith } 44517ab2063SBarry Smith 446416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 44717ab2063SBarry Smith { 448416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 449cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 45017ab2063SBarry Smith return 0; 45117ab2063SBarry Smith } 452416022c9SBarry Smith 45317ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 45417ab2063SBarry Smith { 455416022c9SBarry Smith Mat A = (Mat) obj; 456416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 457d5d45c9bSBarry Smith 45817ab2063SBarry Smith #if defined(PETSC_LOG) 459416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 46017ab2063SBarry Smith #endif 4610452661fSBarry Smith PetscFree(a->a); 4620452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4630452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4640452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4650452661fSBarry Smith if (a->imax) PetscFree(a->imax); 4660452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 46776dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 4680452661fSBarry Smith PetscFree(a); 469416022c9SBarry Smith PLogObjectDestroy(A); 4700452661fSBarry Smith PetscHeaderDestroy(A); 47117ab2063SBarry Smith return 0; 47217ab2063SBarry Smith } 47317ab2063SBarry Smith 474416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 47517ab2063SBarry Smith { 47617ab2063SBarry Smith return 0; 47717ab2063SBarry Smith } 47817ab2063SBarry Smith 479416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 48017ab2063SBarry Smith { 481416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 482416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 4834b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 484416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 485416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 486416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 487e2f28af5SBarry Smith else if (op == ROWS_SORTED || 488e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 489e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 490e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 491df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 492df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4934d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 4946c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_1) a->inode.limit = 1; 4956c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_2) a->inode.limit = 2; 4966c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_3) a->inode.limit = 3; 4976c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_4) a->inode.limit = 4; 4986c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_5) a->inode.limit = 5; 499e2f28af5SBarry Smith else 5004d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 50117ab2063SBarry Smith return 0; 50217ab2063SBarry Smith } 50317ab2063SBarry Smith 504416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 50517ab2063SBarry Smith { 506416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 507416022c9SBarry Smith int i,j, n,shift = a->indexshift; 50817ab2063SBarry Smith Scalar *x, zero = 0.0; 50917ab2063SBarry Smith 51017ab2063SBarry Smith VecSet(&zero,v); 51117ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 512416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 513416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 514416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 515416022c9SBarry Smith if (a->j[j]+shift == i) { 516416022c9SBarry Smith x[i] = a->a[j]; 51717ab2063SBarry Smith break; 51817ab2063SBarry Smith } 51917ab2063SBarry Smith } 52017ab2063SBarry Smith } 52117ab2063SBarry Smith return 0; 52217ab2063SBarry Smith } 52317ab2063SBarry Smith 52417ab2063SBarry Smith /* -------------------------------------------------------*/ 52517ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 52617ab2063SBarry Smith /* -------------------------------------------------------*/ 527416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 52817ab2063SBarry Smith { 529416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 53017ab2063SBarry Smith Scalar *x, *y, *v, alpha; 531416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 53217ab2063SBarry Smith 53317ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 534cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 53517ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 53617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 537416022c9SBarry Smith idx = a->j + a->i[i] + shift; 538416022c9SBarry Smith v = a->a + a->i[i] + shift; 539416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 54017ab2063SBarry Smith alpha = x[i]; 54117ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 54217ab2063SBarry Smith } 543416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 54417ab2063SBarry Smith return 0; 54517ab2063SBarry Smith } 546d5d45c9bSBarry Smith 547416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 54817ab2063SBarry Smith { 549416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 55017ab2063SBarry Smith Scalar *x, *y, *v, alpha; 551416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 55217ab2063SBarry Smith 55317ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 55417ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 55517ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 55617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 557416022c9SBarry Smith idx = a->j + a->i[i] + shift; 558416022c9SBarry Smith v = a->a + a->i[i] + shift; 559416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 56017ab2063SBarry Smith alpha = x[i]; 56117ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 56217ab2063SBarry Smith } 56317ab2063SBarry Smith return 0; 56417ab2063SBarry Smith } 56517ab2063SBarry Smith 566416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 56717ab2063SBarry Smith { 568416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 56917ab2063SBarry Smith Scalar *x, *y, *v, sum; 570416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 57117ab2063SBarry Smith 57217ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 57317ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 574416022c9SBarry Smith idx = a->j; 575416022c9SBarry Smith v = a->a; 576416022c9SBarry Smith ii = a->i; 57717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 578416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 57917ab2063SBarry Smith sum = 0.0; 58017ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 58117ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 582416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 58317ab2063SBarry Smith y[i] = sum; 58417ab2063SBarry Smith } 585416022c9SBarry Smith PLogFlops(2*a->nz - m); 58617ab2063SBarry Smith return 0; 58717ab2063SBarry Smith } 58817ab2063SBarry Smith 589416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 59017ab2063SBarry Smith { 591416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 59217ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 593cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 59417ab2063SBarry Smith 59517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 59617ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 597cddf8d76SBarry Smith idx = a->j; 598cddf8d76SBarry Smith v = a->a; 599cddf8d76SBarry Smith ii = a->i; 60017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 601cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 60217ab2063SBarry Smith sum = y[i]; 603cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 604cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 60517ab2063SBarry Smith z[i] = sum; 60617ab2063SBarry Smith } 607416022c9SBarry Smith PLogFlops(2*a->nz); 60817ab2063SBarry Smith return 0; 60917ab2063SBarry Smith } 61017ab2063SBarry Smith 61117ab2063SBarry Smith /* 61217ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 61317ab2063SBarry Smith */ 61417ab2063SBarry Smith 615416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 61617ab2063SBarry Smith { 617416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 618416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 61917ab2063SBarry Smith 6200452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 621416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 622416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 623416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 624416022c9SBarry Smith if (a->j[j]+shift == i) { 62517ab2063SBarry Smith diag[i] = j - shift; 62617ab2063SBarry Smith break; 62717ab2063SBarry Smith } 62817ab2063SBarry Smith } 62917ab2063SBarry Smith } 630416022c9SBarry Smith a->diag = diag; 63117ab2063SBarry Smith return 0; 63217ab2063SBarry Smith } 63317ab2063SBarry Smith 634416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 63517ab2063SBarry Smith double fshift,int its,Vec xx) 63617ab2063SBarry Smith { 637416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 638416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 639d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 64017ab2063SBarry Smith 64117ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 642416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 643416022c9SBarry Smith diag = a->diag; 644416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 64517ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 64617ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 64717ab2063SBarry Smith bs = b + shift; 64817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 649416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 650416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 651416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 652416022c9SBarry Smith v = a->a + diag[i] + (!shift); 65317ab2063SBarry Smith sum = b[i]*d/omega; 65417ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 65517ab2063SBarry Smith x[i] = sum; 65617ab2063SBarry Smith } 65717ab2063SBarry Smith return 0; 65817ab2063SBarry Smith } 65917ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 66017ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 66117ab2063SBarry Smith } 662416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 66317ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 66417ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 66517ab2063SBarry Smith 66617ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 66717ab2063SBarry Smith 66817ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 66917ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 67017ab2063SBarry Smith is the relaxation factor. 67117ab2063SBarry Smith */ 6720452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 67317ab2063SBarry Smith scale = (2.0/omega) - 1.0; 67417ab2063SBarry Smith 67517ab2063SBarry Smith /* x = (E + U)^{-1} b */ 67617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 677416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 678416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 679416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 680416022c9SBarry Smith v = a->a + diag[i] + (!shift); 68117ab2063SBarry Smith sum = b[i]; 68217ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 68317ab2063SBarry Smith x[i] = omega*(sum/d); 68417ab2063SBarry Smith } 68517ab2063SBarry Smith 68617ab2063SBarry Smith /* t = b - (2*E - D)x */ 687416022c9SBarry Smith v = a->a; 68817ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 68917ab2063SBarry Smith 69017ab2063SBarry Smith /* t = (E + L)^{-1}t */ 691416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 692416022c9SBarry Smith diag = a->diag; 69317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 694416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 695416022c9SBarry Smith n = diag[i] - a->i[i]; 696416022c9SBarry Smith idx = a->j + a->i[i] + shift; 697416022c9SBarry Smith v = a->a + a->i[i] + shift; 69817ab2063SBarry Smith sum = t[i]; 69917ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 70017ab2063SBarry Smith t[i] = omega*(sum/d); 70117ab2063SBarry Smith } 70217ab2063SBarry Smith 70317ab2063SBarry Smith /* x = x + t */ 70417ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7050452661fSBarry Smith PetscFree(t); 70617ab2063SBarry Smith return 0; 70717ab2063SBarry Smith } 70817ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 70917ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 71017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 711416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 712416022c9SBarry Smith n = diag[i] - 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] = omega*(sum/d); 71817ab2063SBarry Smith } 71917ab2063SBarry Smith xb = x; 72017ab2063SBarry Smith } 72117ab2063SBarry Smith else xb = b; 72217ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 72317ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 72417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 725416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 72617ab2063SBarry Smith } 72717ab2063SBarry Smith } 72817ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 72917ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 730416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 731416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 732416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 733416022c9SBarry Smith v = a->a + diag[i] + (!shift); 73417ab2063SBarry Smith sum = xb[i]; 73517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 73617ab2063SBarry Smith x[i] = omega*(sum/d); 73717ab2063SBarry Smith } 73817ab2063SBarry Smith } 73917ab2063SBarry Smith its--; 74017ab2063SBarry Smith } 74117ab2063SBarry Smith while (its--) { 74217ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 74317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 744416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 745416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 746416022c9SBarry Smith idx = a->j + a->i[i] + shift; 747416022c9SBarry Smith v = a->a + a->i[i] + shift; 74817ab2063SBarry Smith sum = b[i]; 74917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 75017ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 75117ab2063SBarry Smith } 75217ab2063SBarry Smith } 75317ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 75417ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 755416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 756416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 757416022c9SBarry Smith idx = a->j + a->i[i] + shift; 758416022c9SBarry Smith v = a->a + a->i[i] + shift; 75917ab2063SBarry Smith sum = b[i]; 76017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 76117ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 76217ab2063SBarry Smith } 76317ab2063SBarry Smith } 76417ab2063SBarry Smith } 76517ab2063SBarry Smith return 0; 76617ab2063SBarry Smith } 76717ab2063SBarry Smith 768d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 76917ab2063SBarry Smith { 770416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 771416022c9SBarry Smith *nz = a->nz; 772416022c9SBarry Smith *nzalloc = a->maxnz; 773416022c9SBarry Smith *mem = (int)A->mem; 77417ab2063SBarry Smith return 0; 77517ab2063SBarry Smith } 77617ab2063SBarry Smith 77717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 77817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 77917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 78017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 78117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 78317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78417ab2063SBarry Smith 78517ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 78617ab2063SBarry Smith { 787416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 788416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 78917ab2063SBarry Smith 79017ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 79117ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 79217ab2063SBarry Smith if (diag) { 79317ab2063SBarry Smith for ( i=0; i<N; i++ ) { 794416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 795416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 796416022c9SBarry Smith a->ilen[rows[i]] = 1; 797416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 798416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 79917ab2063SBarry Smith } 80017ab2063SBarry Smith else { 80117ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 80217ab2063SBarry Smith CHKERRQ(ierr); 80317ab2063SBarry Smith } 80417ab2063SBarry Smith } 80517ab2063SBarry Smith } 80617ab2063SBarry Smith else { 80717ab2063SBarry Smith for ( i=0; i<N; i++ ) { 808416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 809416022c9SBarry Smith a->ilen[rows[i]] = 0; 81017ab2063SBarry Smith } 81117ab2063SBarry Smith } 81217ab2063SBarry Smith ISRestoreIndices(is,&rows); 81317ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 81417ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 81517ab2063SBarry Smith return 0; 81617ab2063SBarry Smith } 81717ab2063SBarry Smith 818416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 81917ab2063SBarry Smith { 820416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 821416022c9SBarry Smith *m = a->m; *n = a->n; 82217ab2063SBarry Smith return 0; 82317ab2063SBarry Smith } 82417ab2063SBarry Smith 825416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 82617ab2063SBarry Smith { 827416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 828416022c9SBarry Smith *m = 0; *n = a->m; 82917ab2063SBarry Smith return 0; 83017ab2063SBarry Smith } 831416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 83217ab2063SBarry Smith { 833416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 834c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 83517ab2063SBarry Smith 836416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 83717ab2063SBarry Smith 838416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 839416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 84017ab2063SBarry Smith if (idx) { 84117ab2063SBarry Smith if (*nz) { 842416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8430452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 84417ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 84517ab2063SBarry Smith } 84617ab2063SBarry Smith else *idx = 0; 84717ab2063SBarry Smith } 84817ab2063SBarry Smith return 0; 84917ab2063SBarry Smith } 85017ab2063SBarry Smith 851416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 85217ab2063SBarry Smith { 8530452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 85417ab2063SBarry Smith return 0; 85517ab2063SBarry Smith } 85617ab2063SBarry Smith 857cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 85817ab2063SBarry Smith { 859416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 860416022c9SBarry Smith Scalar *v = a->a; 86117ab2063SBarry Smith double sum = 0.0; 862416022c9SBarry Smith int i, j,shift = a->indexshift; 86317ab2063SBarry Smith 86417ab2063SBarry Smith if (type == NORM_FROBENIUS) { 865416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 86617ab2063SBarry Smith #if defined(PETSC_COMPLEX) 86717ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 86817ab2063SBarry Smith #else 86917ab2063SBarry Smith sum += (*v)*(*v); v++; 87017ab2063SBarry Smith #endif 87117ab2063SBarry Smith } 87217ab2063SBarry Smith *norm = sqrt(sum); 87317ab2063SBarry Smith } 87417ab2063SBarry Smith else if (type == NORM_1) { 87517ab2063SBarry Smith double *tmp; 876416022c9SBarry Smith int *jj = a->j; 8770452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 878cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 87917ab2063SBarry Smith *norm = 0.0; 880416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 881a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 88217ab2063SBarry Smith } 883416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 88417ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 88517ab2063SBarry Smith } 8860452661fSBarry Smith PetscFree(tmp); 88717ab2063SBarry Smith } 88817ab2063SBarry Smith else if (type == NORM_INFINITY) { 88917ab2063SBarry Smith *norm = 0.0; 890416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 891416022c9SBarry Smith v = a->a + a->i[j] + shift; 89217ab2063SBarry Smith sum = 0.0; 893416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 894cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 89517ab2063SBarry Smith } 89617ab2063SBarry Smith if (sum > *norm) *norm = sum; 89717ab2063SBarry Smith } 89817ab2063SBarry Smith } 89917ab2063SBarry Smith else { 90048d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 90117ab2063SBarry Smith } 90217ab2063SBarry Smith return 0; 90317ab2063SBarry Smith } 90417ab2063SBarry Smith 905416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 90617ab2063SBarry Smith { 907416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 908416022c9SBarry Smith Mat C; 909416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 910416022c9SBarry Smith int shift = a->indexshift; 911d5d45c9bSBarry Smith Scalar *array = a->a; 91217ab2063SBarry Smith 913416022c9SBarry Smith if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 9140452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 915cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 91617ab2063SBarry Smith if (shift) { 91717ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 91817ab2063SBarry Smith } 91917ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 920416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 9210452661fSBarry Smith PetscFree(col); 92217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 92317ab2063SBarry Smith len = ai[i+1]-ai[i]; 924416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 92517ab2063SBarry Smith array += len; aj += len; 92617ab2063SBarry Smith } 92717ab2063SBarry Smith if (shift) { 92817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 92917ab2063SBarry Smith } 93017ab2063SBarry Smith 931416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 932416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 93317ab2063SBarry Smith 934416022c9SBarry Smith if (B) { 935416022c9SBarry Smith *B = C; 93617ab2063SBarry Smith } else { 937416022c9SBarry Smith /* This isn't really an in-place transpose */ 9380452661fSBarry Smith PetscFree(a->a); 9390452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9400452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9410452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9420452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9430452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9440452661fSBarry Smith PetscFree(a); 945416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9460452661fSBarry Smith PetscHeaderDestroy(C); 94717ab2063SBarry Smith } 94817ab2063SBarry Smith return 0; 94917ab2063SBarry Smith } 95017ab2063SBarry Smith 951f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 95217ab2063SBarry Smith { 953416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 95417ab2063SBarry Smith Scalar *l,*r,x,*v; 955d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 95617ab2063SBarry Smith 95717ab2063SBarry Smith if (ll) { 95817ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 959f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 960416022c9SBarry Smith v = a->a; 96117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 96217ab2063SBarry Smith x = l[i]; 963416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 96417ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 96517ab2063SBarry Smith } 96617ab2063SBarry Smith } 96717ab2063SBarry Smith if (rr) { 96817ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 969f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 970416022c9SBarry Smith v = a->a; jj = a->j; 97117ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 97217ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 97317ab2063SBarry Smith } 97417ab2063SBarry Smith } 97517ab2063SBarry Smith return 0; 97617ab2063SBarry Smith } 97717ab2063SBarry Smith 978cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 97917ab2063SBarry Smith { 980db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 98102834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 982a2744918SBarry Smith register int sum,lensi; 98302834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 984a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 985db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 986416022c9SBarry Smith Mat C; 98717ab2063SBarry Smith 98817ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 98917ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 99017ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 99117ab2063SBarry Smith 99202834360SBarry Smith if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 99302834360SBarry Smith /* special case of contiguous rows */ 9940452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 99502834360SBarry Smith starts = lens + ncols; 99602834360SBarry Smith /* loop over new rows determining lens and starting points */ 99702834360SBarry Smith for (i=0; i<nrows; i++) { 998a2744918SBarry Smith kstart = ai[irow[i]]+shift; 999a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 100002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1001d8ced48eSBarry Smith if (aj[k]+shift >= first) { 100202834360SBarry Smith starts[i] = k; 100302834360SBarry Smith break; 100402834360SBarry Smith } 100502834360SBarry Smith } 1006a2744918SBarry Smith sum = 0; 100702834360SBarry Smith while (k < kend) { 1008d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1009a2744918SBarry Smith sum++; 101002834360SBarry Smith } 1011a2744918SBarry Smith lens[i] = sum; 101202834360SBarry Smith } 101302834360SBarry Smith /* create submatrix */ 1014cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 101508480c60SBarry Smith int n_cols,n_rows; 101608480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 101708480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1018d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 101908480c60SBarry Smith C = *B; 102008480c60SBarry Smith } 102108480c60SBarry Smith else { 102202834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 102308480c60SBarry Smith } 1024db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1025db02288aSLois Curfman McInnes 102602834360SBarry Smith /* loop over rows inserting into submatrix */ 1027db02288aSLois Curfman McInnes a_new = c->a; 1028db02288aSLois Curfman McInnes j_new = c->j; 1029db02288aSLois Curfman McInnes i_new = c->i; 1030db02288aSLois Curfman McInnes i_new[0] = -shift; 103102834360SBarry Smith for (i=0; i<nrows; i++) { 1032a2744918SBarry Smith ii = starts[i]; 1033a2744918SBarry Smith lensi = lens[i]; 1034a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1035a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 103602834360SBarry Smith } 1037a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1038a2744918SBarry Smith a_new += lensi; 1039a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1040a2744918SBarry Smith c->ilen[i] = lensi; 104102834360SBarry Smith } 10420452661fSBarry Smith PetscFree(lens); 104302834360SBarry Smith } 104402834360SBarry Smith else { 104502834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 10460452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 104702834360SBarry Smith ssmap = smap + shift; 10487eb43aa7SLois Curfman McInnes cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 104902834360SBarry Smith lens = cwork + ncols; 10500452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1051cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 105217ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 105302834360SBarry Smith /* determine lens of each row */ 105402834360SBarry Smith for (i=0; i<nrows; i++) { 1055d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 105602834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 105702834360SBarry Smith lens[i] = 0; 105802834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1059d8ced48eSBarry Smith if (ssmap[aj[k]]) { 106002834360SBarry Smith lens[i]++; 106102834360SBarry Smith } 106202834360SBarry Smith } 106302834360SBarry Smith } 106417ab2063SBarry Smith /* Create and fill new matrix */ 1065a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 106608480c60SBarry Smith int n_cols,n_rows; 106708480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 106808480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1069d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 107008480c60SBarry Smith C = *B; 107108480c60SBarry Smith } 107208480c60SBarry Smith else { 107302834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 107408480c60SBarry Smith } 107517ab2063SBarry Smith for (i=0; i<nrows; i++) { 107617ab2063SBarry Smith nznew = 0; 1077d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 1078416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 107917ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 108002834360SBarry Smith if (ssmap[a->j[k]]) { 108102834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1082416022c9SBarry Smith vwork[nznew++] = a->a[k]; 108317ab2063SBarry Smith } 108417ab2063SBarry Smith } 1085416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 108617ab2063SBarry Smith } 108702834360SBarry Smith /* Free work space */ 108802834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 10890452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 109002834360SBarry Smith } 1091416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1092416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 109317ab2063SBarry Smith 109417ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1095416022c9SBarry Smith *B = C; 109617ab2063SBarry Smith return 0; 109717ab2063SBarry Smith } 109817ab2063SBarry Smith 1099a871dcd8SBarry Smith /* 110063b91edcSBarry Smith note: This can only work for identity for row and col. It would 110163b91edcSBarry Smith be good to check this and otherwise generate an error. 1102a871dcd8SBarry Smith */ 110363b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1104a871dcd8SBarry Smith { 110563b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 110608480c60SBarry Smith int ierr; 110763b91edcSBarry Smith Mat outA; 110863b91edcSBarry Smith 1109a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1110a871dcd8SBarry Smith 111163b91edcSBarry Smith outA = inA; 111263b91edcSBarry Smith inA->factor = FACTOR_LU; 111363b91edcSBarry Smith a->row = row; 111463b91edcSBarry Smith a->col = col; 111563b91edcSBarry Smith 11160452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 111763b91edcSBarry Smith 111808480c60SBarry Smith if (!a->diag) { 111908480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 112063b91edcSBarry Smith } 112163b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1122a871dcd8SBarry Smith return 0; 1123a871dcd8SBarry Smith } 1124a871dcd8SBarry Smith 1125f0b747eeSBarry Smith #include "pinclude/plapack.h" 1126f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1127f0b747eeSBarry Smith { 1128f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1129f0b747eeSBarry Smith int one = 1; 1130f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1131f0b747eeSBarry Smith PLogFlops(a->nz); 1132f0b747eeSBarry Smith return 0; 1133f0b747eeSBarry Smith } 1134f0b747eeSBarry Smith 1135cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1136cddf8d76SBarry Smith Mat **B) 1137cddf8d76SBarry Smith { 1138cddf8d76SBarry Smith int ierr,i; 1139cddf8d76SBarry Smith 1140cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 11410452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1142cddf8d76SBarry Smith } 1143cddf8d76SBarry Smith 1144cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1145cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1146cddf8d76SBarry Smith } 1147cddf8d76SBarry Smith return 0; 1148cddf8d76SBarry Smith } 1149cddf8d76SBarry Smith 1150e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 11514dcbc457SBarry Smith { 1152e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11538a047759SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *table, *nidx, isz, val; 11548a047759SSatish Balay int start, end, *ai, *aj; 11554dcbc457SBarry Smith 11568a047759SSatish Balay shift = a->indexshift; 1157e4d965acSSatish Balay m = a->m; 1158e4d965acSSatish Balay ai = a->i; 11598a047759SSatish Balay aj = a->j+shift; 11608a047759SSatish Balay 116104a348a9SBarry Smith table = (int *) PetscMalloc(2*m*sizeof(int)); CHKPTRQ(table); nidx = table + m; 11628a047759SSatish Balay 1163e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 1164e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1165e4d965acSSatish Balay /* Initialise the two local arrays */ 1166e4d965acSSatish Balay isz = 0; 1167e4d965acSSatish Balay PetscMemzero(table,m*sizeof(int)); 1168e4d965acSSatish Balay 1169e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 11704dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1171e4d965acSSatish Balay ierr = ISGetLocalSize(is[i],&n); CHKERRQ(ierr); 1172e4d965acSSatish Balay 1173e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1174e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 1175e4d965acSSatish Balay if(!table[idx[j]]++) { nidx[isz++] = idx[j];} 11764dcbc457SBarry Smith } 1177e4d965acSSatish Balay 117804a348a9SBarry Smith k = 0; 117904a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 118004a348a9SBarry Smith n = isz; 118104a348a9SBarry Smith for ( ; k<n ; k++){ /* do only those rows in nidx[k] not yet done */ 1182e4d965acSSatish Balay row = nidx[k]; 1183e4d965acSSatish Balay start = ai[row]; 1184e4d965acSSatish Balay end = ai[row+1]; 118504a348a9SBarry Smith for ( l = start; l<end ; l++){ 11868a047759SSatish Balay val = aj[l] + shift; 1187e4d965acSSatish Balay if (!table[val]++) {nidx[isz++] = val;} 1188e4d965acSSatish Balay } 1189e4d965acSSatish Balay } 1190e4d965acSSatish Balay } 1191e4d965acSSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 1192e4d965acSSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1193e4d965acSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1194e4d965acSSatish Balay } 119504a348a9SBarry Smith PetscFree(table); 1196e4d965acSSatish Balay return 0; 11974dcbc457SBarry Smith } 119817ab2063SBarry Smith 1199682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1200682d7d0cSBarry Smith { 1201682d7d0cSBarry Smith static int called = 0; 1202682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1203682d7d0cSBarry Smith 1204682d7d0cSBarry Smith if (called) return 0; else called = 1; 12052a7368beSLois Curfman McInnes MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 1206682d7d0cSBarry Smith MPIU_printf(comm," -mat_lu_pivotthreshold <threshold>\n"); 12072a7368beSLois Curfman McInnes MPIU_printf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 1208682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_no_inode - Do not use inodes\n"); 1209682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1210682d7d0cSBarry Smith #if defined(HAVE_ESSL) 1211682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1212682d7d0cSBarry Smith #endif 1213682d7d0cSBarry Smith return 0; 1214682d7d0cSBarry Smith } 1215682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 121617ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 121717ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1218416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1219416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 122017ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 122117ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 122217ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 122317ab2063SBarry Smith MatRelax_SeqAIJ, 122417ab2063SBarry Smith MatTranspose_SeqAIJ, 122517ab2063SBarry Smith MatGetInfo_SeqAIJ,0, 1226f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 122717ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 122817ab2063SBarry Smith MatCompress_SeqAIJ, 122917ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 123017ab2063SBarry Smith MatGetReordering_SeqAIJ, 123117ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 123217ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 123317ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 123417ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1235416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 12363d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1237cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 12387eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1239682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1240f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1241f0b747eeSBarry Smith MatScale_SeqAIJ}; 124217ab2063SBarry Smith 124317ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 124417ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 124517ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 124617ab2063SBarry Smith 124717ab2063SBarry Smith /*@C 1248682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 124917ab2063SBarry Smith (the default uniprocessor PETSc format). 125017ab2063SBarry Smith 125117ab2063SBarry Smith Input Parameters: 125217ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 125317ab2063SBarry Smith . m - number of rows 125417ab2063SBarry Smith . n - number of columns 125517ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 125617ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 125717ab2063SBarry Smith 125817ab2063SBarry Smith Output Parameter: 1259416022c9SBarry Smith . A - the matrix 126017ab2063SBarry Smith 126117ab2063SBarry Smith Notes: 126217ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 126317ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 12640002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 12650002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 126617ab2063SBarry Smith 126717ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1268a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 126917ab2063SBarry Smith allocation. 127017ab2063SBarry Smith 1271682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1272682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1273682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 12746c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 12756c7ebb05SLois Curfman McInnes 12766c7ebb05SLois Curfman McInnes Options Database Keys: 12776c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 12786c7ebb05SLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit (max limit=5) 12794b0e389bSBarry Smith $ -mat_aij_oneindex - Internally using indexing starting at 1 rather then zero. 12804b0e389bSBarry Smith . Note: You still index entries starting at 0! 128117ab2063SBarry Smith 128217ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 128317ab2063SBarry Smith @*/ 1284416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 128517ab2063SBarry Smith { 1286416022c9SBarry Smith Mat B; 1287416022c9SBarry Smith Mat_SeqAIJ *b; 128869957df2SSatish Balay int i,len,ierr, flg; 1289d5d45c9bSBarry Smith 1290416022c9SBarry Smith *A = 0; 12910452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1292416022c9SBarry Smith PLogObjectCreate(B); 12930452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1294416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1295416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1296416022c9SBarry Smith B->view = MatView_SeqAIJ; 1297416022c9SBarry Smith B->factor = 0; 1298416022c9SBarry Smith B->lupivotthreshold = 1.0; 129969957df2SSatish Balay ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, \ 130069957df2SSatish Balay &flg); CHKERRQ(ierr); 1301416022c9SBarry Smith b->row = 0; 1302416022c9SBarry Smith b->col = 0; 1303416022c9SBarry Smith b->indexshift = 0; 1304*b810aeb4SBarry Smith b->reallocs = 0; 130569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 130669957df2SSatish Balay if (flg) b->indexshift = -1; 130717ab2063SBarry Smith 1308416022c9SBarry Smith b->m = m; 1309416022c9SBarry Smith b->n = n; 13100452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1311b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 13127b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 13137b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1314416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 131517ab2063SBarry Smith nz = nz*m; 131617ab2063SBarry Smith } 131717ab2063SBarry Smith else { 131817ab2063SBarry Smith nz = 0; 1319416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 132017ab2063SBarry Smith } 132117ab2063SBarry Smith 132217ab2063SBarry Smith /* allocate the matrix space */ 1323416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 13240452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1325416022c9SBarry Smith b->j = (int *) (b->a + nz); 1326cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1327416022c9SBarry Smith b->i = b->j + nz; 1328416022c9SBarry Smith b->singlemalloc = 1; 132917ab2063SBarry Smith 1330416022c9SBarry Smith b->i[0] = -b->indexshift; 133117ab2063SBarry Smith for (i=1; i<m+1; i++) { 1332416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 133317ab2063SBarry Smith } 133417ab2063SBarry Smith 1335416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 13360452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1337416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1338416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 133917ab2063SBarry Smith 1340416022c9SBarry Smith b->nz = 0; 1341416022c9SBarry Smith b->maxnz = nz; 1342416022c9SBarry Smith b->sorted = 0; 1343416022c9SBarry Smith b->roworiented = 1; 1344416022c9SBarry Smith b->nonew = 0; 1345416022c9SBarry Smith b->diag = 0; 1346416022c9SBarry Smith b->solve_work = 0; 134771bd300dSLois Curfman McInnes b->spptr = 0; 1348754ec7b1SSatish Balay b->inode.node_count = 0; 1349754ec7b1SSatish Balay b->inode.size = 0; 13506c7ebb05SLois Curfman McInnes b->inode.limit = 5; 13516c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 135217ab2063SBarry Smith 1353416022c9SBarry Smith *A = B; 13544b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 13554b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 135669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 135769957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 13584b14c69eSBarry Smith #endif 135969957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 136069957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 136169957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 136269957df2SSatish Balay if (flg) { 1363416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1364416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 136517ab2063SBarry Smith } 136669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 136769957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 136817ab2063SBarry Smith return 0; 136917ab2063SBarry Smith } 137017ab2063SBarry Smith 13713d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 137217ab2063SBarry Smith { 1373416022c9SBarry Smith Mat C; 1374416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 137508480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 137617ab2063SBarry Smith 13774043dd9cSLois Curfman McInnes *B = 0; 13780452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1379416022c9SBarry Smith PLogObjectCreate(C); 13800452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 138141c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1382416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1383416022c9SBarry Smith C->view = MatView_SeqAIJ; 1384416022c9SBarry Smith C->factor = A->factor; 1385416022c9SBarry Smith c->row = 0; 1386416022c9SBarry Smith c->col = 0; 1387416022c9SBarry Smith c->indexshift = shift; 1388c456f294SBarry Smith C->assembled = PETSC_TRUE; 138917ab2063SBarry Smith 1390416022c9SBarry Smith c->m = a->m; 1391416022c9SBarry Smith c->n = a->n; 139217ab2063SBarry Smith 13930452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 13940452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 139517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1396416022c9SBarry Smith c->imax[i] = a->imax[i]; 1397416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 139817ab2063SBarry Smith } 139917ab2063SBarry Smith 140017ab2063SBarry Smith /* allocate the matrix space */ 1401416022c9SBarry Smith c->singlemalloc = 1; 1402416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 14030452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1404416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1405416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1406416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 140717ab2063SBarry Smith if (m > 0) { 1408416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 140908480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1410416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 141117ab2063SBarry Smith } 141208480c60SBarry Smith } 141317ab2063SBarry Smith 1414416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1415416022c9SBarry Smith c->sorted = a->sorted; 1416416022c9SBarry Smith c->roworiented = a->roworiented; 1417416022c9SBarry Smith c->nonew = a->nonew; 141817ab2063SBarry Smith 1419416022c9SBarry Smith if (a->diag) { 14200452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1421416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 142217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1423416022c9SBarry Smith c->diag[i] = a->diag[i]; 142417ab2063SBarry Smith } 142517ab2063SBarry Smith } 1426416022c9SBarry Smith else c->diag = 0; 14276c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 14286c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1429754ec7b1SSatish Balay if (a->inode.size){ 1430754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1431754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1432754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1433754ec7b1SSatish Balay } else { 1434754ec7b1SSatish Balay c->inode.size = 0; 1435754ec7b1SSatish Balay c->inode.node_count = 0; 1436754ec7b1SSatish Balay } 1437416022c9SBarry Smith c->nz = a->nz; 1438416022c9SBarry Smith c->maxnz = a->maxnz; 1439416022c9SBarry Smith c->solve_work = 0; 144076dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1441754ec7b1SSatish Balay 1442416022c9SBarry Smith *B = C; 144317ab2063SBarry Smith return 0; 144417ab2063SBarry Smith } 144517ab2063SBarry Smith 1446416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 144717ab2063SBarry Smith { 1448416022c9SBarry Smith Mat_SeqAIJ *a; 1449416022c9SBarry Smith Mat B; 145017699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 145117ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 145217ab2063SBarry Smith MPI_Comm comm = vobj->comm; 145317ab2063SBarry Smith 145417699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 145517699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 145617ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1457416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 145848d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 145917ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 146017ab2063SBarry Smith 146117ab2063SBarry Smith /* read in row lengths */ 14620452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1463416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 146417ab2063SBarry Smith 146517ab2063SBarry Smith /* create our matrix */ 1466416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1467416022c9SBarry Smith B = *A; 1468416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1469416022c9SBarry Smith shift = a->indexshift; 147017ab2063SBarry Smith 147117ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1472416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 147317ab2063SBarry Smith if (shift) { 147417ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1475416022c9SBarry Smith a->j[i] += 1; 147617ab2063SBarry Smith } 147717ab2063SBarry Smith } 147817ab2063SBarry Smith 147917ab2063SBarry Smith /* read in nonzero values */ 1480416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 148117ab2063SBarry Smith 148217ab2063SBarry Smith /* set matrix "i" values */ 1483416022c9SBarry Smith a->i[0] = -shift; 148417ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1485416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1486416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 148717ab2063SBarry Smith } 14880452661fSBarry Smith PetscFree(rowlengths); 148917ab2063SBarry Smith 1490416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1491416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 149217ab2063SBarry Smith return 0; 149317ab2063SBarry Smith } 149417ab2063SBarry Smith 149517ab2063SBarry Smith 149617ab2063SBarry Smith 1497