16d84be18SBarry Smith 217ab2063SBarry Smith #ifndef lint 3*3d323bbdSBarry Smith static char vcid[] = "$Id: aij.c,v 1.172 1996/06/26 18:11:42 curfman Exp bsmith $"; 417ab2063SBarry Smith #endif 517ab2063SBarry Smith 6d5d45c9bSBarry Smith /* 7d5d45c9bSBarry Smith Defines the basic matrix operations for the AIJ (compressed row) 8d5d45c9bSBarry Smith matrix storage format. 9d5d45c9bSBarry Smith */ 1017ab2063SBarry Smith #include "aij.h" 11f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 12f5eb4b81SSatish Balay #include "src/inline/spops.h" 13e4d965acSSatish Balay #include "petsc.h" 14f5eb4b81SSatish Balay #include "src/inline/bitarray.h" 1517ab2063SBarry Smith 16bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 1717ab2063SBarry Smith 18416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1917ab2063SBarry Smith { 20416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 21bcd2baecSBarry Smith int ierr, *ia, *ja,n,*idx,i,oshift,ishift; 2217ab2063SBarry Smith 23a2744918SBarry Smith /* 24a2744918SBarry Smith this is tacky: In the future when we have written special factorization 25a2744918SBarry Smith and solve routines for the identity permutation we should use a 26a2744918SBarry Smith stride index set instead of the general one. 27a2744918SBarry Smith */ 28a2744918SBarry Smith if (type == ORDER_NATURAL) { 29a2744918SBarry Smith n = a->n; 30a2744918SBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 31a2744918SBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 32a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 33a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 34a2744918SBarry Smith PetscFree(idx); 35a2744918SBarry Smith ISSetPermutation(*rperm); 36a2744918SBarry Smith ISSetPermutation(*cperm); 37a2744918SBarry Smith ISSetIdentity(*rperm); 38a2744918SBarry Smith ISSetIdentity(*cperm); 39a2744918SBarry Smith return 0; 40a2744918SBarry Smith } 41a2744918SBarry Smith 42bcd2baecSBarry Smith MatReorderingRegisterAll(); 43bcd2baecSBarry Smith ishift = a->indexshift; 44bcd2baecSBarry Smith oshift = -MatReorderingIndexShift[(int)type]; 45bcd2baecSBarry Smith if (MatReorderingRequiresSymmetric[(int)type]) { 46bcd2baecSBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja); 47bcd2baecSBarry Smith CHKERRQ(ierr); 48416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 490452661fSBarry Smith PetscFree(ia); PetscFree(ja); 50bcd2baecSBarry Smith } else { 51bcd2baecSBarry Smith if (ishift == oshift) { 52bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 53bcd2baecSBarry Smith } 54bcd2baecSBarry Smith else if (ishift == -1) { 55bcd2baecSBarry Smith /* temporarily subtract 1 from i and j indices */ 56bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 57bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 58bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 59bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 60bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 61bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 62bcd2baecSBarry Smith } else { 63bcd2baecSBarry Smith /* temporarily add 1 to i and j indices */ 64bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 65bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 66bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 67bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 68bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 69bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 70bcd2baecSBarry Smith } 71bcd2baecSBarry Smith } 7217ab2063SBarry Smith return 0; 7317ab2063SBarry Smith } 7417ab2063SBarry Smith 75227d817aSBarry Smith #define CHUNKSIZE 15 7617ab2063SBarry Smith 7717ab2063SBarry Smith /* This version has row oriented v */ 78416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 7917ab2063SBarry Smith { 80416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 81416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 824b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 83d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 84416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 8517ab2063SBarry Smith 8617ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 87416022c9SBarry Smith row = im[k]; 8817ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 89416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 9017ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 9117ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 92416022c9SBarry Smith low = 0; 9317ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 94416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 95416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 964b0e389bSBarry Smith col = in[l] - shift; 974b0e389bSBarry Smith if (roworiented) { 984b0e389bSBarry Smith value = *v++; 994b0e389bSBarry Smith } 1004b0e389bSBarry Smith else { 1014b0e389bSBarry Smith value = v[k + l*m]; 1024b0e389bSBarry Smith } 103416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 104416022c9SBarry Smith while (high-low > 5) { 105416022c9SBarry Smith t = (low+high)/2; 106416022c9SBarry Smith if (rp[t] > col) high = t; 107416022c9SBarry Smith else low = t; 10817ab2063SBarry Smith } 109416022c9SBarry Smith for ( i=low; i<high; i++ ) { 11017ab2063SBarry Smith if (rp[i] > col) break; 11117ab2063SBarry Smith if (rp[i] == col) { 112416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 11317ab2063SBarry Smith else ap[i] = value; 11417ab2063SBarry Smith goto noinsert; 11517ab2063SBarry Smith } 11617ab2063SBarry Smith } 11717ab2063SBarry Smith if (nonew) goto noinsert; 11817ab2063SBarry Smith if (nrow >= rmax) { 11917ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 120416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 12117ab2063SBarry Smith Scalar *new_a; 12217ab2063SBarry Smith 12317ab2063SBarry Smith /* malloc new storage space */ 124416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1250452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 12617ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 12717ab2063SBarry Smith new_i = new_j + new_nz; 12817ab2063SBarry Smith 12917ab2063SBarry Smith /* copy over old data into new slots */ 13017ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 131416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 132416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 133416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 134416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 13517ab2063SBarry Smith len*sizeof(int)); 136416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 137416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 13817ab2063SBarry Smith len*sizeof(Scalar)); 13917ab2063SBarry Smith /* free up old matrix storage */ 1400452661fSBarry Smith PetscFree(a->a); 1410452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 142416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 143416022c9SBarry Smith a->singlemalloc = 1; 14417ab2063SBarry Smith 14517ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 146416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 147416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 148416022c9SBarry Smith a->maxnz += CHUNKSIZE; 149b810aeb4SBarry Smith a->reallocs++; 15017ab2063SBarry Smith } 151416022c9SBarry Smith N = nrow++ - 1; a->nz++; 152416022c9SBarry Smith /* shift up all the later entries in this row */ 153416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 15417ab2063SBarry Smith rp[ii+1] = rp[ii]; 15517ab2063SBarry Smith ap[ii+1] = ap[ii]; 15617ab2063SBarry Smith } 15717ab2063SBarry Smith rp[i] = col; 15817ab2063SBarry Smith ap[i] = value; 15917ab2063SBarry Smith noinsert:; 160416022c9SBarry Smith low = i + 1; 16117ab2063SBarry Smith } 16217ab2063SBarry Smith ailen[row] = nrow; 16317ab2063SBarry Smith } 16417ab2063SBarry Smith return 0; 16517ab2063SBarry Smith } 16617ab2063SBarry Smith 1677eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 1687eb43aa7SLois Curfman McInnes { 1697eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 170b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 1717eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 1727eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 1737eb43aa7SLois Curfman McInnes 1747eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 1757eb43aa7SLois Curfman McInnes row = im[k]; 1767eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 1777eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 1787eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 1797eb43aa7SLois Curfman McInnes nrow = ailen[row]; 1807eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 1817eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 1827eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 1837eb43aa7SLois Curfman McInnes col = in[l] - shift; 1847eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 1857eb43aa7SLois Curfman McInnes while (high-low > 5) { 1867eb43aa7SLois Curfman McInnes t = (low+high)/2; 1877eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 1887eb43aa7SLois Curfman McInnes else low = t; 1897eb43aa7SLois Curfman McInnes } 1907eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 1917eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 1927eb43aa7SLois Curfman McInnes if (rp[i] == col) { 193b49de8d1SLois Curfman McInnes *v++ = ap[i]; 1947eb43aa7SLois Curfman McInnes goto finished; 1957eb43aa7SLois Curfman McInnes } 1967eb43aa7SLois Curfman McInnes } 197b49de8d1SLois Curfman McInnes *v++ = zero; 1987eb43aa7SLois Curfman McInnes finished:; 1997eb43aa7SLois Curfman McInnes } 2007eb43aa7SLois Curfman McInnes } 2017eb43aa7SLois Curfman McInnes return 0; 2027eb43aa7SLois Curfman McInnes } 2037eb43aa7SLois Curfman McInnes 20417ab2063SBarry Smith #include "draw.h" 20517ab2063SBarry Smith #include "pinclude/pviewer.h" 20677c4ece6SBarry Smith #include "sys.h" 20717ab2063SBarry Smith 208416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 20917ab2063SBarry Smith { 210416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 211416022c9SBarry Smith int i, fd, *col_lens, ierr; 21217ab2063SBarry Smith 21390ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2140452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 215416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 216416022c9SBarry Smith col_lens[1] = a->m; 217416022c9SBarry Smith col_lens[2] = a->n; 218416022c9SBarry Smith col_lens[3] = a->nz; 219416022c9SBarry Smith 220416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 221416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 222416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 22317ab2063SBarry Smith } 22477c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 2250452661fSBarry Smith PetscFree(col_lens); 226416022c9SBarry Smith 227416022c9SBarry Smith /* store column indices (zero start index) */ 228416022c9SBarry Smith if (a->indexshift) { 229416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 23017ab2063SBarry Smith } 23177c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 232416022c9SBarry Smith if (a->indexshift) { 233416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 23417ab2063SBarry Smith } 235416022c9SBarry Smith 236416022c9SBarry Smith /* store nonzero values */ 23777c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 23817ab2063SBarry Smith return 0; 23917ab2063SBarry Smith } 240416022c9SBarry Smith 241416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 242416022c9SBarry Smith { 243416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 244496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 24517ab2063SBarry Smith FILE *fd; 24617ab2063SBarry Smith char *outputname; 24717ab2063SBarry Smith 24890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 249416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 25090ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 25190ace30eSBarry Smith if (format == ASCII_FORMAT_INFO) { 25295e01e2fSLois Curfman McInnes return 0; 25395e01e2fSLois Curfman McInnes } 25490ace30eSBarry Smith else if (format == ASCII_FORMAT_INFO_DETAILED) { 255496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 256496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 257496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 25895e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 25995e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 26017ab2063SBarry Smith } 26190ace30eSBarry Smith else if (format == ASCII_FORMAT_MATLAB) { 26217ab2063SBarry Smith int nz, nzalloc, mem; 263416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 264416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 26517ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 26617ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 26717ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 26817ab2063SBarry Smith 26917ab2063SBarry Smith for (i=0; i<m; i++) { 270416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 27117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 2727a743949SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]), 273416022c9SBarry Smith imag(a->a[j])); 27417ab2063SBarry Smith #else 2757a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 27617ab2063SBarry Smith #endif 27717ab2063SBarry Smith } 27817ab2063SBarry Smith } 27917ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 28017ab2063SBarry Smith } 28144cd7ae7SLois Curfman McInnes else if (format == ASCII_FORMAT_COMMON) { 28244cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 28344cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 28444cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 28544cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 28644cd7ae7SLois Curfman McInnes if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0) 28744cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 28844cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 28944cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 29044cd7ae7SLois Curfman McInnes #else 29144cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 29244cd7ae7SLois Curfman McInnes #endif 29344cd7ae7SLois Curfman McInnes } 29444cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 29544cd7ae7SLois Curfman McInnes } 29644cd7ae7SLois Curfman McInnes } 29717ab2063SBarry Smith else { 29817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 29917ab2063SBarry Smith fprintf(fd,"row %d:",i); 300416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 30117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 302416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 303416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 30417ab2063SBarry Smith } 30517ab2063SBarry Smith else { 306416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 30717ab2063SBarry Smith } 30817ab2063SBarry Smith #else 309416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 31017ab2063SBarry Smith #endif 31117ab2063SBarry Smith } 31217ab2063SBarry Smith fprintf(fd,"\n"); 31317ab2063SBarry Smith } 31417ab2063SBarry Smith } 31517ab2063SBarry Smith fflush(fd); 316416022c9SBarry Smith return 0; 317416022c9SBarry Smith } 318416022c9SBarry Smith 319416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 320416022c9SBarry Smith { 321416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 322cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 323cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 324bcd2baecSBarry Smith Draw draw; 325cddf8d76SBarry Smith DrawButton button; 32619bcc07fSBarry Smith PetscTruth isnull; 327cddf8d76SBarry Smith 328bcd2baecSBarry Smith ViewerDrawGetDraw(viewer,&draw); 32919bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 33019bcc07fSBarry Smith 331416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 332416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 333416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 334416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 335cddf8d76SBarry Smith color = DRAW_BLUE; 336416022c9SBarry Smith for ( i=0; i<m; i++ ) { 337cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 338416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 339cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 340cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 341cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 342cddf8d76SBarry Smith #else 343cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 344cddf8d76SBarry Smith #endif 345cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 346cddf8d76SBarry Smith } 347cddf8d76SBarry Smith } 348cddf8d76SBarry Smith color = DRAW_CYAN; 349cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 350cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 351cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 352cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 353cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 354cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 355cddf8d76SBarry Smith } 356cddf8d76SBarry Smith } 357cddf8d76SBarry Smith color = DRAW_RED; 358cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 359cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 360cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 361cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 362cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 363cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 364cddf8d76SBarry Smith #else 365cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 366cddf8d76SBarry Smith #endif 367cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 368416022c9SBarry Smith } 369416022c9SBarry Smith } 370416022c9SBarry Smith DrawFlush(draw); 371cddf8d76SBarry Smith DrawGetPause(draw,&pause); 372cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 373cddf8d76SBarry Smith 374cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 375cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 376cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 377cddf8d76SBarry Smith DrawClear(draw); 378cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 379cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 380cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 381cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 382cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 383cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 384cddf8d76SBarry Smith w *= scale; h *= scale; 385cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 386cddf8d76SBarry Smith color = DRAW_BLUE; 387cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 388cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 389cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 390cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 391cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 392cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 393cddf8d76SBarry Smith #else 394cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 395cddf8d76SBarry Smith #endif 396cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 397cddf8d76SBarry Smith } 398cddf8d76SBarry Smith } 399cddf8d76SBarry Smith color = DRAW_CYAN; 400cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 401cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 402cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 403cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 404cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 405cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 406cddf8d76SBarry Smith } 407cddf8d76SBarry Smith } 408cddf8d76SBarry Smith color = DRAW_RED; 409cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 410cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 411cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 412cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 413cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 414cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 415cddf8d76SBarry Smith #else 416cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 417cddf8d76SBarry Smith #endif 418cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 419cddf8d76SBarry Smith } 420cddf8d76SBarry Smith } 421cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 422cddf8d76SBarry Smith } 423416022c9SBarry Smith return 0; 424416022c9SBarry Smith } 425416022c9SBarry Smith 426416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 427416022c9SBarry Smith { 428416022c9SBarry Smith Mat A = (Mat) obj; 429416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 430bcd2baecSBarry Smith ViewerType vtype; 431bcd2baecSBarry Smith int ierr; 432416022c9SBarry Smith 433416022c9SBarry Smith if (!viewer) { 434bcd2baecSBarry Smith viewer = STDOUT_VIEWER_SELF; 435416022c9SBarry Smith } 436bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 437bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 438416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 439416022c9SBarry Smith } 440bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 441416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 442416022c9SBarry Smith } 443bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 444416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 445416022c9SBarry Smith } 446bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 447bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 44817ab2063SBarry Smith } 44917ab2063SBarry Smith return 0; 45017ab2063SBarry Smith } 45119bcc07fSBarry Smith 452c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 453416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 45417ab2063SBarry Smith { 455416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 45641c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 457416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 458416022c9SBarry Smith Scalar *aa = a->a, *ap; 45917ab2063SBarry Smith 46017ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 46117ab2063SBarry Smith 46217ab2063SBarry Smith for ( i=1; i<m; i++ ) { 463416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 46417ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 46517ab2063SBarry Smith if (fshift) { 466416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 46717ab2063SBarry Smith N = ailen[i]; 46817ab2063SBarry Smith for ( j=0; j<N; j++ ) { 46917ab2063SBarry Smith ip[j-fshift] = ip[j]; 47017ab2063SBarry Smith ap[j-fshift] = ap[j]; 47117ab2063SBarry Smith } 47217ab2063SBarry Smith } 47317ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 47417ab2063SBarry Smith } 47517ab2063SBarry Smith if (m) { 47617ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 47717ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 47817ab2063SBarry Smith } 47917ab2063SBarry Smith /* reset ilen and imax for each row */ 48017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 48117ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 48217ab2063SBarry Smith } 483416022c9SBarry Smith a->nz = ai[m] + shift; 48417ab2063SBarry Smith 48517ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 486416022c9SBarry Smith if (fshift && a->diag) { 4870452661fSBarry Smith PetscFree(a->diag); 488416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 489416022c9SBarry Smith a->diag = 0; 49017ab2063SBarry Smith } 49194a424c1SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n", 492b810aeb4SBarry Smith fshift,a->nz,m); 49394a424c1SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n", 494b810aeb4SBarry Smith a->reallocs); 49576dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 49641c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 49717ab2063SBarry Smith return 0; 49817ab2063SBarry Smith } 49917ab2063SBarry Smith 500416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 50117ab2063SBarry Smith { 502416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 503cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 50417ab2063SBarry Smith return 0; 50517ab2063SBarry Smith } 506416022c9SBarry Smith 50717ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 50817ab2063SBarry Smith { 509416022c9SBarry Smith Mat A = (Mat) obj; 510416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 511d5d45c9bSBarry Smith 51217ab2063SBarry Smith #if defined(PETSC_LOG) 513416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 51417ab2063SBarry Smith #endif 5150452661fSBarry Smith PetscFree(a->a); 5160452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 5170452661fSBarry Smith if (a->diag) PetscFree(a->diag); 5180452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 5190452661fSBarry Smith if (a->imax) PetscFree(a->imax); 5200452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 52176dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 5220452661fSBarry Smith PetscFree(a); 523f2655603SLois Curfman McInnes PLogObjectDestroy(A); 524f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 52517ab2063SBarry Smith return 0; 52617ab2063SBarry Smith } 52717ab2063SBarry Smith 528416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 52917ab2063SBarry Smith { 53017ab2063SBarry Smith return 0; 53117ab2063SBarry Smith } 53217ab2063SBarry Smith 533416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 53417ab2063SBarry Smith { 535416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 536416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 5374b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 538416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 539416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 540416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 541e2f28af5SBarry Smith else if (op == ROWS_SORTED || 542e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 543e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 544e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 54594a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 546df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 5474d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 5486c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_1) a->inode.limit = 1; 5496c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_2) a->inode.limit = 2; 5506c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_3) a->inode.limit = 3; 5516c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_4) a->inode.limit = 4; 5526c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_5) a->inode.limit = 5; 553e2f28af5SBarry Smith else 5544d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 55517ab2063SBarry Smith return 0; 55617ab2063SBarry Smith } 55717ab2063SBarry Smith 558416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 55917ab2063SBarry Smith { 560416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 561416022c9SBarry Smith int i,j, n,shift = a->indexshift; 56217ab2063SBarry Smith Scalar *x, zero = 0.0; 56317ab2063SBarry Smith 56417ab2063SBarry Smith VecSet(&zero,v); 56517ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 566416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 567416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 568416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 569416022c9SBarry Smith if (a->j[j]+shift == i) { 570416022c9SBarry Smith x[i] = a->a[j]; 57117ab2063SBarry Smith break; 57217ab2063SBarry Smith } 57317ab2063SBarry Smith } 57417ab2063SBarry Smith } 57517ab2063SBarry Smith return 0; 57617ab2063SBarry Smith } 57717ab2063SBarry Smith 57817ab2063SBarry Smith /* -------------------------------------------------------*/ 57917ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 58017ab2063SBarry Smith /* -------------------------------------------------------*/ 58144cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 58217ab2063SBarry Smith { 583416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 58417ab2063SBarry Smith Scalar *x, *y, *v, alpha; 585416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 58617ab2063SBarry Smith 58717ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 588cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 58917ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 59017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 591416022c9SBarry Smith idx = a->j + a->i[i] + shift; 592416022c9SBarry Smith v = a->a + a->i[i] + shift; 593416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 59417ab2063SBarry Smith alpha = x[i]; 59517ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 59617ab2063SBarry Smith } 597416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 59817ab2063SBarry Smith return 0; 59917ab2063SBarry Smith } 600d5d45c9bSBarry Smith 60144cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 60217ab2063SBarry Smith { 603416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 60417ab2063SBarry Smith Scalar *x, *y, *v, alpha; 605416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 60617ab2063SBarry Smith 60717ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 60817ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 60917ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 61017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 611416022c9SBarry Smith idx = a->j + a->i[i] + shift; 612416022c9SBarry Smith v = a->a + a->i[i] + shift; 613416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 61417ab2063SBarry Smith alpha = x[i]; 61517ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 61617ab2063SBarry Smith } 61717ab2063SBarry Smith return 0; 61817ab2063SBarry Smith } 61917ab2063SBarry Smith 62044cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 62117ab2063SBarry Smith { 622416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 62317ab2063SBarry Smith Scalar *x, *y, *v, sum; 624416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 62517ab2063SBarry Smith 62617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 62717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 628416022c9SBarry Smith idx = a->j; 629416022c9SBarry Smith v = a->a; 630416022c9SBarry Smith ii = a->i; 63117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 632416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 63317ab2063SBarry Smith sum = 0.0; 63417ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 63517ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 636416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 63717ab2063SBarry Smith y[i] = sum; 63817ab2063SBarry Smith } 639416022c9SBarry Smith PLogFlops(2*a->nz - m); 64017ab2063SBarry Smith return 0; 64117ab2063SBarry Smith } 64217ab2063SBarry Smith 64344cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 64417ab2063SBarry Smith { 645416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 64617ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 647cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 64817ab2063SBarry Smith 64917ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 65017ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 651cddf8d76SBarry Smith idx = a->j; 652cddf8d76SBarry Smith v = a->a; 653cddf8d76SBarry Smith ii = a->i; 65417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 655cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 65617ab2063SBarry Smith sum = y[i]; 657cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 658cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 65917ab2063SBarry Smith z[i] = sum; 66017ab2063SBarry Smith } 661416022c9SBarry Smith PLogFlops(2*a->nz); 66217ab2063SBarry Smith return 0; 66317ab2063SBarry Smith } 66417ab2063SBarry Smith 66517ab2063SBarry Smith /* 66617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 66717ab2063SBarry Smith */ 66817ab2063SBarry Smith 669416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 67017ab2063SBarry Smith { 671416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 672416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 67317ab2063SBarry Smith 6740452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 675416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 676416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 677416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 678416022c9SBarry Smith if (a->j[j]+shift == i) { 67917ab2063SBarry Smith diag[i] = j - shift; 68017ab2063SBarry Smith break; 68117ab2063SBarry Smith } 68217ab2063SBarry Smith } 68317ab2063SBarry Smith } 684416022c9SBarry Smith a->diag = diag; 68517ab2063SBarry Smith return 0; 68617ab2063SBarry Smith } 68717ab2063SBarry Smith 68844cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 68917ab2063SBarry Smith double fshift,int its,Vec xx) 69017ab2063SBarry Smith { 691416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 692416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 693d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 69417ab2063SBarry Smith 69517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 696416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 697416022c9SBarry Smith diag = a->diag; 698416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 69917ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 70017ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 70117ab2063SBarry Smith bs = b + shift; 70217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 703416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 704416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 705416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 706416022c9SBarry Smith v = a->a + diag[i] + (!shift); 70717ab2063SBarry Smith sum = b[i]*d/omega; 70817ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 70917ab2063SBarry Smith x[i] = sum; 71017ab2063SBarry Smith } 71117ab2063SBarry Smith return 0; 71217ab2063SBarry Smith } 71317ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 71417ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 71517ab2063SBarry Smith } 716416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 71717ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 71817ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 71917ab2063SBarry Smith 72017ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 72117ab2063SBarry Smith 72217ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 72317ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 72417ab2063SBarry Smith is the relaxation factor. 72517ab2063SBarry Smith */ 7260452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 72717ab2063SBarry Smith scale = (2.0/omega) - 1.0; 72817ab2063SBarry Smith 72917ab2063SBarry Smith /* x = (E + U)^{-1} b */ 73017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 731416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 732416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 733416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 734416022c9SBarry Smith v = a->a + diag[i] + (!shift); 73517ab2063SBarry Smith sum = b[i]; 73617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 73717ab2063SBarry Smith x[i] = omega*(sum/d); 73817ab2063SBarry Smith } 73917ab2063SBarry Smith 74017ab2063SBarry Smith /* t = b - (2*E - D)x */ 741416022c9SBarry Smith v = a->a; 74217ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 74317ab2063SBarry Smith 74417ab2063SBarry Smith /* t = (E + L)^{-1}t */ 745416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 746416022c9SBarry Smith diag = a->diag; 74717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 748416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 749416022c9SBarry Smith n = diag[i] - a->i[i]; 750416022c9SBarry Smith idx = a->j + a->i[i] + shift; 751416022c9SBarry Smith v = a->a + a->i[i] + shift; 75217ab2063SBarry Smith sum = t[i]; 75317ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 75417ab2063SBarry Smith t[i] = omega*(sum/d); 75517ab2063SBarry Smith } 75617ab2063SBarry Smith 75717ab2063SBarry Smith /* x = x + t */ 75817ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7590452661fSBarry Smith PetscFree(t); 76017ab2063SBarry Smith return 0; 76117ab2063SBarry Smith } 76217ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 76317ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 76417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 765416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 766416022c9SBarry Smith n = diag[i] - a->i[i]; 767416022c9SBarry Smith idx = a->j + a->i[i] + shift; 768416022c9SBarry Smith v = a->a + a->i[i] + shift; 76917ab2063SBarry Smith sum = b[i]; 77017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 77117ab2063SBarry Smith x[i] = omega*(sum/d); 77217ab2063SBarry Smith } 77317ab2063SBarry Smith xb = x; 77417ab2063SBarry Smith } 77517ab2063SBarry Smith else xb = b; 77617ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 77717ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 77817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 779416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 78017ab2063SBarry Smith } 78117ab2063SBarry Smith } 78217ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 78317ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 784416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 785416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 786416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 787416022c9SBarry Smith v = a->a + diag[i] + (!shift); 78817ab2063SBarry Smith sum = xb[i]; 78917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 79017ab2063SBarry Smith x[i] = omega*(sum/d); 79117ab2063SBarry Smith } 79217ab2063SBarry Smith } 79317ab2063SBarry Smith its--; 79417ab2063SBarry Smith } 79517ab2063SBarry Smith while (its--) { 79617ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 79717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 798416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 799416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 800416022c9SBarry Smith idx = a->j + a->i[i] + shift; 801416022c9SBarry Smith v = a->a + a->i[i] + shift; 80217ab2063SBarry Smith sum = b[i]; 80317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 80417ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 80517ab2063SBarry Smith } 80617ab2063SBarry Smith } 80717ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 80817ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 809416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 810416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 811416022c9SBarry Smith idx = a->j + a->i[i] + shift; 812416022c9SBarry Smith v = a->a + a->i[i] + shift; 81317ab2063SBarry Smith sum = b[i]; 81417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 81517ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 81617ab2063SBarry Smith } 81717ab2063SBarry Smith } 81817ab2063SBarry Smith } 81917ab2063SBarry Smith return 0; 82017ab2063SBarry Smith } 82117ab2063SBarry Smith 822d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 82317ab2063SBarry Smith { 824416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 825bcd2baecSBarry Smith if (nz) *nz = a->nz; 826bcd2baecSBarry Smith if (nzalloc) *nzalloc = a->maxnz; 827bcd2baecSBarry Smith if (mem) *mem = (int)A->mem; 82817ab2063SBarry Smith return 0; 82917ab2063SBarry Smith } 83017ab2063SBarry Smith 83117ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 83217ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 83317ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 83417ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 83517ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 83617ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 83717ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 83817ab2063SBarry Smith 83917ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 84017ab2063SBarry Smith { 841416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 842416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 84317ab2063SBarry Smith 84477c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 84517ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 84617ab2063SBarry Smith if (diag) { 84717ab2063SBarry Smith for ( i=0; i<N; i++ ) { 848416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 849416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 850416022c9SBarry Smith a->ilen[rows[i]] = 1; 851416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 852416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 85317ab2063SBarry Smith } 85417ab2063SBarry Smith else { 85517ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 85617ab2063SBarry Smith CHKERRQ(ierr); 85717ab2063SBarry Smith } 85817ab2063SBarry Smith } 85917ab2063SBarry Smith } 86017ab2063SBarry Smith else { 86117ab2063SBarry Smith for ( i=0; i<N; i++ ) { 862416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 863416022c9SBarry Smith a->ilen[rows[i]] = 0; 86417ab2063SBarry Smith } 86517ab2063SBarry Smith } 86617ab2063SBarry Smith ISRestoreIndices(is,&rows); 86717ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 86817ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 86917ab2063SBarry Smith return 0; 87017ab2063SBarry Smith } 87117ab2063SBarry Smith 872416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 87317ab2063SBarry Smith { 874416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 875416022c9SBarry Smith *m = a->m; *n = a->n; 87617ab2063SBarry Smith return 0; 87717ab2063SBarry Smith } 87817ab2063SBarry Smith 879416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 88017ab2063SBarry Smith { 881416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 882416022c9SBarry Smith *m = 0; *n = a->m; 88317ab2063SBarry Smith return 0; 88417ab2063SBarry Smith } 8854e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 88617ab2063SBarry Smith { 887416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 888c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 88917ab2063SBarry Smith 890416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 89117ab2063SBarry Smith 892416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 893416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 89417ab2063SBarry Smith if (idx) { 895416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8964e093b46SBarry Smith if (*nz && shift) { 8970452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 89817ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 8994e093b46SBarry Smith } else if (*nz) { 9004e093b46SBarry Smith *idx = itmp; 90117ab2063SBarry Smith } 90217ab2063SBarry Smith else *idx = 0; 90317ab2063SBarry Smith } 90417ab2063SBarry Smith return 0; 90517ab2063SBarry Smith } 90617ab2063SBarry Smith 9074e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 90817ab2063SBarry Smith { 9094e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 9104e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 91117ab2063SBarry Smith return 0; 91217ab2063SBarry Smith } 91317ab2063SBarry Smith 914cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 91517ab2063SBarry Smith { 916416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 917416022c9SBarry Smith Scalar *v = a->a; 91817ab2063SBarry Smith double sum = 0.0; 919416022c9SBarry Smith int i, j,shift = a->indexshift; 92017ab2063SBarry Smith 92117ab2063SBarry Smith if (type == NORM_FROBENIUS) { 922416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 92317ab2063SBarry Smith #if defined(PETSC_COMPLEX) 92417ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 92517ab2063SBarry Smith #else 92617ab2063SBarry Smith sum += (*v)*(*v); v++; 92717ab2063SBarry Smith #endif 92817ab2063SBarry Smith } 92917ab2063SBarry Smith *norm = sqrt(sum); 93017ab2063SBarry Smith } 93117ab2063SBarry Smith else if (type == NORM_1) { 93217ab2063SBarry Smith double *tmp; 933416022c9SBarry Smith int *jj = a->j; 9340452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 935cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 93617ab2063SBarry Smith *norm = 0.0; 937416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 938a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 93917ab2063SBarry Smith } 940416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 94117ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 94217ab2063SBarry Smith } 9430452661fSBarry Smith PetscFree(tmp); 94417ab2063SBarry Smith } 94517ab2063SBarry Smith else if (type == NORM_INFINITY) { 94617ab2063SBarry Smith *norm = 0.0; 947416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 948416022c9SBarry Smith v = a->a + a->i[j] + shift; 94917ab2063SBarry Smith sum = 0.0; 950416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 951cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 95217ab2063SBarry Smith } 95317ab2063SBarry Smith if (sum > *norm) *norm = sum; 95417ab2063SBarry Smith } 95517ab2063SBarry Smith } 95617ab2063SBarry Smith else { 95748d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 95817ab2063SBarry Smith } 95917ab2063SBarry Smith return 0; 96017ab2063SBarry Smith } 96117ab2063SBarry Smith 962416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 96317ab2063SBarry Smith { 964416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 965416022c9SBarry Smith Mat C; 966416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 967416022c9SBarry Smith int shift = a->indexshift; 968d5d45c9bSBarry Smith Scalar *array = a->a; 96917ab2063SBarry Smith 9703638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 971dfa27b74SSatish Balay SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place"); 9720452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 973cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 97417ab2063SBarry Smith if (shift) { 97517ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 97617ab2063SBarry Smith } 97717ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 978416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 9790452661fSBarry Smith PetscFree(col); 98017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 98117ab2063SBarry Smith len = ai[i+1]-ai[i]; 982416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 98317ab2063SBarry Smith array += len; aj += len; 98417ab2063SBarry Smith } 98517ab2063SBarry Smith if (shift) { 98617ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 98717ab2063SBarry Smith } 98817ab2063SBarry Smith 989416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 990416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 99117ab2063SBarry Smith 9923638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 993416022c9SBarry Smith *B = C; 99417ab2063SBarry Smith } else { 995416022c9SBarry Smith /* This isn't really an in-place transpose */ 9960452661fSBarry Smith PetscFree(a->a); 9970452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9980452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9990452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 10000452661fSBarry Smith if (a->imax) PetscFree(a->imax); 10010452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 10021073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 10030452661fSBarry Smith PetscFree(a); 1004416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 10050452661fSBarry Smith PetscHeaderDestroy(C); 100617ab2063SBarry Smith } 100717ab2063SBarry Smith return 0; 100817ab2063SBarry Smith } 100917ab2063SBarry Smith 1010f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 101117ab2063SBarry Smith { 1012416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 101317ab2063SBarry Smith Scalar *l,*r,x,*v; 1014d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 101517ab2063SBarry Smith 101617ab2063SBarry Smith if (ll) { 101717ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 1018f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 1019416022c9SBarry Smith v = a->a; 102017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 102117ab2063SBarry Smith x = l[i]; 1022416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 102317ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 102417ab2063SBarry Smith } 102544cd7ae7SLois Curfman McInnes PLogFlops(nz); 102617ab2063SBarry Smith } 102717ab2063SBarry Smith if (rr) { 102817ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 1029f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 1030416022c9SBarry Smith v = a->a; jj = a->j; 103117ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 103217ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 103317ab2063SBarry Smith } 103444cd7ae7SLois Curfman McInnes PLogFlops(nz); 103517ab2063SBarry Smith } 103617ab2063SBarry Smith return 0; 103717ab2063SBarry Smith } 103817ab2063SBarry Smith 1039cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 104017ab2063SBarry Smith { 1041db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 104202834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 104399141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1044a2744918SBarry Smith register int sum,lensi; 104599141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 104699141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 104799141d43SSatish Balay Scalar *a_new,*mat_a; 1048416022c9SBarry Smith Mat C; 104917ab2063SBarry Smith 105099141d43SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i); 105199141d43SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IS is not sorted"); 105299141d43SSatish Balay 105317ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 105417ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 105517ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 105617ab2063SBarry Smith 10577264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 105802834360SBarry Smith /* special case of contiguous rows */ 105957faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 106002834360SBarry Smith starts = lens + ncols; 106102834360SBarry Smith /* loop over new rows determining lens and starting points */ 106202834360SBarry Smith for (i=0; i<nrows; i++) { 1063a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1064a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 106502834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1066d8ced48eSBarry Smith if (aj[k]+shift >= first) { 106702834360SBarry Smith starts[i] = k; 106802834360SBarry Smith break; 106902834360SBarry Smith } 107002834360SBarry Smith } 1071a2744918SBarry Smith sum = 0; 107202834360SBarry Smith while (k < kend) { 1073d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1074a2744918SBarry Smith sum++; 107502834360SBarry Smith } 1076a2744918SBarry Smith lens[i] = sum; 107702834360SBarry Smith } 107802834360SBarry Smith /* create submatrix */ 1079cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 108008480c60SBarry Smith int n_cols,n_rows; 108108480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 108208480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1083d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 108408480c60SBarry Smith C = *B; 108508480c60SBarry Smith } 108608480c60SBarry Smith else { 108702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 108808480c60SBarry Smith } 1089db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1090db02288aSLois Curfman McInnes 109102834360SBarry Smith /* loop over rows inserting into submatrix */ 1092db02288aSLois Curfman McInnes a_new = c->a; 1093db02288aSLois Curfman McInnes j_new = c->j; 1094db02288aSLois Curfman McInnes i_new = c->i; 1095db02288aSLois Curfman McInnes i_new[0] = -shift; 109602834360SBarry Smith for (i=0; i<nrows; i++) { 1097a2744918SBarry Smith ii = starts[i]; 1098a2744918SBarry Smith lensi = lens[i]; 1099a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1100a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 110102834360SBarry Smith } 1102a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1103a2744918SBarry Smith a_new += lensi; 1104a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1105a2744918SBarry Smith c->ilen[i] = lensi; 110602834360SBarry Smith } 11070452661fSBarry Smith PetscFree(lens); 110802834360SBarry Smith } 110902834360SBarry Smith else { 111002834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 11110452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 111202834360SBarry Smith ssmap = smap + shift; 111399141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1114cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 111517ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 111602834360SBarry Smith /* determine lens of each row */ 111702834360SBarry Smith for (i=0; i<nrows; i++) { 1118d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 111902834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 112002834360SBarry Smith lens[i] = 0; 112102834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1122d8ced48eSBarry Smith if (ssmap[aj[k]]) { 112302834360SBarry Smith lens[i]++; 112402834360SBarry Smith } 112502834360SBarry Smith } 112602834360SBarry Smith } 112717ab2063SBarry Smith /* Create and fill new matrix */ 1128a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 112999141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 113099141d43SSatish Balay 113199141d43SSatish Balay if (c->m != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 113299141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 113399141d43SSatish Balay SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros"); 113499141d43SSatish Balay } 113599141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 113608480c60SBarry Smith C = *B; 113708480c60SBarry Smith } 113808480c60SBarry Smith else { 113902834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 114008480c60SBarry Smith } 114199141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 114217ab2063SBarry Smith for (i=0; i<nrows; i++) { 114399141d43SSatish Balay row = irow[i]; 114417ab2063SBarry Smith nznew = 0; 114599141d43SSatish Balay kstart = ai[row]+shift; 114699141d43SSatish Balay kend = kstart + a->ilen[row]; 114799141d43SSatish Balay mat_i = c->i[i]+shift; 114899141d43SSatish Balay mat_j = c->j + mat_i; 114999141d43SSatish Balay mat_a = c->a + mat_i; 115099141d43SSatish Balay mat_ilen = c->ilen + i; 115117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 115299141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 115399141d43SSatish Balay *mat_j++ = tcol - (!shift); 115499141d43SSatish Balay *mat_a++ = a->a[k]; 115599141d43SSatish Balay (*mat_ilen)++; 115699141d43SSatish Balay 115717ab2063SBarry Smith } 115817ab2063SBarry Smith } 115917ab2063SBarry Smith } 116002834360SBarry Smith /* Free work space */ 116102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 116299141d43SSatish Balay PetscFree(smap); PetscFree(lens); 116302834360SBarry Smith } 1164416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1165416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 116617ab2063SBarry Smith 116717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1168416022c9SBarry Smith *B = C; 116917ab2063SBarry Smith return 0; 117017ab2063SBarry Smith } 117117ab2063SBarry Smith 1172a871dcd8SBarry Smith /* 117363b91edcSBarry Smith note: This can only work for identity for row and col. It would 117463b91edcSBarry Smith be good to check this and otherwise generate an error. 1175a871dcd8SBarry Smith */ 117663b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1177a871dcd8SBarry Smith { 117863b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 117908480c60SBarry Smith int ierr; 118063b91edcSBarry Smith Mat outA; 118163b91edcSBarry Smith 1182a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1183a871dcd8SBarry Smith 118463b91edcSBarry Smith outA = inA; 118563b91edcSBarry Smith inA->factor = FACTOR_LU; 118663b91edcSBarry Smith a->row = row; 118763b91edcSBarry Smith a->col = col; 118863b91edcSBarry Smith 11890452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 119063b91edcSBarry Smith 119108480c60SBarry Smith if (!a->diag) { 119208480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 119363b91edcSBarry Smith } 119463b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1195a871dcd8SBarry Smith return 0; 1196a871dcd8SBarry Smith } 1197a871dcd8SBarry Smith 1198f0b747eeSBarry Smith #include "pinclude/plapack.h" 1199f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1200f0b747eeSBarry Smith { 1201f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1202f0b747eeSBarry Smith int one = 1; 1203f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1204f0b747eeSBarry Smith PLogFlops(a->nz); 1205f0b747eeSBarry Smith return 0; 1206f0b747eeSBarry Smith } 1207f0b747eeSBarry Smith 1208cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1209cddf8d76SBarry Smith Mat **B) 1210cddf8d76SBarry Smith { 1211cddf8d76SBarry Smith int ierr,i; 1212cddf8d76SBarry Smith 1213cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 12140452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1215cddf8d76SBarry Smith } 1216cddf8d76SBarry Smith 1217cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1218cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1219cddf8d76SBarry Smith } 1220cddf8d76SBarry Smith return 0; 1221cddf8d76SBarry Smith } 1222cddf8d76SBarry Smith 1223e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 12244dcbc457SBarry Smith { 1225e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 122606763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 12278a047759SSatish Balay int start, end, *ai, *aj; 122806763907SSatish Balay char *table; 12298a047759SSatish Balay shift = a->indexshift; 1230e4d965acSSatish Balay m = a->m; 1231e4d965acSSatish Balay ai = a->i; 12328a047759SSatish Balay aj = a->j+shift; 12338a047759SSatish Balay 1234e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 123506763907SSatish Balay 123606763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 123706763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 123806763907SSatish Balay 1239e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1240e4d965acSSatish Balay /* Initialise the two local arrays */ 1241e4d965acSSatish Balay isz = 0; 124206763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1243e4d965acSSatish Balay 1244e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 12454dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 124677c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1247e4d965acSSatish Balay 1248e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1249e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 125006763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 12514dcbc457SBarry Smith } 125206763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 125306763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1254e4d965acSSatish Balay 125504a348a9SBarry Smith k = 0; 125604a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 125704a348a9SBarry Smith n = isz; 125806763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1259e4d965acSSatish Balay row = nidx[k]; 1260e4d965acSSatish Balay start = ai[row]; 1261e4d965acSSatish Balay end = ai[row+1]; 126204a348a9SBarry Smith for ( l = start; l<end ; l++){ 12638a047759SSatish Balay val = aj[l] + shift; 126406763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1265e4d965acSSatish Balay } 1266e4d965acSSatish Balay } 1267e4d965acSSatish Balay } 1268e4d965acSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1269e4d965acSSatish Balay } 127004a348a9SBarry Smith PetscFree(table); 127106763907SSatish Balay PetscFree(nidx); 1272e4d965acSSatish Balay return 0; 12734dcbc457SBarry Smith } 127417ab2063SBarry Smith 1275682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1276682d7d0cSBarry Smith { 1277682d7d0cSBarry Smith static int called = 0; 1278682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1279682d7d0cSBarry Smith 1280682d7d0cSBarry Smith if (called) return 0; else called = 1; 128177c4ece6SBarry Smith PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 128277c4ece6SBarry Smith PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>\n"); 128377c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 128477c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_no_inode - Do not use inodes\n"); 128577c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1286682d7d0cSBarry Smith #if defined(HAVE_ESSL) 128777c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1288682d7d0cSBarry Smith #endif 1289682d7d0cSBarry Smith return 0; 1290682d7d0cSBarry Smith } 1291205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1292682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 129317ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 129417ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1295416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1296416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 129717ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 129817ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 129917ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 130017ab2063SBarry Smith MatRelax_SeqAIJ, 130117ab2063SBarry Smith MatTranspose_SeqAIJ, 13027264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1303f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 130417ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 130517ab2063SBarry Smith MatCompress_SeqAIJ, 130617ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 130717ab2063SBarry Smith MatGetReordering_SeqAIJ, 130817ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 130917ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 131017ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 131117ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1312416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 13133d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1314cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 13157eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1316682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1317f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1318f0b747eeSBarry Smith MatScale_SeqAIJ}; 131917ab2063SBarry Smith 132017ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 132117ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 132217ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 132317ab2063SBarry Smith 132417ab2063SBarry Smith /*@C 1325682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 13260d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 13276e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 13282bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 13292bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 133017ab2063SBarry Smith 133117ab2063SBarry Smith Input Parameters: 133217ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 133317ab2063SBarry Smith . m - number of rows 133417ab2063SBarry Smith . n - number of columns 133517ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 13362bd5e0b2SLois Curfman McInnes . nzz - array containing the number of nonzeros in the various rows 13372bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 133817ab2063SBarry Smith 133917ab2063SBarry Smith Output Parameter: 1340416022c9SBarry Smith . A - the matrix 134117ab2063SBarry Smith 134217ab2063SBarry Smith Notes: 134317ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 134417ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 13450002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 134644cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 134717ab2063SBarry Smith 134817ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1349a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1350*3d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 1351*3d323bbdSBarry Smith will get TERRIBLE performance, see the users' manual chapter on 13520d15e28bSLois Curfman McInnes matrices and the file $(PETSC_DIR)/Performance. 135317ab2063SBarry Smith 1354682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1355682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1356682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 13576c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 13586c7ebb05SLois Curfman McInnes 13596c7ebb05SLois Curfman McInnes Options Database Keys: 13606c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 13616e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 13626e62573dSLois Curfman McInnes $ (max limit=5) 13636e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 13646e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 13656e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 136617ab2063SBarry Smith 136717ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 136817ab2063SBarry Smith @*/ 1369416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 137017ab2063SBarry Smith { 1371416022c9SBarry Smith Mat B; 1372416022c9SBarry Smith Mat_SeqAIJ *b; 137369957df2SSatish Balay int i, len, ierr, flg; 1374d5d45c9bSBarry Smith 1375416022c9SBarry Smith *A = 0; 13760452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1377416022c9SBarry Smith PLogObjectCreate(B); 13780452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 137944cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1380416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1381416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1382416022c9SBarry Smith B->view = MatView_SeqAIJ; 1383416022c9SBarry Smith B->factor = 0; 1384416022c9SBarry Smith B->lupivotthreshold = 1.0; 13857a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 138669957df2SSatish Balay &flg); CHKERRQ(ierr); 13877a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 13887a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 13897a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1390416022c9SBarry Smith b->row = 0; 1391416022c9SBarry Smith b->col = 0; 1392416022c9SBarry Smith b->indexshift = 0; 1393b810aeb4SBarry Smith b->reallocs = 0; 139469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 139569957df2SSatish Balay if (flg) b->indexshift = -1; 139617ab2063SBarry Smith 139744cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 139844cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 13990452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1400b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 14017b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 14027b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1403416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 140417ab2063SBarry Smith nz = nz*m; 140517ab2063SBarry Smith } 140617ab2063SBarry Smith else { 140717ab2063SBarry Smith nz = 0; 1408416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 140917ab2063SBarry Smith } 141017ab2063SBarry Smith 141117ab2063SBarry Smith /* allocate the matrix space */ 1412416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 14130452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1414416022c9SBarry Smith b->j = (int *) (b->a + nz); 1415cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1416416022c9SBarry Smith b->i = b->j + nz; 1417416022c9SBarry Smith b->singlemalloc = 1; 141817ab2063SBarry Smith 1419416022c9SBarry Smith b->i[0] = -b->indexshift; 142017ab2063SBarry Smith for (i=1; i<m+1; i++) { 1421416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 142217ab2063SBarry Smith } 142317ab2063SBarry Smith 1424416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 14250452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1426416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1427416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 142817ab2063SBarry Smith 1429416022c9SBarry Smith b->nz = 0; 1430416022c9SBarry Smith b->maxnz = nz; 1431416022c9SBarry Smith b->sorted = 0; 1432416022c9SBarry Smith b->roworiented = 1; 1433416022c9SBarry Smith b->nonew = 0; 1434416022c9SBarry Smith b->diag = 0; 1435416022c9SBarry Smith b->solve_work = 0; 143671bd300dSLois Curfman McInnes b->spptr = 0; 1437754ec7b1SSatish Balay b->inode.node_count = 0; 1438754ec7b1SSatish Balay b->inode.size = 0; 14396c7ebb05SLois Curfman McInnes b->inode.limit = 5; 14406c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 144117ab2063SBarry Smith 1442416022c9SBarry Smith *A = B; 14434b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 14444b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 144569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 144669957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 14474b14c69eSBarry Smith #endif 144869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 144969957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 145069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 145169957df2SSatish Balay if (flg) { 1452416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1453416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 145417ab2063SBarry Smith } 145569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 145669957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 145717ab2063SBarry Smith return 0; 145817ab2063SBarry Smith } 145917ab2063SBarry Smith 14603d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 146117ab2063SBarry Smith { 1462416022c9SBarry Smith Mat C; 1463416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 146408480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 146517ab2063SBarry Smith 14664043dd9cSLois Curfman McInnes *B = 0; 14670452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1468416022c9SBarry Smith PLogObjectCreate(C); 14690452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 147041c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1471416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1472416022c9SBarry Smith C->view = MatView_SeqAIJ; 1473416022c9SBarry Smith C->factor = A->factor; 1474416022c9SBarry Smith c->row = 0; 1475416022c9SBarry Smith c->col = 0; 1476416022c9SBarry Smith c->indexshift = shift; 1477c456f294SBarry Smith C->assembled = PETSC_TRUE; 147817ab2063SBarry Smith 147944cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 148044cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 148144cd7ae7SLois Curfman McInnes C->M = a->m; 148244cd7ae7SLois Curfman McInnes C->N = a->n; 148317ab2063SBarry Smith 14840452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 14850452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 148617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1487416022c9SBarry Smith c->imax[i] = a->imax[i]; 1488416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 148917ab2063SBarry Smith } 149017ab2063SBarry Smith 149117ab2063SBarry Smith /* allocate the matrix space */ 1492416022c9SBarry Smith c->singlemalloc = 1; 1493416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 14940452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1495416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1496416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1497416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 149817ab2063SBarry Smith if (m > 0) { 1499416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 150008480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1501416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 150217ab2063SBarry Smith } 150308480c60SBarry Smith } 150417ab2063SBarry Smith 1505416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1506416022c9SBarry Smith c->sorted = a->sorted; 1507416022c9SBarry Smith c->roworiented = a->roworiented; 1508416022c9SBarry Smith c->nonew = a->nonew; 15097a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 151017ab2063SBarry Smith 1511416022c9SBarry Smith if (a->diag) { 15120452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1513416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 151417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1515416022c9SBarry Smith c->diag[i] = a->diag[i]; 151617ab2063SBarry Smith } 151717ab2063SBarry Smith } 1518416022c9SBarry Smith else c->diag = 0; 15196c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 15206c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1521754ec7b1SSatish Balay if (a->inode.size){ 1522754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1523754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1524754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1525754ec7b1SSatish Balay } else { 1526754ec7b1SSatish Balay c->inode.size = 0; 1527754ec7b1SSatish Balay c->inode.node_count = 0; 1528754ec7b1SSatish Balay } 1529416022c9SBarry Smith c->nz = a->nz; 1530416022c9SBarry Smith c->maxnz = a->maxnz; 1531416022c9SBarry Smith c->solve_work = 0; 153276dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1533754ec7b1SSatish Balay 1534416022c9SBarry Smith *B = C; 153517ab2063SBarry Smith return 0; 153617ab2063SBarry Smith } 153717ab2063SBarry Smith 153819bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 153917ab2063SBarry Smith { 1540416022c9SBarry Smith Mat_SeqAIJ *a; 1541416022c9SBarry Smith Mat B; 154217699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1543bcd2baecSBarry Smith MPI_Comm comm; 154417ab2063SBarry Smith 154519bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 154617699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 154717699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 154890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 154977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 155048d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 155117ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 155217ab2063SBarry Smith 155317ab2063SBarry Smith /* read in row lengths */ 15540452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 155577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 155617ab2063SBarry Smith 155717ab2063SBarry Smith /* create our matrix */ 1558416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1559416022c9SBarry Smith B = *A; 1560416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1561416022c9SBarry Smith shift = a->indexshift; 156217ab2063SBarry Smith 156317ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 156477c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 156517ab2063SBarry Smith if (shift) { 156617ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1567416022c9SBarry Smith a->j[i] += 1; 156817ab2063SBarry Smith } 156917ab2063SBarry Smith } 157017ab2063SBarry Smith 157117ab2063SBarry Smith /* read in nonzero values */ 157277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 157317ab2063SBarry Smith 157417ab2063SBarry Smith /* set matrix "i" values */ 1575416022c9SBarry Smith a->i[0] = -shift; 157617ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1577416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1578416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 157917ab2063SBarry Smith } 15800452661fSBarry Smith PetscFree(rowlengths); 158117ab2063SBarry Smith 1582416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1583416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 158417ab2063SBarry Smith return 0; 158517ab2063SBarry Smith } 158617ab2063SBarry Smith 1587686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 15887264ac53SSatish Balay { 15897264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 15907264ac53SSatish Balay 1591bcd2baecSBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type"); 15927264ac53SSatish Balay 15937264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 15947264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1595bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 159677c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1597bcd2baecSBarry Smith } 15987264ac53SSatish Balay 15997264ac53SSatish Balay /* if the a->i are the same */ 1600bcd2baecSBarry Smith if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 160177c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 16027264ac53SSatish Balay } 16037264ac53SSatish Balay 16047264ac53SSatish Balay /* if a->j are the same */ 1605bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 160677c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1607bcd2baecSBarry Smith } 1608bcd2baecSBarry Smith 1609bcd2baecSBarry Smith /* if a->a are the same */ 161019bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 161177c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 16127264ac53SSatish Balay } 161377c4ece6SBarry Smith *flg = PETSC_TRUE; 16147264ac53SSatish Balay return 0; 16157264ac53SSatish Balay 16167264ac53SSatish Balay } 1617