16d84be18SBarry Smith 2*6945ee14SBarry Smith 317ab2063SBarry Smith #ifndef lint 4*6945ee14SBarry Smith static char vcid[] = "$Id: aij.c,v 1.183 1996/08/22 19:52:29 curfman Exp bsmith $"; 517ab2063SBarry Smith #endif 617ab2063SBarry Smith 7d5d45c9bSBarry Smith /* 85a838052SSatish Balay B Defines the basic matrix operations for the AIJ (compressed row) 9d5d45c9bSBarry Smith matrix storage format. 10d5d45c9bSBarry Smith */ 1170f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 12f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 13f5eb4b81SSatish Balay #include "src/inline/spops.h" 14e4d965acSSatish Balay #include "petsc.h" 15f5eb4b81SSatish Balay #include "src/inline/bitarray.h" 1617ab2063SBarry Smith 17a2ce50c7SBarry Smith /* 18a2ce50c7SBarry Smith Basic AIJ format ILU based on drop tolerance 19a2ce50c7SBarry Smith */ 20a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact) 21a2ce50c7SBarry Smith { 22a2ce50c7SBarry Smith /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */ 23a2ce50c7SBarry Smith int ierr = 1; 24a2ce50c7SBarry Smith 25a2ce50c7SBarry Smith SETERRQ(ierr,"MatILUDTFactor_SeqAIJ:Not implemented"); 26a2ce50c7SBarry Smith } 27a2ce50c7SBarry Smith 28bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 2917ab2063SBarry Smith 30*6945ee14SBarry Smith static int MatGetIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 31*6945ee14SBarry Smith PetscTruth *done) 3217ab2063SBarry Smith { 33416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 34*6945ee14SBarry Smith int ierr,i,ishift; 3517ab2063SBarry Smith 36*6945ee14SBarry Smith *n = A->n; 37*6945ee14SBarry Smith if (!ia) return 0; 38*6945ee14SBarry Smith ishift = a->indexshift; 39*6945ee14SBarry Smith if (symmetric) { 40*6945ee14SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 41*6945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 42*6945ee14SBarry Smith /* temporarily subtract 1 from i and j indices */ 43*6945ee14SBarry Smith int nz = a->i[a->n]; 44*6945ee14SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 45*6945ee14SBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 46*6945ee14SBarry Smith *ia = a->i; *ja = a->j; 47*6945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 48*6945ee14SBarry Smith /* temporarily add 1 to i and j indices */ 49*6945ee14SBarry Smith int nz = a->i[a->n] + 1; 50*6945ee14SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 51*6945ee14SBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 52*6945ee14SBarry Smith *ia = a->i; *ja = a->j; 53*6945ee14SBarry Smith } else { 54*6945ee14SBarry Smith *ia = a->i; *ja = a->j; 55a2ce50c7SBarry Smith } 56a2ce50c7SBarry Smith 57a2744918SBarry Smith return 0; 58a2744918SBarry Smith } 59a2744918SBarry Smith 60*6945ee14SBarry Smith static int MatRestoreIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 61*6945ee14SBarry Smith PetscTruth *done) 62*6945ee14SBarry Smith { 63*6945ee14SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 64*6945ee14SBarry Smith int i,ishift; 65*6945ee14SBarry Smith 66*6945ee14SBarry Smith if (!ia) return 0; 67bcd2baecSBarry Smith ishift = a->indexshift; 68*6945ee14SBarry Smith if (symmetric) { 69*6945ee14SBarry Smith PetscFree(*ia); 70*6945ee14SBarry Smith PetscFree(*ja); 71*6945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 72*6945ee14SBarry Smith int nz = a->i[a->n] + 1; 73*6945ee14SBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 74*6945ee14SBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 75*6945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 76bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 77bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 78bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 79bcd2baecSBarry Smith } 8017ab2063SBarry Smith return 0; 8117ab2063SBarry Smith } 8217ab2063SBarry Smith 83227d817aSBarry Smith #define CHUNKSIZE 15 8417ab2063SBarry Smith 8517ab2063SBarry Smith /* This version has row oriented v */ 86416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 8717ab2063SBarry Smith { 88416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 89416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 904b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 91d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 92416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 9317ab2063SBarry Smith 9417ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 95416022c9SBarry Smith row = im[k]; 9617ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 97416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 9817ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 9917ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 100416022c9SBarry Smith low = 0; 10117ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 102416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 103416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 1044b0e389bSBarry Smith col = in[l] - shift; 1054b0e389bSBarry Smith if (roworiented) { 1064b0e389bSBarry Smith value = *v++; 1074b0e389bSBarry Smith } 1084b0e389bSBarry Smith else { 1094b0e389bSBarry Smith value = v[k + l*m]; 1104b0e389bSBarry Smith } 111416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 112416022c9SBarry Smith while (high-low > 5) { 113416022c9SBarry Smith t = (low+high)/2; 114416022c9SBarry Smith if (rp[t] > col) high = t; 115416022c9SBarry Smith else low = t; 11617ab2063SBarry Smith } 117416022c9SBarry Smith for ( i=low; i<high; i++ ) { 11817ab2063SBarry Smith if (rp[i] > col) break; 11917ab2063SBarry Smith if (rp[i] == col) { 120416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 12117ab2063SBarry Smith else ap[i] = value; 12217ab2063SBarry Smith goto noinsert; 12317ab2063SBarry Smith } 12417ab2063SBarry Smith } 12517ab2063SBarry Smith if (nonew) goto noinsert; 12617ab2063SBarry Smith if (nrow >= rmax) { 12717ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 128416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 12917ab2063SBarry Smith Scalar *new_a; 13017ab2063SBarry Smith 13117ab2063SBarry Smith /* malloc new storage space */ 132416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1330452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 13417ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 13517ab2063SBarry Smith new_i = new_j + new_nz; 13617ab2063SBarry Smith 13717ab2063SBarry Smith /* copy over old data into new slots */ 13817ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 139416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 140416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 141416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 142416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 14317ab2063SBarry Smith len*sizeof(int)); 144416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 145416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 14617ab2063SBarry Smith len*sizeof(Scalar)); 14717ab2063SBarry Smith /* free up old matrix storage */ 1480452661fSBarry Smith PetscFree(a->a); 1490452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 150416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 151416022c9SBarry Smith a->singlemalloc = 1; 15217ab2063SBarry Smith 15317ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 154416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 155416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 156416022c9SBarry Smith a->maxnz += CHUNKSIZE; 157b810aeb4SBarry Smith a->reallocs++; 15817ab2063SBarry Smith } 159416022c9SBarry Smith N = nrow++ - 1; a->nz++; 160416022c9SBarry Smith /* shift up all the later entries in this row */ 161416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 16217ab2063SBarry Smith rp[ii+1] = rp[ii]; 16317ab2063SBarry Smith ap[ii+1] = ap[ii]; 16417ab2063SBarry Smith } 16517ab2063SBarry Smith rp[i] = col; 16617ab2063SBarry Smith ap[i] = value; 16717ab2063SBarry Smith noinsert:; 168416022c9SBarry Smith low = i + 1; 16917ab2063SBarry Smith } 17017ab2063SBarry Smith ailen[row] = nrow; 17117ab2063SBarry Smith } 17217ab2063SBarry Smith return 0; 17317ab2063SBarry Smith } 17417ab2063SBarry Smith 1757eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 1767eb43aa7SLois Curfman McInnes { 1777eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 178b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 1797eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 1807eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 1817eb43aa7SLois Curfman McInnes 1827eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 1837eb43aa7SLois Curfman McInnes row = im[k]; 1847eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 1857eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 1867eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 1877eb43aa7SLois Curfman McInnes nrow = ailen[row]; 1887eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 1897eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 1907eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 1917eb43aa7SLois Curfman McInnes col = in[l] - shift; 1927eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 1937eb43aa7SLois Curfman McInnes while (high-low > 5) { 1947eb43aa7SLois Curfman McInnes t = (low+high)/2; 1957eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 1967eb43aa7SLois Curfman McInnes else low = t; 1977eb43aa7SLois Curfman McInnes } 1987eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 1997eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2007eb43aa7SLois Curfman McInnes if (rp[i] == col) { 201b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2027eb43aa7SLois Curfman McInnes goto finished; 2037eb43aa7SLois Curfman McInnes } 2047eb43aa7SLois Curfman McInnes } 205b49de8d1SLois Curfman McInnes *v++ = zero; 2067eb43aa7SLois Curfman McInnes finished:; 2077eb43aa7SLois Curfman McInnes } 2087eb43aa7SLois Curfman McInnes } 2097eb43aa7SLois Curfman McInnes return 0; 2107eb43aa7SLois Curfman McInnes } 2117eb43aa7SLois Curfman McInnes 21217ab2063SBarry Smith #include "draw.h" 21317ab2063SBarry Smith #include "pinclude/pviewer.h" 21477c4ece6SBarry Smith #include "sys.h" 21517ab2063SBarry Smith 216416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 21717ab2063SBarry Smith { 218416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 219416022c9SBarry Smith int i, fd, *col_lens, ierr; 22017ab2063SBarry Smith 22190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2220452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 223416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 224416022c9SBarry Smith col_lens[1] = a->m; 225416022c9SBarry Smith col_lens[2] = a->n; 226416022c9SBarry Smith col_lens[3] = a->nz; 227416022c9SBarry Smith 228416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 229416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 230416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 23117ab2063SBarry Smith } 23277c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 2330452661fSBarry Smith PetscFree(col_lens); 234416022c9SBarry Smith 235416022c9SBarry Smith /* store column indices (zero start index) */ 236416022c9SBarry Smith if (a->indexshift) { 237416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 23817ab2063SBarry Smith } 23977c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 240416022c9SBarry Smith if (a->indexshift) { 241416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 24217ab2063SBarry Smith } 243416022c9SBarry Smith 244416022c9SBarry Smith /* store nonzero values */ 24577c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 24617ab2063SBarry Smith return 0; 24717ab2063SBarry Smith } 248416022c9SBarry Smith 249416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 250416022c9SBarry Smith { 251416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 252496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 25317ab2063SBarry Smith FILE *fd; 25417ab2063SBarry Smith char *outputname; 25517ab2063SBarry Smith 25690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 257416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 25890ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 25990ace30eSBarry Smith if (format == ASCII_FORMAT_INFO) { 26095e01e2fSLois Curfman McInnes return 0; 26195e01e2fSLois Curfman McInnes } 26290ace30eSBarry Smith else if (format == ASCII_FORMAT_INFO_DETAILED) { 263496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 264496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 265496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 26695e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 26795e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 26817ab2063SBarry Smith } 26990ace30eSBarry Smith else if (format == ASCII_FORMAT_MATLAB) { 270416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 2714e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 2724e220ebcSLois Curfman McInnes fprintf(fd,"zzz = zeros(%d,3);\n",a->nz); 27317ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 27417ab2063SBarry Smith 27517ab2063SBarry Smith for (i=0; i<m; i++) { 276416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 27717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 278*6945ee14SBarry Smith fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]), 279416022c9SBarry Smith imag(a->a[j])); 28017ab2063SBarry Smith #else 2817a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 28217ab2063SBarry Smith #endif 28317ab2063SBarry Smith } 28417ab2063SBarry Smith } 28517ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 28617ab2063SBarry Smith } 28744cd7ae7SLois Curfman McInnes else if (format == ASCII_FORMAT_COMMON) { 28844cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 28944cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 29044cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 29144cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 29244cd7ae7SLois Curfman McInnes if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0) 29344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 29444cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 29544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 29644cd7ae7SLois Curfman McInnes #else 29744cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 29844cd7ae7SLois Curfman McInnes #endif 29944cd7ae7SLois Curfman McInnes } 30044cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 30144cd7ae7SLois Curfman McInnes } 30244cd7ae7SLois Curfman McInnes } 30317ab2063SBarry Smith else { 30417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 30517ab2063SBarry Smith fprintf(fd,"row %d:",i); 306416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 30717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 308416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 309416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 31017ab2063SBarry Smith } 31117ab2063SBarry Smith else { 312416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 31317ab2063SBarry Smith } 31417ab2063SBarry Smith #else 315416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 31617ab2063SBarry Smith #endif 31717ab2063SBarry Smith } 31817ab2063SBarry Smith fprintf(fd,"\n"); 31917ab2063SBarry Smith } 32017ab2063SBarry Smith } 32117ab2063SBarry Smith fflush(fd); 322416022c9SBarry Smith return 0; 323416022c9SBarry Smith } 324416022c9SBarry Smith 325416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 326416022c9SBarry Smith { 327416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 328cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 329cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 330bcd2baecSBarry Smith Draw draw; 331cddf8d76SBarry Smith DrawButton button; 33219bcc07fSBarry Smith PetscTruth isnull; 333cddf8d76SBarry Smith 334bcd2baecSBarry Smith ViewerDrawGetDraw(viewer,&draw); 33519bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 33619bcc07fSBarry Smith 337416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 338416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 339416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 340416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 341cddf8d76SBarry Smith color = DRAW_BLUE; 342416022c9SBarry Smith for ( i=0; i<m; i++ ) { 343cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 344416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 345cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 346cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 347cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 348cddf8d76SBarry Smith #else 349cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 350cddf8d76SBarry Smith #endif 351cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 352cddf8d76SBarry Smith } 353cddf8d76SBarry Smith } 354cddf8d76SBarry Smith color = DRAW_CYAN; 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 (a->a[j] != 0.) continue; 360cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 361cddf8d76SBarry Smith } 362cddf8d76SBarry Smith } 363cddf8d76SBarry Smith color = DRAW_RED; 364cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 365cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 366cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 367cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 368cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 369cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 370cddf8d76SBarry Smith #else 371cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 372cddf8d76SBarry Smith #endif 373cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 374416022c9SBarry Smith } 375416022c9SBarry Smith } 376416022c9SBarry Smith DrawFlush(draw); 377cddf8d76SBarry Smith DrawGetPause(draw,&pause); 378cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 379cddf8d76SBarry Smith 380cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 381*6945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 382cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 383cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 384cddf8d76SBarry Smith DrawClear(draw); 385cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 386cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 387cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 388cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 389cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 390cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 391cddf8d76SBarry Smith w *= scale; h *= scale; 392cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 393cddf8d76SBarry Smith color = DRAW_BLUE; 394cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 395cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 396cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 397cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 398cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 399cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 400cddf8d76SBarry Smith #else 401cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 402cddf8d76SBarry Smith #endif 403cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 404cddf8d76SBarry Smith } 405cddf8d76SBarry Smith } 406cddf8d76SBarry Smith color = DRAW_CYAN; 407cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 408cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 409cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 410cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 411cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 412cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 413cddf8d76SBarry Smith } 414cddf8d76SBarry Smith } 415cddf8d76SBarry Smith color = DRAW_RED; 416cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 417cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 418cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 419cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 420cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 421cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 422cddf8d76SBarry Smith #else 423cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 424cddf8d76SBarry Smith #endif 425cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 426cddf8d76SBarry Smith } 427cddf8d76SBarry Smith } 428*6945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 429cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 430cddf8d76SBarry Smith } 431416022c9SBarry Smith return 0; 432416022c9SBarry Smith } 433416022c9SBarry Smith 434416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 435416022c9SBarry Smith { 436416022c9SBarry Smith Mat A = (Mat) obj; 437416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 438bcd2baecSBarry Smith ViewerType vtype; 439bcd2baecSBarry Smith int ierr; 440416022c9SBarry Smith 441bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 442bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 443416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 444416022c9SBarry Smith } 445bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 446416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 447416022c9SBarry Smith } 448bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 449416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 450416022c9SBarry Smith } 451bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 452bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 45317ab2063SBarry Smith } 45417ab2063SBarry Smith return 0; 45517ab2063SBarry Smith } 45619bcc07fSBarry Smith 457c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 458416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 45917ab2063SBarry Smith { 460416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 46141c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 462416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 463416022c9SBarry Smith Scalar *aa = a->a, *ap; 46417ab2063SBarry Smith 4656d4a8577SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) return 0; 46617ab2063SBarry Smith 46717ab2063SBarry Smith for ( i=1; i<m; i++ ) { 468416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 46917ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 47017ab2063SBarry Smith if (fshift) { 471416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 47217ab2063SBarry Smith N = ailen[i]; 47317ab2063SBarry Smith for ( j=0; j<N; j++ ) { 47417ab2063SBarry Smith ip[j-fshift] = ip[j]; 47517ab2063SBarry Smith ap[j-fshift] = ap[j]; 47617ab2063SBarry Smith } 47717ab2063SBarry Smith } 47817ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 47917ab2063SBarry Smith } 48017ab2063SBarry Smith if (m) { 48117ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 48217ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 48317ab2063SBarry Smith } 48417ab2063SBarry Smith /* reset ilen and imax for each row */ 48517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 48617ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 48717ab2063SBarry Smith } 488416022c9SBarry Smith a->nz = ai[m] + shift; 48917ab2063SBarry Smith 49017ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 491416022c9SBarry Smith if (fshift && a->diag) { 4920452661fSBarry Smith PetscFree(a->diag); 493416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 494416022c9SBarry Smith a->diag = 0; 49517ab2063SBarry Smith } 4964e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 4974e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 4984e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 499b810aeb4SBarry Smith a->reallocs); 5004e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 5014e220ebcSLois Curfman McInnes 50276dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 50341c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 50417ab2063SBarry Smith return 0; 50517ab2063SBarry Smith } 50617ab2063SBarry Smith 507416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 50817ab2063SBarry Smith { 509416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 510cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 51117ab2063SBarry Smith return 0; 51217ab2063SBarry Smith } 513416022c9SBarry Smith 51417ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 51517ab2063SBarry Smith { 516416022c9SBarry Smith Mat A = (Mat) obj; 517416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 518d5d45c9bSBarry Smith 51917ab2063SBarry Smith #if defined(PETSC_LOG) 520416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 52117ab2063SBarry Smith #endif 5220452661fSBarry Smith PetscFree(a->a); 5230452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 5240452661fSBarry Smith if (a->diag) PetscFree(a->diag); 5250452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 5260452661fSBarry Smith if (a->imax) PetscFree(a->imax); 5270452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 52876dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 5290452661fSBarry Smith PetscFree(a); 530f2655603SLois Curfman McInnes PLogObjectDestroy(A); 531f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 53217ab2063SBarry Smith return 0; 53317ab2063SBarry Smith } 53417ab2063SBarry Smith 535416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 53617ab2063SBarry Smith { 53717ab2063SBarry Smith return 0; 53817ab2063SBarry Smith } 53917ab2063SBarry Smith 540416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 54117ab2063SBarry Smith { 542416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 5436d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 5446d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 5456d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 5466d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 5476d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 5486d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 5496d4a8577SBarry Smith op == MAT_SYMMETRIC || 5506d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 5516d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 55294a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 5536d4a8577SBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) 5546d4a8577SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:MAT_NO_NEW_DIAGONALS");} 5556d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 5566d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 5576d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 5586d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 5596d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 560e2f28af5SBarry Smith else 5614d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 56217ab2063SBarry Smith return 0; 56317ab2063SBarry Smith } 56417ab2063SBarry Smith 565416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 56617ab2063SBarry Smith { 567416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 568416022c9SBarry Smith int i,j, n,shift = a->indexshift; 56917ab2063SBarry Smith Scalar *x, zero = 0.0; 57017ab2063SBarry Smith 57117ab2063SBarry Smith VecSet(&zero,v); 57217ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 573416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 574416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 575416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 576416022c9SBarry Smith if (a->j[j]+shift == i) { 577416022c9SBarry Smith x[i] = a->a[j]; 57817ab2063SBarry Smith break; 57917ab2063SBarry Smith } 58017ab2063SBarry Smith } 58117ab2063SBarry Smith } 58217ab2063SBarry Smith return 0; 58317ab2063SBarry Smith } 58417ab2063SBarry Smith 58517ab2063SBarry Smith /* -------------------------------------------------------*/ 58617ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 58717ab2063SBarry Smith /* -------------------------------------------------------*/ 58844cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 58917ab2063SBarry Smith { 590416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 59117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 592416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 59317ab2063SBarry Smith 59417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 595cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 59617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 59717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 598416022c9SBarry Smith idx = a->j + a->i[i] + shift; 599416022c9SBarry Smith v = a->a + a->i[i] + shift; 600416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 60117ab2063SBarry Smith alpha = x[i]; 60217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 60317ab2063SBarry Smith } 604416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 60517ab2063SBarry Smith return 0; 60617ab2063SBarry Smith } 607d5d45c9bSBarry Smith 60844cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 60917ab2063SBarry Smith { 610416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 61117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 612416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 61317ab2063SBarry Smith 61417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 61517ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 61617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 61717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 618416022c9SBarry Smith idx = a->j + a->i[i] + shift; 619416022c9SBarry Smith v = a->a + a->i[i] + shift; 620416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 62117ab2063SBarry Smith alpha = x[i]; 62217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 62317ab2063SBarry Smith } 62417ab2063SBarry Smith return 0; 62517ab2063SBarry Smith } 62617ab2063SBarry Smith 62744cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 62817ab2063SBarry Smith { 629416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 63017ab2063SBarry Smith Scalar *x, *y, *v, sum; 631416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 63217ab2063SBarry Smith 63317ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 63417ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 635416022c9SBarry Smith idx = a->j; 636416022c9SBarry Smith v = a->a; 637416022c9SBarry Smith ii = a->i; 63817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 639416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 64017ab2063SBarry Smith sum = 0.0; 64117ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 64217ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 643416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 64417ab2063SBarry Smith y[i] = sum; 64517ab2063SBarry Smith } 646416022c9SBarry Smith PLogFlops(2*a->nz - m); 64717ab2063SBarry Smith return 0; 64817ab2063SBarry Smith } 64917ab2063SBarry Smith 65044cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 65117ab2063SBarry Smith { 652416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 65317ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 654cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 65517ab2063SBarry Smith 65617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 65717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 658cddf8d76SBarry Smith idx = a->j; 659cddf8d76SBarry Smith v = a->a; 660cddf8d76SBarry Smith ii = a->i; 66117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 662cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 66317ab2063SBarry Smith sum = y[i]; 664cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 665cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 66617ab2063SBarry Smith z[i] = sum; 66717ab2063SBarry Smith } 668416022c9SBarry Smith PLogFlops(2*a->nz); 66917ab2063SBarry Smith return 0; 67017ab2063SBarry Smith } 67117ab2063SBarry Smith 67217ab2063SBarry Smith /* 67317ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 67417ab2063SBarry Smith */ 67517ab2063SBarry Smith 676416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 67717ab2063SBarry Smith { 678416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 679416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 68017ab2063SBarry Smith 6810452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 682416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 683416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 684416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 685416022c9SBarry Smith if (a->j[j]+shift == i) { 68617ab2063SBarry Smith diag[i] = j - shift; 68717ab2063SBarry Smith break; 68817ab2063SBarry Smith } 68917ab2063SBarry Smith } 69017ab2063SBarry Smith } 691416022c9SBarry Smith a->diag = diag; 69217ab2063SBarry Smith return 0; 69317ab2063SBarry Smith } 69417ab2063SBarry Smith 69544cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 69617ab2063SBarry Smith double fshift,int its,Vec xx) 69717ab2063SBarry Smith { 698416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 699416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 700d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 70117ab2063SBarry Smith 70217ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 703416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 704416022c9SBarry Smith diag = a->diag; 705416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 70617ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 70717ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 70817ab2063SBarry Smith bs = b + shift; 70917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 710416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 711416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 712416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 713416022c9SBarry Smith v = a->a + diag[i] + (!shift); 71417ab2063SBarry Smith sum = b[i]*d/omega; 71517ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 71617ab2063SBarry Smith x[i] = sum; 71717ab2063SBarry Smith } 71817ab2063SBarry Smith return 0; 71917ab2063SBarry Smith } 72017ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 72117ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 72217ab2063SBarry Smith } 723416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 72417ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 72517ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 72617ab2063SBarry Smith 72717ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 72817ab2063SBarry Smith 72917ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 73017ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 73117ab2063SBarry Smith is the relaxation factor. 73217ab2063SBarry Smith */ 7330452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 73417ab2063SBarry Smith scale = (2.0/omega) - 1.0; 73517ab2063SBarry Smith 73617ab2063SBarry Smith /* x = (E + U)^{-1} b */ 73717ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 738416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 739416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 740416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 741416022c9SBarry Smith v = a->a + diag[i] + (!shift); 74217ab2063SBarry Smith sum = b[i]; 74317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 74417ab2063SBarry Smith x[i] = omega*(sum/d); 74517ab2063SBarry Smith } 74617ab2063SBarry Smith 74717ab2063SBarry Smith /* t = b - (2*E - D)x */ 748416022c9SBarry Smith v = a->a; 74917ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 75017ab2063SBarry Smith 75117ab2063SBarry Smith /* t = (E + L)^{-1}t */ 752416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 753416022c9SBarry Smith diag = a->diag; 75417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 755416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 756416022c9SBarry Smith n = diag[i] - a->i[i]; 757416022c9SBarry Smith idx = a->j + a->i[i] + shift; 758416022c9SBarry Smith v = a->a + a->i[i] + shift; 75917ab2063SBarry Smith sum = t[i]; 76017ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 76117ab2063SBarry Smith t[i] = omega*(sum/d); 76217ab2063SBarry Smith } 76317ab2063SBarry Smith 76417ab2063SBarry Smith /* x = x + t */ 76517ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7660452661fSBarry Smith PetscFree(t); 76717ab2063SBarry Smith return 0; 76817ab2063SBarry Smith } 76917ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 77017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 77117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 772416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 773416022c9SBarry Smith n = diag[i] - a->i[i]; 774416022c9SBarry Smith idx = a->j + a->i[i] + shift; 775416022c9SBarry Smith v = a->a + a->i[i] + shift; 77617ab2063SBarry Smith sum = b[i]; 77717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 77817ab2063SBarry Smith x[i] = omega*(sum/d); 77917ab2063SBarry Smith } 78017ab2063SBarry Smith xb = x; 78117ab2063SBarry Smith } 78217ab2063SBarry Smith else xb = b; 78317ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 78417ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 78517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 786416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 78717ab2063SBarry Smith } 78817ab2063SBarry Smith } 78917ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 79017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 791416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 792416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 793416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 794416022c9SBarry Smith v = a->a + diag[i] + (!shift); 79517ab2063SBarry Smith sum = xb[i]; 79617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 79717ab2063SBarry Smith x[i] = omega*(sum/d); 79817ab2063SBarry Smith } 79917ab2063SBarry Smith } 80017ab2063SBarry Smith its--; 80117ab2063SBarry Smith } 80217ab2063SBarry Smith while (its--) { 80317ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 80417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 805416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 806416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 807416022c9SBarry Smith idx = a->j + a->i[i] + shift; 808416022c9SBarry Smith v = a->a + a->i[i] + shift; 80917ab2063SBarry Smith sum = b[i]; 81017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 81117ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 81217ab2063SBarry Smith } 81317ab2063SBarry Smith } 81417ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 81517ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 816416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 817416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 818416022c9SBarry Smith idx = a->j + a->i[i] + shift; 819416022c9SBarry Smith v = a->a + a->i[i] + shift; 82017ab2063SBarry Smith sum = b[i]; 82117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 82217ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 82317ab2063SBarry Smith } 82417ab2063SBarry Smith } 82517ab2063SBarry Smith } 82617ab2063SBarry Smith return 0; 82717ab2063SBarry Smith } 82817ab2063SBarry Smith 8294e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 83017ab2063SBarry Smith { 831416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 8324e220ebcSLois Curfman McInnes 8334e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 8344e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 8354e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 8364e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 8374e220ebcSLois Curfman McInnes info->block_size = 1.0; 8384e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 8394e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 8404e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 8414e220ebcSLois Curfman McInnes /* if (info->nz_unneeded != A->info.nz_unneeded) 8424e220ebcSLois Curfman McInnes printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 8434e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 8444e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 8454e220ebcSLois Curfman McInnes info->memory = A->mem; 8464e220ebcSLois Curfman McInnes if (A->factor) { 8474e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 8484e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 8494e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 8504e220ebcSLois Curfman McInnes } else { 8514e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 8524e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 8534e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 8544e220ebcSLois Curfman McInnes } 85517ab2063SBarry Smith return 0; 85617ab2063SBarry Smith } 85717ab2063SBarry Smith 85817ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 85917ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 86017ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 86117ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 86217ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 86317ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 86417ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 86517ab2063SBarry Smith 86617ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 86717ab2063SBarry Smith { 868416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 869416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 87017ab2063SBarry Smith 87177c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 87217ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 87317ab2063SBarry Smith if (diag) { 87417ab2063SBarry Smith for ( i=0; i<N; i++ ) { 875416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 876416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 877416022c9SBarry Smith a->ilen[rows[i]] = 1; 878416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 879416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 88017ab2063SBarry Smith } 88117ab2063SBarry Smith else { 88217ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 88317ab2063SBarry Smith CHKERRQ(ierr); 88417ab2063SBarry Smith } 88517ab2063SBarry Smith } 88617ab2063SBarry Smith } 88717ab2063SBarry Smith else { 88817ab2063SBarry Smith for ( i=0; i<N; i++ ) { 889416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 890416022c9SBarry Smith a->ilen[rows[i]] = 0; 89117ab2063SBarry Smith } 89217ab2063SBarry Smith } 89317ab2063SBarry Smith ISRestoreIndices(is,&rows); 8946d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8956d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 89617ab2063SBarry Smith return 0; 89717ab2063SBarry Smith } 89817ab2063SBarry Smith 899416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 90017ab2063SBarry Smith { 901416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 902416022c9SBarry Smith *m = a->m; *n = a->n; 90317ab2063SBarry Smith return 0; 90417ab2063SBarry Smith } 90517ab2063SBarry Smith 906416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 90717ab2063SBarry Smith { 908416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 909416022c9SBarry Smith *m = 0; *n = a->m; 91017ab2063SBarry Smith return 0; 91117ab2063SBarry Smith } 9124e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 91317ab2063SBarry Smith { 914416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 915c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 91617ab2063SBarry Smith 917416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 91817ab2063SBarry Smith 919416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 920416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 92117ab2063SBarry Smith if (idx) { 922416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 9234e093b46SBarry Smith if (*nz && shift) { 9240452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 92517ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 9264e093b46SBarry Smith } else if (*nz) { 9274e093b46SBarry Smith *idx = itmp; 92817ab2063SBarry Smith } 92917ab2063SBarry Smith else *idx = 0; 93017ab2063SBarry Smith } 93117ab2063SBarry Smith return 0; 93217ab2063SBarry Smith } 93317ab2063SBarry Smith 9344e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 93517ab2063SBarry Smith { 9364e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 9374e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 93817ab2063SBarry Smith return 0; 93917ab2063SBarry Smith } 94017ab2063SBarry Smith 941cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 94217ab2063SBarry Smith { 943416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 944416022c9SBarry Smith Scalar *v = a->a; 94517ab2063SBarry Smith double sum = 0.0; 946416022c9SBarry Smith int i, j,shift = a->indexshift; 94717ab2063SBarry Smith 94817ab2063SBarry Smith if (type == NORM_FROBENIUS) { 949416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 95017ab2063SBarry Smith #if defined(PETSC_COMPLEX) 95117ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 95217ab2063SBarry Smith #else 95317ab2063SBarry Smith sum += (*v)*(*v); v++; 95417ab2063SBarry Smith #endif 95517ab2063SBarry Smith } 95617ab2063SBarry Smith *norm = sqrt(sum); 95717ab2063SBarry Smith } 95817ab2063SBarry Smith else if (type == NORM_1) { 95917ab2063SBarry Smith double *tmp; 960416022c9SBarry Smith int *jj = a->j; 9610452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 962cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 96317ab2063SBarry Smith *norm = 0.0; 964416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 965a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 96617ab2063SBarry Smith } 967416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 96817ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 96917ab2063SBarry Smith } 9700452661fSBarry Smith PetscFree(tmp); 97117ab2063SBarry Smith } 97217ab2063SBarry Smith else if (type == NORM_INFINITY) { 97317ab2063SBarry Smith *norm = 0.0; 974416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 975416022c9SBarry Smith v = a->a + a->i[j] + shift; 97617ab2063SBarry Smith sum = 0.0; 977416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 978cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 97917ab2063SBarry Smith } 98017ab2063SBarry Smith if (sum > *norm) *norm = sum; 98117ab2063SBarry Smith } 98217ab2063SBarry Smith } 98317ab2063SBarry Smith else { 98448d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 98517ab2063SBarry Smith } 98617ab2063SBarry Smith return 0; 98717ab2063SBarry Smith } 98817ab2063SBarry Smith 989416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 99017ab2063SBarry Smith { 991416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 992416022c9SBarry Smith Mat C; 993416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 994416022c9SBarry Smith int shift = a->indexshift; 995d5d45c9bSBarry Smith Scalar *array = a->a; 99617ab2063SBarry Smith 9973638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 998dfa27b74SSatish Balay SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place"); 9990452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1000cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 100117ab2063SBarry Smith if (shift) { 100217ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 100317ab2063SBarry Smith } 100417ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1005416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 10060452661fSBarry Smith PetscFree(col); 100717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 100817ab2063SBarry Smith len = ai[i+1]-ai[i]; 1009416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 101017ab2063SBarry Smith array += len; aj += len; 101117ab2063SBarry Smith } 101217ab2063SBarry Smith if (shift) { 101317ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 101417ab2063SBarry Smith } 101517ab2063SBarry Smith 10166d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10176d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 101817ab2063SBarry Smith 10193638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1020416022c9SBarry Smith *B = C; 102117ab2063SBarry Smith } else { 1022416022c9SBarry Smith /* This isn't really an in-place transpose */ 10230452661fSBarry Smith PetscFree(a->a); 10240452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 10250452661fSBarry Smith if (a->diag) PetscFree(a->diag); 10260452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 10270452661fSBarry Smith if (a->imax) PetscFree(a->imax); 10280452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 10291073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 10300452661fSBarry Smith PetscFree(a); 1031416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 10320452661fSBarry Smith PetscHeaderDestroy(C); 103317ab2063SBarry Smith } 103417ab2063SBarry Smith return 0; 103517ab2063SBarry Smith } 103617ab2063SBarry Smith 1037f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 103817ab2063SBarry Smith { 1039416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 104017ab2063SBarry Smith Scalar *l,*r,x,*v; 1041d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 104217ab2063SBarry Smith 104317ab2063SBarry Smith if (ll) { 104417ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 1045f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 1046416022c9SBarry Smith v = a->a; 104717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 104817ab2063SBarry Smith x = l[i]; 1049416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 105017ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 105117ab2063SBarry Smith } 105244cd7ae7SLois Curfman McInnes PLogFlops(nz); 105317ab2063SBarry Smith } 105417ab2063SBarry Smith if (rr) { 105517ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 1056f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 1057416022c9SBarry Smith v = a->a; jj = a->j; 105817ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 105917ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 106017ab2063SBarry Smith } 106144cd7ae7SLois Curfman McInnes PLogFlops(nz); 106217ab2063SBarry Smith } 106317ab2063SBarry Smith return 0; 106417ab2063SBarry Smith } 106517ab2063SBarry Smith 1066cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 106717ab2063SBarry Smith { 1068db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 106902834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 107099141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1071a2744918SBarry Smith register int sum,lensi; 107299141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 107399141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 107499141d43SSatish Balay Scalar *a_new,*mat_a; 1075416022c9SBarry Smith Mat C; 107617ab2063SBarry Smith 1077b48a1e75SSatish Balay ierr = ISSorted(isrow,(PetscTruth*)&i); 1078b48a1e75SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:ISrow is not sorted"); 107999141d43SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i); 1080b48a1e75SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IScol is not sorted"); 108199141d43SSatish Balay 108217ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 108317ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 108417ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 108517ab2063SBarry Smith 10867264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 108702834360SBarry Smith /* special case of contiguous rows */ 108857faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 108902834360SBarry Smith starts = lens + ncols; 109002834360SBarry Smith /* loop over new rows determining lens and starting points */ 109102834360SBarry Smith for (i=0; i<nrows; i++) { 1092a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1093a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 109402834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1095d8ced48eSBarry Smith if (aj[k]+shift >= first) { 109602834360SBarry Smith starts[i] = k; 109702834360SBarry Smith break; 109802834360SBarry Smith } 109902834360SBarry Smith } 1100a2744918SBarry Smith sum = 0; 110102834360SBarry Smith while (k < kend) { 1102d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1103a2744918SBarry Smith sum++; 110402834360SBarry Smith } 1105a2744918SBarry Smith lens[i] = sum; 110602834360SBarry Smith } 110702834360SBarry Smith /* create submatrix */ 1108cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 110908480c60SBarry Smith int n_cols,n_rows; 111008480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 111108480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1112d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 111308480c60SBarry Smith C = *B; 111408480c60SBarry Smith } 111508480c60SBarry Smith else { 111602834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 111708480c60SBarry Smith } 1118db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1119db02288aSLois Curfman McInnes 112002834360SBarry Smith /* loop over rows inserting into submatrix */ 1121db02288aSLois Curfman McInnes a_new = c->a; 1122db02288aSLois Curfman McInnes j_new = c->j; 1123db02288aSLois Curfman McInnes i_new = c->i; 1124db02288aSLois Curfman McInnes i_new[0] = -shift; 112502834360SBarry Smith for (i=0; i<nrows; i++) { 1126a2744918SBarry Smith ii = starts[i]; 1127a2744918SBarry Smith lensi = lens[i]; 1128a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1129a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 113002834360SBarry Smith } 1131a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1132a2744918SBarry Smith a_new += lensi; 1133a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1134a2744918SBarry Smith c->ilen[i] = lensi; 113502834360SBarry Smith } 11360452661fSBarry Smith PetscFree(lens); 113702834360SBarry Smith } 113802834360SBarry Smith else { 113902834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 11400452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 114102834360SBarry Smith ssmap = smap + shift; 114299141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1143cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 114417ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 114502834360SBarry Smith /* determine lens of each row */ 114602834360SBarry Smith for (i=0; i<nrows; i++) { 1147d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 114802834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 114902834360SBarry Smith lens[i] = 0; 115002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1151d8ced48eSBarry Smith if (ssmap[aj[k]]) { 115202834360SBarry Smith lens[i]++; 115302834360SBarry Smith } 115402834360SBarry Smith } 115502834360SBarry Smith } 115617ab2063SBarry Smith /* Create and fill new matrix */ 1157a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 115899141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 115999141d43SSatish Balay 116099141d43SSatish Balay if (c->m != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 116199141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 116299141d43SSatish Balay SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros"); 116399141d43SSatish Balay } 116499141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 116508480c60SBarry Smith C = *B; 116608480c60SBarry Smith } 116708480c60SBarry Smith else { 116802834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 116908480c60SBarry Smith } 117099141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 117117ab2063SBarry Smith for (i=0; i<nrows; i++) { 117299141d43SSatish Balay row = irow[i]; 117317ab2063SBarry Smith nznew = 0; 117499141d43SSatish Balay kstart = ai[row]+shift; 117599141d43SSatish Balay kend = kstart + a->ilen[row]; 117699141d43SSatish Balay mat_i = c->i[i]+shift; 117799141d43SSatish Balay mat_j = c->j + mat_i; 117899141d43SSatish Balay mat_a = c->a + mat_i; 117999141d43SSatish Balay mat_ilen = c->ilen + i; 118017ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 118199141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 118299141d43SSatish Balay *mat_j++ = tcol - (!shift); 118399141d43SSatish Balay *mat_a++ = a->a[k]; 118499141d43SSatish Balay (*mat_ilen)++; 118599141d43SSatish Balay 118617ab2063SBarry Smith } 118717ab2063SBarry Smith } 118817ab2063SBarry Smith } 118902834360SBarry Smith /* Free work space */ 119002834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 119199141d43SSatish Balay PetscFree(smap); PetscFree(lens); 119202834360SBarry Smith } 11936d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11946d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 119517ab2063SBarry Smith 119617ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1197416022c9SBarry Smith *B = C; 119817ab2063SBarry Smith return 0; 119917ab2063SBarry Smith } 120017ab2063SBarry Smith 1201a871dcd8SBarry Smith /* 120263b91edcSBarry Smith note: This can only work for identity for row and col. It would 120363b91edcSBarry Smith be good to check this and otherwise generate an error. 1204a871dcd8SBarry Smith */ 120563b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1206a871dcd8SBarry Smith { 120763b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 120808480c60SBarry Smith int ierr; 120963b91edcSBarry Smith Mat outA; 121063b91edcSBarry Smith 1211a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1212a871dcd8SBarry Smith 121363b91edcSBarry Smith outA = inA; 121463b91edcSBarry Smith inA->factor = FACTOR_LU; 121563b91edcSBarry Smith a->row = row; 121663b91edcSBarry Smith a->col = col; 121763b91edcSBarry Smith 12180452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 121963b91edcSBarry Smith 122008480c60SBarry Smith if (!a->diag) { 122108480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 122263b91edcSBarry Smith } 122363b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1224a871dcd8SBarry Smith return 0; 1225a871dcd8SBarry Smith } 1226a871dcd8SBarry Smith 1227f0b747eeSBarry Smith #include "pinclude/plapack.h" 1228f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1229f0b747eeSBarry Smith { 1230f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1231f0b747eeSBarry Smith int one = 1; 1232f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1233f0b747eeSBarry Smith PLogFlops(a->nz); 1234f0b747eeSBarry Smith return 0; 1235f0b747eeSBarry Smith } 1236f0b747eeSBarry Smith 1237cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1238cddf8d76SBarry Smith Mat **B) 1239cddf8d76SBarry Smith { 1240cddf8d76SBarry Smith int ierr,i; 1241cddf8d76SBarry Smith 1242cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 12430452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1244cddf8d76SBarry Smith } 1245cddf8d76SBarry Smith 1246cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1247905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1248cddf8d76SBarry Smith } 1249cddf8d76SBarry Smith return 0; 1250cddf8d76SBarry Smith } 1251cddf8d76SBarry Smith 12525a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 12535a838052SSatish Balay { 12545a838052SSatish Balay *bs = 1; 12555a838052SSatish Balay return 0; 12565a838052SSatish Balay } 12575a838052SSatish Balay 1258e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 12594dcbc457SBarry Smith { 1260e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 126106763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 12628a047759SSatish Balay int start, end, *ai, *aj; 126306763907SSatish Balay char *table; 12648a047759SSatish Balay shift = a->indexshift; 1265e4d965acSSatish Balay m = a->m; 1266e4d965acSSatish Balay ai = a->i; 12678a047759SSatish Balay aj = a->j+shift; 12688a047759SSatish Balay 1269e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 127006763907SSatish Balay 127106763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 127206763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 127306763907SSatish Balay 1274e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1275e4d965acSSatish Balay /* Initialise the two local arrays */ 1276e4d965acSSatish Balay isz = 0; 127706763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1278e4d965acSSatish Balay 1279e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 12804dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 128177c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1282e4d965acSSatish Balay 1283e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1284e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 128506763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 12864dcbc457SBarry Smith } 128706763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 128806763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1289e4d965acSSatish Balay 129004a348a9SBarry Smith k = 0; 129104a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 129204a348a9SBarry Smith n = isz; 129306763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1294e4d965acSSatish Balay row = nidx[k]; 1295e4d965acSSatish Balay start = ai[row]; 1296e4d965acSSatish Balay end = ai[row+1]; 129704a348a9SBarry Smith for ( l = start; l<end ; l++){ 12988a047759SSatish Balay val = aj[l] + shift; 129906763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1300e4d965acSSatish Balay } 1301e4d965acSSatish Balay } 1302e4d965acSSatish Balay } 1303537820f0SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1304e4d965acSSatish Balay } 130504a348a9SBarry Smith PetscFree(table); 130606763907SSatish Balay PetscFree(nidx); 1307e4d965acSSatish Balay return 0; 13084dcbc457SBarry Smith } 130917ab2063SBarry Smith 1310682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1311682d7d0cSBarry Smith { 1312682d7d0cSBarry Smith static int called = 0; 1313682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1314682d7d0cSBarry Smith 1315682d7d0cSBarry Smith if (called) return 0; else called = 1; 131677c4ece6SBarry Smith PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 131777c4ece6SBarry Smith PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>\n"); 131877c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 131977c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_no_inode - Do not use inodes\n"); 132077c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1321682d7d0cSBarry Smith #if defined(HAVE_ESSL) 132277c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1323682d7d0cSBarry Smith #endif 1324682d7d0cSBarry Smith return 0; 1325682d7d0cSBarry Smith } 1326205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1327682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 132817ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 132917ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1330416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1331416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 133217ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 133317ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 133417ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 133517ab2063SBarry Smith MatRelax_SeqAIJ, 133617ab2063SBarry Smith MatTranspose_SeqAIJ, 13377264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1338f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 133917ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 134017ab2063SBarry Smith MatCompress_SeqAIJ, 134117ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 134217ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 134317ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 134417ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 134517ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 13463d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1347cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 13487eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1349682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1350f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 13515a838052SSatish Balay MatScale_SeqAIJ,0,0, 1352*6945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 1353*6945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 1354*6945ee14SBarry Smith MatGetIJ_SeqAIJ, 1355*6945ee14SBarry Smith MatRestoreIJ_SeqAIJ}; 135617ab2063SBarry Smith 135717ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 135817ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 135917ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 136017ab2063SBarry Smith 136117ab2063SBarry Smith /*@C 1362682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 13630d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 13646e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 13652bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 13662bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 136717ab2063SBarry Smith 136817ab2063SBarry Smith Input Parameters: 136917ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 137017ab2063SBarry Smith . m - number of rows 137117ab2063SBarry Smith . n - number of columns 137217ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 13732bd5e0b2SLois Curfman McInnes . nzz - array containing the number of nonzeros in the various rows 13742bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 137517ab2063SBarry Smith 137617ab2063SBarry Smith Output Parameter: 1377416022c9SBarry Smith . A - the matrix 137817ab2063SBarry Smith 137917ab2063SBarry Smith Notes: 138017ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 138117ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 13820002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 138344cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 138417ab2063SBarry Smith 138517ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1386a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 13873d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 13883d323bbdSBarry Smith will get TERRIBLE performance, see the users' manual chapter on 13890d15e28bSLois Curfman McInnes matrices and the file $(PETSC_DIR)/Performance. 139017ab2063SBarry Smith 1391682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1392682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1393682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 13946c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 13956c7ebb05SLois Curfman McInnes 13966c7ebb05SLois Curfman McInnes Options Database Keys: 13976c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 13986e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 13996e62573dSLois Curfman McInnes $ (max limit=5) 14006e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 14016e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 14026e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 140317ab2063SBarry Smith 140417ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 140517ab2063SBarry Smith @*/ 1406416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 140717ab2063SBarry Smith { 1408416022c9SBarry Smith Mat B; 1409416022c9SBarry Smith Mat_SeqAIJ *b; 1410*6945ee14SBarry Smith int i, len, ierr, flg,size; 1411*6945ee14SBarry Smith 1412*6945ee14SBarry Smith MPI_Comm_size(comm,&size); 1413*6945ee14SBarry Smith if (size > 1) SETERRQ(1,"MatCreateSeqAIJ:Comm must be of size 1"); 1414d5d45c9bSBarry Smith 1415416022c9SBarry Smith *A = 0; 14160452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1417416022c9SBarry Smith PLogObjectCreate(B); 14180452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 141944cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1420416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1421416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1422416022c9SBarry Smith B->view = MatView_SeqAIJ; 1423416022c9SBarry Smith B->factor = 0; 1424416022c9SBarry Smith B->lupivotthreshold = 1.0; 14257a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 142669957df2SSatish Balay &flg); CHKERRQ(ierr); 14277a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 14287a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 14297a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1430416022c9SBarry Smith b->row = 0; 1431416022c9SBarry Smith b->col = 0; 1432416022c9SBarry Smith b->indexshift = 0; 1433b810aeb4SBarry Smith b->reallocs = 0; 143469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 143569957df2SSatish Balay if (flg) b->indexshift = -1; 143617ab2063SBarry Smith 143744cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 143844cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 14390452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1440b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 14417b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 14427b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1443416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 144417ab2063SBarry Smith nz = nz*m; 144517ab2063SBarry Smith } 144617ab2063SBarry Smith else { 144717ab2063SBarry Smith nz = 0; 1448416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 144917ab2063SBarry Smith } 145017ab2063SBarry Smith 145117ab2063SBarry Smith /* allocate the matrix space */ 1452416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 14530452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1454416022c9SBarry Smith b->j = (int *) (b->a + nz); 1455cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1456416022c9SBarry Smith b->i = b->j + nz; 1457416022c9SBarry Smith b->singlemalloc = 1; 145817ab2063SBarry Smith 1459416022c9SBarry Smith b->i[0] = -b->indexshift; 146017ab2063SBarry Smith for (i=1; i<m+1; i++) { 1461416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 146217ab2063SBarry Smith } 146317ab2063SBarry Smith 1464416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 14650452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1466416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1467416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 146817ab2063SBarry Smith 1469416022c9SBarry Smith b->nz = 0; 1470416022c9SBarry Smith b->maxnz = nz; 1471416022c9SBarry Smith b->sorted = 0; 1472416022c9SBarry Smith b->roworiented = 1; 1473416022c9SBarry Smith b->nonew = 0; 1474416022c9SBarry Smith b->diag = 0; 1475416022c9SBarry Smith b->solve_work = 0; 147671bd300dSLois Curfman McInnes b->spptr = 0; 1477754ec7b1SSatish Balay b->inode.node_count = 0; 1478754ec7b1SSatish Balay b->inode.size = 0; 14796c7ebb05SLois Curfman McInnes b->inode.limit = 5; 14806c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 14814e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 148217ab2063SBarry Smith 1483416022c9SBarry Smith *A = B; 14844e220ebcSLois Curfman McInnes 14854b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 14864b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 148769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 148869957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 14894b14c69eSBarry Smith #endif 149069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 149169957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 149269957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 149369957df2SSatish Balay if (flg) { 1494416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1495416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 149617ab2063SBarry Smith } 149769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 149869957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 149917ab2063SBarry Smith return 0; 150017ab2063SBarry Smith } 150117ab2063SBarry Smith 15023d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 150317ab2063SBarry Smith { 1504416022c9SBarry Smith Mat C; 1505416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 150608480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 150717ab2063SBarry Smith 15084043dd9cSLois Curfman McInnes *B = 0; 15090452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1510416022c9SBarry Smith PLogObjectCreate(C); 15110452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 151241c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1513416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1514416022c9SBarry Smith C->view = MatView_SeqAIJ; 1515416022c9SBarry Smith C->factor = A->factor; 1516416022c9SBarry Smith c->row = 0; 1517416022c9SBarry Smith c->col = 0; 1518416022c9SBarry Smith c->indexshift = shift; 1519c456f294SBarry Smith C->assembled = PETSC_TRUE; 152017ab2063SBarry Smith 152144cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 152244cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 152344cd7ae7SLois Curfman McInnes C->M = a->m; 152444cd7ae7SLois Curfman McInnes C->N = a->n; 152517ab2063SBarry Smith 15260452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 15270452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 152817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1529416022c9SBarry Smith c->imax[i] = a->imax[i]; 1530416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 153117ab2063SBarry Smith } 153217ab2063SBarry Smith 153317ab2063SBarry Smith /* allocate the matrix space */ 1534416022c9SBarry Smith c->singlemalloc = 1; 1535416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 15360452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1537416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1538416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1539416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 154017ab2063SBarry Smith if (m > 0) { 1541416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 154208480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1543416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 154417ab2063SBarry Smith } 154508480c60SBarry Smith } 154617ab2063SBarry Smith 1547416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1548416022c9SBarry Smith c->sorted = a->sorted; 1549416022c9SBarry Smith c->roworiented = a->roworiented; 1550416022c9SBarry Smith c->nonew = a->nonew; 15517a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 155217ab2063SBarry Smith 1553416022c9SBarry Smith if (a->diag) { 15540452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1555416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 155617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1557416022c9SBarry Smith c->diag[i] = a->diag[i]; 155817ab2063SBarry Smith } 155917ab2063SBarry Smith } 1560416022c9SBarry Smith else c->diag = 0; 15616c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 15626c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1563754ec7b1SSatish Balay if (a->inode.size){ 1564754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1565754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1566754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1567754ec7b1SSatish Balay } else { 1568754ec7b1SSatish Balay c->inode.size = 0; 1569754ec7b1SSatish Balay c->inode.node_count = 0; 1570754ec7b1SSatish Balay } 1571416022c9SBarry Smith c->nz = a->nz; 1572416022c9SBarry Smith c->maxnz = a->maxnz; 1573416022c9SBarry Smith c->solve_work = 0; 157476dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1575754ec7b1SSatish Balay 1576416022c9SBarry Smith *B = C; 157717ab2063SBarry Smith return 0; 157817ab2063SBarry Smith } 157917ab2063SBarry Smith 158019bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 158117ab2063SBarry Smith { 1582416022c9SBarry Smith Mat_SeqAIJ *a; 1583416022c9SBarry Smith Mat B; 158417699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1585bcd2baecSBarry Smith MPI_Comm comm; 158617ab2063SBarry Smith 158719bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 158817699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 158917699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 159090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 159177c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 159248d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 159317ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 159417ab2063SBarry Smith 159517ab2063SBarry Smith /* read in row lengths */ 15960452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 159777c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 159817ab2063SBarry Smith 159917ab2063SBarry Smith /* create our matrix */ 1600416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1601416022c9SBarry Smith B = *A; 1602416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1603416022c9SBarry Smith shift = a->indexshift; 160417ab2063SBarry Smith 160517ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 160677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 160717ab2063SBarry Smith if (shift) { 160817ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1609416022c9SBarry Smith a->j[i] += 1; 161017ab2063SBarry Smith } 161117ab2063SBarry Smith } 161217ab2063SBarry Smith 161317ab2063SBarry Smith /* read in nonzero values */ 161477c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 161517ab2063SBarry Smith 161617ab2063SBarry Smith /* set matrix "i" values */ 1617416022c9SBarry Smith a->i[0] = -shift; 161817ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1619416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1620416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 162117ab2063SBarry Smith } 16220452661fSBarry Smith PetscFree(rowlengths); 162317ab2063SBarry Smith 16246d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 16256d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 162617ab2063SBarry Smith return 0; 162717ab2063SBarry Smith } 162817ab2063SBarry Smith 1629686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 16307264ac53SSatish Balay { 16317264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 16327264ac53SSatish Balay 1633bcd2baecSBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type"); 16347264ac53SSatish Balay 16357264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 16367264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1637bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 163877c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1639bcd2baecSBarry Smith } 16407264ac53SSatish Balay 16417264ac53SSatish Balay /* if the a->i are the same */ 1642bcd2baecSBarry Smith if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 164377c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 16447264ac53SSatish Balay } 16457264ac53SSatish Balay 16467264ac53SSatish Balay /* if a->j are the same */ 1647bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 164877c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1649bcd2baecSBarry Smith } 1650bcd2baecSBarry Smith 1651bcd2baecSBarry Smith /* if a->a are the same */ 165219bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 165377c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 16547264ac53SSatish Balay } 165577c4ece6SBarry Smith *flg = PETSC_TRUE; 16567264ac53SSatish Balay return 0; 16577264ac53SSatish Balay 16587264ac53SSatish Balay } 1659