16d84be18SBarry Smith 217ab2063SBarry Smith #ifndef lint 3*bcd2baecSBarry Smith static char vcid[] = "$Id: aij.c,v 1.153 1996/03/04 05:15:52 bsmith 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" 1117ab2063SBarry Smith #include "vec/vecimpl.h" 1217ab2063SBarry Smith #include "inline/spops.h" 13e4d965acSSatish Balay #include "petsc.h" 1406763907SSatish Balay #include "inline/bitarray.h" 1517ab2063SBarry Smith 16*bcd2baecSBarry 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; 21*bcd2baecSBarry 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 42*bcd2baecSBarry Smith MatReorderingRegisterAll(); 43*bcd2baecSBarry Smith ishift = a->indexshift; 44*bcd2baecSBarry Smith oshift = -MatReorderingIndexShift[(int)type]; 45*bcd2baecSBarry Smith if (MatReorderingRequiresSymmetric[(int)type]) { 46*bcd2baecSBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja); 47*bcd2baecSBarry Smith CHKERRQ(ierr); 48416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 490452661fSBarry Smith PetscFree(ia); PetscFree(ja); 50*bcd2baecSBarry Smith } else { 51*bcd2baecSBarry Smith if (ishift == oshift) { 52*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 53*bcd2baecSBarry Smith } 54*bcd2baecSBarry Smith else if (ishift == -1) { 55*bcd2baecSBarry Smith /* temporarily subtract 1 from i and j indices */ 56*bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 57*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 58*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 59*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 60*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 61*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 62*bcd2baecSBarry Smith } else { 63*bcd2baecSBarry Smith /* temporarily add 1 to i and j indices */ 64*bcd2baecSBarry Smith int nz = a->i[a->n] - 1; 65*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]++; 66*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]++; 67*bcd2baecSBarry Smith ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm); CHKERRQ(ierr); 68*bcd2baecSBarry Smith for ( i=0; i<nz; i++ ) a->j[i]--; 69*bcd2baecSBarry Smith for ( i=0; i<a->n+1; i++ ) a->i[i]--; 70*bcd2baecSBarry Smith } 71*bcd2baecSBarry 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" 206416022c9SBarry Smith #include "sysio.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 213416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(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 } 224416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,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 } 231416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,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 */ 237416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,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; 24495e01e2fSLois Curfman McInnes int ierr, i,j, m = a->m, shift = a->indexshift, format, flg; 24517ab2063SBarry Smith FILE *fd; 24617ab2063SBarry Smith char *outputname; 24717ab2063SBarry Smith 248227d817aSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 249416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 250416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 25117ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 25295e01e2fSLois Curfman McInnes return 0; 25395e01e2fSLois Curfman McInnes } 25495e01e2fSLois Curfman McInnes else if (format == FILE_FORMAT_INFO_DETAILED) { 25595e01e2fSLois Curfman McInnes ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 25695e01e2fSLois Curfman McInnes if (flg) fprintf(fd," not using I-node routines\n"); 25795e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 25895e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 25917ab2063SBarry Smith } 26017ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 26117ab2063SBarry Smith int nz, nzalloc, mem; 262416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 263416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 26417ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 26517ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 26617ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 26717ab2063SBarry Smith 26817ab2063SBarry Smith for (i=0; i<m; i++) { 269416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 27017ab2063SBarry Smith #if defined(PETSC_COMPLEX) 2717a743949SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]), 272416022c9SBarry Smith imag(a->a[j])); 27317ab2063SBarry Smith #else 2747a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 27517ab2063SBarry Smith #endif 27617ab2063SBarry Smith } 27717ab2063SBarry Smith } 27817ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 27917ab2063SBarry Smith } 28017ab2063SBarry Smith else { 28117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 28217ab2063SBarry Smith fprintf(fd,"row %d:",i); 283416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 28417ab2063SBarry Smith #if defined(PETSC_COMPLEX) 285416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 286416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 28717ab2063SBarry Smith } 28817ab2063SBarry Smith else { 289416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 29017ab2063SBarry Smith } 29117ab2063SBarry Smith #else 292416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 29317ab2063SBarry Smith #endif 29417ab2063SBarry Smith } 29517ab2063SBarry Smith fprintf(fd,"\n"); 29617ab2063SBarry Smith } 29717ab2063SBarry Smith } 29817ab2063SBarry Smith fflush(fd); 299416022c9SBarry Smith return 0; 300416022c9SBarry Smith } 301416022c9SBarry Smith 302416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 303416022c9SBarry Smith { 304416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 305cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 306cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 307*bcd2baecSBarry Smith Draw draw; 308cddf8d76SBarry Smith DrawButton button; 309cddf8d76SBarry Smith 310*bcd2baecSBarry Smith ViewerDrawGetDraw(viewer,&draw); 311416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 312416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 313416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 314416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 315cddf8d76SBarry Smith color = DRAW_BLUE; 316416022c9SBarry Smith for ( i=0; i<m; i++ ) { 317cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 318416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 319cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 320cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 321cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 322cddf8d76SBarry Smith #else 323cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 324cddf8d76SBarry Smith #endif 325cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 326cddf8d76SBarry Smith } 327cddf8d76SBarry Smith } 328cddf8d76SBarry Smith color = DRAW_CYAN; 329cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 330cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 331cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 332cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 333cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 334cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 335cddf8d76SBarry Smith } 336cddf8d76SBarry Smith } 337cddf8d76SBarry Smith color = DRAW_RED; 338cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 339cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 340cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 341cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 342cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 343cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 344cddf8d76SBarry Smith #else 345cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 346cddf8d76SBarry Smith #endif 347cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 348416022c9SBarry Smith } 349416022c9SBarry Smith } 350416022c9SBarry Smith DrawFlush(draw); 351cddf8d76SBarry Smith DrawGetPause(draw,&pause); 352cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 353cddf8d76SBarry Smith 354cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 355cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 356cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 357cddf8d76SBarry Smith DrawClear(draw); 358cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 359cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 360cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 361cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 362cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 363cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 364cddf8d76SBarry Smith w *= scale; h *= scale; 365cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 366cddf8d76SBarry Smith color = DRAW_BLUE; 367cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 368cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 369cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 370cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 371cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 372cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 373cddf8d76SBarry Smith #else 374cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 375cddf8d76SBarry Smith #endif 376cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 377cddf8d76SBarry Smith } 378cddf8d76SBarry Smith } 379cddf8d76SBarry Smith color = DRAW_CYAN; 380cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 381cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 382cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 383cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 384cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 385cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 386cddf8d76SBarry Smith } 387cddf8d76SBarry Smith } 388cddf8d76SBarry Smith color = DRAW_RED; 389cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 390cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 391cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 392cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 393cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 394cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 395cddf8d76SBarry Smith #else 396cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 397cddf8d76SBarry Smith #endif 398cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 399cddf8d76SBarry Smith } 400cddf8d76SBarry Smith } 401cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 402cddf8d76SBarry Smith } 403416022c9SBarry Smith return 0; 404416022c9SBarry Smith } 405416022c9SBarry Smith 406*bcd2baecSBarry Smith 407416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 408416022c9SBarry Smith { 409416022c9SBarry Smith Mat A = (Mat) obj; 410416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 411*bcd2baecSBarry Smith ViewerType vtype; 412*bcd2baecSBarry Smith int ierr; 413416022c9SBarry Smith 414416022c9SBarry Smith if (!viewer) { 415*bcd2baecSBarry Smith viewer = STDOUT_VIEWER_SELF; 416416022c9SBarry Smith } 417*bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 418*bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 419416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 420416022c9SBarry Smith } 421*bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 422416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 423416022c9SBarry Smith } 424*bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 425416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 426416022c9SBarry Smith } 427*bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 428*bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 42917ab2063SBarry Smith } 43017ab2063SBarry Smith return 0; 43117ab2063SBarry Smith } 432c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 433416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 43417ab2063SBarry Smith { 435416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 43641c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 437416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 438416022c9SBarry Smith Scalar *aa = a->a, *ap; 43917ab2063SBarry Smith 44017ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 44117ab2063SBarry Smith 44217ab2063SBarry Smith for ( i=1; i<m; i++ ) { 443416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 44417ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 44517ab2063SBarry Smith if (fshift) { 446416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 44717ab2063SBarry Smith N = ailen[i]; 44817ab2063SBarry Smith for ( j=0; j<N; j++ ) { 44917ab2063SBarry Smith ip[j-fshift] = ip[j]; 45017ab2063SBarry Smith ap[j-fshift] = ap[j]; 45117ab2063SBarry Smith } 45217ab2063SBarry Smith } 45317ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 45417ab2063SBarry Smith } 45517ab2063SBarry Smith if (m) { 45617ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 45717ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 45817ab2063SBarry Smith } 45917ab2063SBarry Smith /* reset ilen and imax for each row */ 46017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 46117ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 46217ab2063SBarry Smith } 463416022c9SBarry Smith a->nz = ai[m] + shift; 46417ab2063SBarry Smith 46517ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 466416022c9SBarry Smith if (fshift && a->diag) { 4670452661fSBarry Smith PetscFree(a->diag); 468416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 469416022c9SBarry Smith a->diag = 0; 47017ab2063SBarry Smith } 471b810aeb4SBarry Smith PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n", 472b810aeb4SBarry Smith fshift,a->nz,m); 473b810aeb4SBarry Smith PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n", 474b810aeb4SBarry Smith a->reallocs); 47576dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 47641c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 47717ab2063SBarry Smith return 0; 47817ab2063SBarry Smith } 47917ab2063SBarry Smith 480416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 48117ab2063SBarry Smith { 482416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 483cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 48417ab2063SBarry Smith return 0; 48517ab2063SBarry Smith } 486416022c9SBarry Smith 48717ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 48817ab2063SBarry Smith { 489416022c9SBarry Smith Mat A = (Mat) obj; 490416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 491d5d45c9bSBarry Smith 49217ab2063SBarry Smith #if defined(PETSC_LOG) 493416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 49417ab2063SBarry Smith #endif 4950452661fSBarry Smith PetscFree(a->a); 4960452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4970452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4980452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4990452661fSBarry Smith if (a->imax) PetscFree(a->imax); 5000452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 50176dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 5020452661fSBarry Smith PetscFree(a); 503416022c9SBarry Smith PLogObjectDestroy(A); 5040452661fSBarry Smith PetscHeaderDestroy(A); 50517ab2063SBarry Smith return 0; 50617ab2063SBarry Smith } 50717ab2063SBarry Smith 508416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 50917ab2063SBarry Smith { 51017ab2063SBarry Smith return 0; 51117ab2063SBarry Smith } 51217ab2063SBarry Smith 513416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 51417ab2063SBarry Smith { 515416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 516416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 5174b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 518416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 519416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 520416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 521e2f28af5SBarry Smith else if (op == ROWS_SORTED || 522e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 523e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 524e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 525df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 526df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 5274d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 5286c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_1) a->inode.limit = 1; 5296c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_2) a->inode.limit = 2; 5306c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_3) a->inode.limit = 3; 5316c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_4) a->inode.limit = 4; 5326c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_5) a->inode.limit = 5; 533e2f28af5SBarry Smith else 5344d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 53517ab2063SBarry Smith return 0; 53617ab2063SBarry Smith } 53717ab2063SBarry Smith 538416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 53917ab2063SBarry Smith { 540416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 541416022c9SBarry Smith int i,j, n,shift = a->indexshift; 54217ab2063SBarry Smith Scalar *x, zero = 0.0; 54317ab2063SBarry Smith 54417ab2063SBarry Smith VecSet(&zero,v); 54517ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 546416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 547416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 548416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 549416022c9SBarry Smith if (a->j[j]+shift == i) { 550416022c9SBarry Smith x[i] = a->a[j]; 55117ab2063SBarry Smith break; 55217ab2063SBarry Smith } 55317ab2063SBarry Smith } 55417ab2063SBarry Smith } 55517ab2063SBarry Smith return 0; 55617ab2063SBarry Smith } 55717ab2063SBarry Smith 55817ab2063SBarry Smith /* -------------------------------------------------------*/ 55917ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 56017ab2063SBarry Smith /* -------------------------------------------------------*/ 561416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 56217ab2063SBarry Smith { 563416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 56417ab2063SBarry Smith Scalar *x, *y, *v, alpha; 565416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 56617ab2063SBarry Smith 56717ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 568cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 56917ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 57017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 571416022c9SBarry Smith idx = a->j + a->i[i] + shift; 572416022c9SBarry Smith v = a->a + a->i[i] + shift; 573416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 57417ab2063SBarry Smith alpha = x[i]; 57517ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 57617ab2063SBarry Smith } 577416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 57817ab2063SBarry Smith return 0; 57917ab2063SBarry Smith } 580d5d45c9bSBarry Smith 581416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,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); 58817ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 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 } 59717ab2063SBarry Smith return 0; 59817ab2063SBarry Smith } 59917ab2063SBarry Smith 600416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 60117ab2063SBarry Smith { 602416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 60317ab2063SBarry Smith Scalar *x, *y, *v, sum; 604416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 60517ab2063SBarry Smith 60617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 60717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 608416022c9SBarry Smith idx = a->j; 609416022c9SBarry Smith v = a->a; 610416022c9SBarry Smith ii = a->i; 61117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 612416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 61317ab2063SBarry Smith sum = 0.0; 61417ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 61517ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 616416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 61717ab2063SBarry Smith y[i] = sum; 61817ab2063SBarry Smith } 619416022c9SBarry Smith PLogFlops(2*a->nz - m); 62017ab2063SBarry Smith return 0; 62117ab2063SBarry Smith } 62217ab2063SBarry Smith 623416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 62417ab2063SBarry Smith { 625416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 62617ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 627cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 62817ab2063SBarry Smith 62917ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 63017ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 631cddf8d76SBarry Smith idx = a->j; 632cddf8d76SBarry Smith v = a->a; 633cddf8d76SBarry Smith ii = a->i; 63417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 635cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 63617ab2063SBarry Smith sum = y[i]; 637cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 638cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 63917ab2063SBarry Smith z[i] = sum; 64017ab2063SBarry Smith } 641416022c9SBarry Smith PLogFlops(2*a->nz); 64217ab2063SBarry Smith return 0; 64317ab2063SBarry Smith } 64417ab2063SBarry Smith 64517ab2063SBarry Smith /* 64617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 64717ab2063SBarry Smith */ 64817ab2063SBarry Smith 649416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 65017ab2063SBarry Smith { 651416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 652416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 65317ab2063SBarry Smith 6540452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 655416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 656416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 657416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 658416022c9SBarry Smith if (a->j[j]+shift == i) { 65917ab2063SBarry Smith diag[i] = j - shift; 66017ab2063SBarry Smith break; 66117ab2063SBarry Smith } 66217ab2063SBarry Smith } 66317ab2063SBarry Smith } 664416022c9SBarry Smith a->diag = diag; 66517ab2063SBarry Smith return 0; 66617ab2063SBarry Smith } 66717ab2063SBarry Smith 668416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 66917ab2063SBarry Smith double fshift,int its,Vec xx) 67017ab2063SBarry Smith { 671416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 672416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 673d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 67417ab2063SBarry Smith 67517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 676416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 677416022c9SBarry Smith diag = a->diag; 678416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 67917ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 68017ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 68117ab2063SBarry Smith bs = b + shift; 68217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 683416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 684416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 685416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 686416022c9SBarry Smith v = a->a + diag[i] + (!shift); 68717ab2063SBarry Smith sum = b[i]*d/omega; 68817ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 68917ab2063SBarry Smith x[i] = sum; 69017ab2063SBarry Smith } 69117ab2063SBarry Smith return 0; 69217ab2063SBarry Smith } 69317ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 69417ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 69517ab2063SBarry Smith } 696416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 69717ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 69817ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 69917ab2063SBarry Smith 70017ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 70117ab2063SBarry Smith 70217ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 70317ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 70417ab2063SBarry Smith is the relaxation factor. 70517ab2063SBarry Smith */ 7060452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 70717ab2063SBarry Smith scale = (2.0/omega) - 1.0; 70817ab2063SBarry Smith 70917ab2063SBarry Smith /* x = (E + U)^{-1} b */ 71017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 711416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 712416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 713416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 714416022c9SBarry Smith v = a->a + diag[i] + (!shift); 71517ab2063SBarry Smith sum = b[i]; 71617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 71717ab2063SBarry Smith x[i] = omega*(sum/d); 71817ab2063SBarry Smith } 71917ab2063SBarry Smith 72017ab2063SBarry Smith /* t = b - (2*E - D)x */ 721416022c9SBarry Smith v = a->a; 72217ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 72317ab2063SBarry Smith 72417ab2063SBarry Smith /* t = (E + L)^{-1}t */ 725416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 726416022c9SBarry Smith diag = a->diag; 72717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 728416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 729416022c9SBarry Smith n = diag[i] - a->i[i]; 730416022c9SBarry Smith idx = a->j + a->i[i] + shift; 731416022c9SBarry Smith v = a->a + a->i[i] + shift; 73217ab2063SBarry Smith sum = t[i]; 73317ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 73417ab2063SBarry Smith t[i] = omega*(sum/d); 73517ab2063SBarry Smith } 73617ab2063SBarry Smith 73717ab2063SBarry Smith /* x = x + t */ 73817ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7390452661fSBarry Smith PetscFree(t); 74017ab2063SBarry Smith return 0; 74117ab2063SBarry Smith } 74217ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 74317ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 74417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 745416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 746416022c9SBarry Smith n = diag[i] - a->i[i]; 747416022c9SBarry Smith idx = a->j + a->i[i] + shift; 748416022c9SBarry Smith v = a->a + a->i[i] + shift; 74917ab2063SBarry Smith sum = b[i]; 75017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 75117ab2063SBarry Smith x[i] = omega*(sum/d); 75217ab2063SBarry Smith } 75317ab2063SBarry Smith xb = x; 75417ab2063SBarry Smith } 75517ab2063SBarry Smith else xb = b; 75617ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 75717ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 75817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 759416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 76017ab2063SBarry Smith } 76117ab2063SBarry Smith } 76217ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 76317ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 764416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 765416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 766416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 767416022c9SBarry Smith v = a->a + diag[i] + (!shift); 76817ab2063SBarry Smith sum = xb[i]; 76917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 77017ab2063SBarry Smith x[i] = omega*(sum/d); 77117ab2063SBarry Smith } 77217ab2063SBarry Smith } 77317ab2063SBarry Smith its--; 77417ab2063SBarry Smith } 77517ab2063SBarry Smith while (its--) { 77617ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 77717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 778416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 779416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 780416022c9SBarry Smith idx = a->j + a->i[i] + shift; 781416022c9SBarry Smith v = a->a + a->i[i] + shift; 78217ab2063SBarry Smith sum = b[i]; 78317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 78417ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 78517ab2063SBarry Smith } 78617ab2063SBarry Smith } 78717ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 78817ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 789416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 790416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 791416022c9SBarry Smith idx = a->j + a->i[i] + shift; 792416022c9SBarry Smith v = a->a + a->i[i] + shift; 79317ab2063SBarry Smith sum = b[i]; 79417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 79517ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 79617ab2063SBarry Smith } 79717ab2063SBarry Smith } 79817ab2063SBarry Smith } 79917ab2063SBarry Smith return 0; 80017ab2063SBarry Smith } 80117ab2063SBarry Smith 802d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 80317ab2063SBarry Smith { 804416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 805*bcd2baecSBarry Smith if (nz) *nz = a->nz; 806*bcd2baecSBarry Smith if (nzalloc) *nzalloc = a->maxnz; 807*bcd2baecSBarry Smith if (mem) *mem = (int)A->mem; 80817ab2063SBarry Smith return 0; 80917ab2063SBarry Smith } 81017ab2063SBarry Smith 81117ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 81217ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 81317ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 81417ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 81517ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 81617ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 81717ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 81817ab2063SBarry Smith 81917ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 82017ab2063SBarry Smith { 821416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 822416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 82317ab2063SBarry Smith 82417ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 82517ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 82617ab2063SBarry Smith if (diag) { 82717ab2063SBarry Smith for ( i=0; i<N; i++ ) { 828416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 829416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 830416022c9SBarry Smith a->ilen[rows[i]] = 1; 831416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 832416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 83317ab2063SBarry Smith } 83417ab2063SBarry Smith else { 83517ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 83617ab2063SBarry Smith CHKERRQ(ierr); 83717ab2063SBarry Smith } 83817ab2063SBarry Smith } 83917ab2063SBarry Smith } 84017ab2063SBarry Smith else { 84117ab2063SBarry Smith for ( i=0; i<N; i++ ) { 842416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 843416022c9SBarry Smith a->ilen[rows[i]] = 0; 84417ab2063SBarry Smith } 84517ab2063SBarry Smith } 84617ab2063SBarry Smith ISRestoreIndices(is,&rows); 84717ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 84817ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 84917ab2063SBarry Smith return 0; 85017ab2063SBarry Smith } 85117ab2063SBarry Smith 852416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 85317ab2063SBarry Smith { 854416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 855416022c9SBarry Smith *m = a->m; *n = a->n; 85617ab2063SBarry Smith return 0; 85717ab2063SBarry Smith } 85817ab2063SBarry Smith 859416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 86017ab2063SBarry Smith { 861416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 862416022c9SBarry Smith *m = 0; *n = a->m; 86317ab2063SBarry Smith return 0; 86417ab2063SBarry Smith } 8654e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 86617ab2063SBarry Smith { 867416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 868c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 86917ab2063SBarry Smith 870416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 87117ab2063SBarry Smith 872416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 873416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 87417ab2063SBarry Smith if (idx) { 875416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8764e093b46SBarry Smith if (*nz && shift) { 8770452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 87817ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 8794e093b46SBarry Smith } else if (*nz) { 8804e093b46SBarry Smith *idx = itmp; 88117ab2063SBarry Smith } 88217ab2063SBarry Smith else *idx = 0; 88317ab2063SBarry Smith } 88417ab2063SBarry Smith return 0; 88517ab2063SBarry Smith } 88617ab2063SBarry Smith 8874e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 88817ab2063SBarry Smith { 8894e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 8904e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 89117ab2063SBarry Smith return 0; 89217ab2063SBarry Smith } 89317ab2063SBarry Smith 894cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 89517ab2063SBarry Smith { 896416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 897416022c9SBarry Smith Scalar *v = a->a; 89817ab2063SBarry Smith double sum = 0.0; 899416022c9SBarry Smith int i, j,shift = a->indexshift; 90017ab2063SBarry Smith 90117ab2063SBarry Smith if (type == NORM_FROBENIUS) { 902416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 90317ab2063SBarry Smith #if defined(PETSC_COMPLEX) 90417ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 90517ab2063SBarry Smith #else 90617ab2063SBarry Smith sum += (*v)*(*v); v++; 90717ab2063SBarry Smith #endif 90817ab2063SBarry Smith } 90917ab2063SBarry Smith *norm = sqrt(sum); 91017ab2063SBarry Smith } 91117ab2063SBarry Smith else if (type == NORM_1) { 91217ab2063SBarry Smith double *tmp; 913416022c9SBarry Smith int *jj = a->j; 9140452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 915cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 91617ab2063SBarry Smith *norm = 0.0; 917416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 918a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 91917ab2063SBarry Smith } 920416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 92117ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 92217ab2063SBarry Smith } 9230452661fSBarry Smith PetscFree(tmp); 92417ab2063SBarry Smith } 92517ab2063SBarry Smith else if (type == NORM_INFINITY) { 92617ab2063SBarry Smith *norm = 0.0; 927416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 928416022c9SBarry Smith v = a->a + a->i[j] + shift; 92917ab2063SBarry Smith sum = 0.0; 930416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 931cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 93217ab2063SBarry Smith } 93317ab2063SBarry Smith if (sum > *norm) *norm = sum; 93417ab2063SBarry Smith } 93517ab2063SBarry Smith } 93617ab2063SBarry Smith else { 93748d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 93817ab2063SBarry Smith } 93917ab2063SBarry Smith return 0; 94017ab2063SBarry Smith } 94117ab2063SBarry Smith 942416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 94317ab2063SBarry Smith { 944416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 945416022c9SBarry Smith Mat C; 946416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 947416022c9SBarry Smith int shift = a->indexshift; 948d5d45c9bSBarry Smith Scalar *array = a->a; 94917ab2063SBarry Smith 9503638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 9513638b69dSLois Curfman McInnes SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place"); 9520452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 953cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 95417ab2063SBarry Smith if (shift) { 95517ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 95617ab2063SBarry Smith } 95717ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 958416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 9590452661fSBarry Smith PetscFree(col); 96017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 96117ab2063SBarry Smith len = ai[i+1]-ai[i]; 962416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 96317ab2063SBarry Smith array += len; aj += len; 96417ab2063SBarry Smith } 96517ab2063SBarry Smith if (shift) { 96617ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 96717ab2063SBarry Smith } 96817ab2063SBarry Smith 969416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 970416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 97117ab2063SBarry Smith 9723638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 973416022c9SBarry Smith *B = C; 97417ab2063SBarry Smith } else { 975416022c9SBarry Smith /* This isn't really an in-place transpose */ 9760452661fSBarry Smith PetscFree(a->a); 9770452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9780452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9790452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9800452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9810452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9820452661fSBarry Smith PetscFree(a); 983416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9840452661fSBarry Smith PetscHeaderDestroy(C); 98517ab2063SBarry Smith } 98617ab2063SBarry Smith return 0; 98717ab2063SBarry Smith } 98817ab2063SBarry Smith 989f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 99017ab2063SBarry Smith { 991416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 99217ab2063SBarry Smith Scalar *l,*r,x,*v; 993d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 99417ab2063SBarry Smith 99517ab2063SBarry Smith if (ll) { 99617ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 997f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 998416022c9SBarry Smith v = a->a; 99917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 100017ab2063SBarry Smith x = l[i]; 1001416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 100217ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 100317ab2063SBarry Smith } 100417ab2063SBarry Smith } 100517ab2063SBarry Smith if (rr) { 100617ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 1007f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 1008416022c9SBarry Smith v = a->a; jj = a->j; 100917ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 101017ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 101117ab2063SBarry Smith } 101217ab2063SBarry Smith } 101317ab2063SBarry Smith return 0; 101417ab2063SBarry Smith } 101517ab2063SBarry Smith 1016cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 101717ab2063SBarry Smith { 1018db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 101902834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 1020a2744918SBarry Smith register int sum,lensi; 102102834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 1022a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 1023db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 1024416022c9SBarry Smith Mat C; 102517ab2063SBarry Smith 102617ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 102717ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 102817ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 102917ab2063SBarry Smith 10307264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 103102834360SBarry Smith /* special case of contiguous rows */ 10320452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 103302834360SBarry Smith starts = lens + ncols; 103402834360SBarry Smith /* loop over new rows determining lens and starting points */ 103502834360SBarry Smith for (i=0; i<nrows; i++) { 1036a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1037a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 103802834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1039d8ced48eSBarry Smith if (aj[k]+shift >= first) { 104002834360SBarry Smith starts[i] = k; 104102834360SBarry Smith break; 104202834360SBarry Smith } 104302834360SBarry Smith } 1044a2744918SBarry Smith sum = 0; 104502834360SBarry Smith while (k < kend) { 1046d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1047a2744918SBarry Smith sum++; 104802834360SBarry Smith } 1049a2744918SBarry Smith lens[i] = sum; 105002834360SBarry Smith } 105102834360SBarry Smith /* create submatrix */ 1052cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 105308480c60SBarry Smith int n_cols,n_rows; 105408480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 105508480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1056d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 105708480c60SBarry Smith C = *B; 105808480c60SBarry Smith } 105908480c60SBarry Smith else { 106002834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 106108480c60SBarry Smith } 1062db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1063db02288aSLois Curfman McInnes 106402834360SBarry Smith /* loop over rows inserting into submatrix */ 1065db02288aSLois Curfman McInnes a_new = c->a; 1066db02288aSLois Curfman McInnes j_new = c->j; 1067db02288aSLois Curfman McInnes i_new = c->i; 1068db02288aSLois Curfman McInnes i_new[0] = -shift; 106902834360SBarry Smith for (i=0; i<nrows; i++) { 1070a2744918SBarry Smith ii = starts[i]; 1071a2744918SBarry Smith lensi = lens[i]; 1072a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1073a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 107402834360SBarry Smith } 1075a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1076a2744918SBarry Smith a_new += lensi; 1077a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1078a2744918SBarry Smith c->ilen[i] = lensi; 107902834360SBarry Smith } 10800452661fSBarry Smith PetscFree(lens); 108102834360SBarry Smith } 108202834360SBarry Smith else { 108302834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 10847264ac53SSatish Balay ierr = SYIsort(ncols,icol); CHKERRQ(ierr); 10850452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 108602834360SBarry Smith ssmap = smap + shift; 10877eb43aa7SLois Curfman McInnes cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 108802834360SBarry Smith lens = cwork + ncols; 10890452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1090cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 109117ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 109202834360SBarry Smith /* determine lens of each row */ 109302834360SBarry Smith for (i=0; i<nrows; i++) { 1094d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 109502834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 109602834360SBarry Smith lens[i] = 0; 109702834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1098d8ced48eSBarry Smith if (ssmap[aj[k]]) { 109902834360SBarry Smith lens[i]++; 110002834360SBarry Smith } 110102834360SBarry Smith } 110202834360SBarry Smith } 110317ab2063SBarry Smith /* Create and fill new matrix */ 1104a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 110508480c60SBarry Smith int n_cols,n_rows; 110608480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 110708480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1108d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 110908480c60SBarry Smith C = *B; 111008480c60SBarry Smith } 111108480c60SBarry Smith else { 111202834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 111308480c60SBarry Smith } 111417ab2063SBarry Smith for (i=0; i<nrows; i++) { 111517ab2063SBarry Smith nznew = 0; 1116d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 1117416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 111817ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 111902834360SBarry Smith if (ssmap[a->j[k]]) { 112002834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1121416022c9SBarry Smith vwork[nznew++] = a->a[k]; 112217ab2063SBarry Smith } 112317ab2063SBarry Smith } 1124416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 112517ab2063SBarry Smith } 112602834360SBarry Smith /* Free work space */ 112702834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 11280452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 112902834360SBarry Smith } 1130416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1131416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 113217ab2063SBarry Smith 113317ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1134416022c9SBarry Smith *B = C; 113517ab2063SBarry Smith return 0; 113617ab2063SBarry Smith } 113717ab2063SBarry Smith 1138a871dcd8SBarry Smith /* 113963b91edcSBarry Smith note: This can only work for identity for row and col. It would 114063b91edcSBarry Smith be good to check this and otherwise generate an error. 1141a871dcd8SBarry Smith */ 114263b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1143a871dcd8SBarry Smith { 114463b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 114508480c60SBarry Smith int ierr; 114663b91edcSBarry Smith Mat outA; 114763b91edcSBarry Smith 1148a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1149a871dcd8SBarry Smith 115063b91edcSBarry Smith outA = inA; 115163b91edcSBarry Smith inA->factor = FACTOR_LU; 115263b91edcSBarry Smith a->row = row; 115363b91edcSBarry Smith a->col = col; 115463b91edcSBarry Smith 11550452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 115663b91edcSBarry Smith 115708480c60SBarry Smith if (!a->diag) { 115808480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 115963b91edcSBarry Smith } 116063b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1161a871dcd8SBarry Smith return 0; 1162a871dcd8SBarry Smith } 1163a871dcd8SBarry Smith 1164f0b747eeSBarry Smith #include "pinclude/plapack.h" 1165f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1166f0b747eeSBarry Smith { 1167f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1168f0b747eeSBarry Smith int one = 1; 1169f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1170f0b747eeSBarry Smith PLogFlops(a->nz); 1171f0b747eeSBarry Smith return 0; 1172f0b747eeSBarry Smith } 1173f0b747eeSBarry Smith 1174cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1175cddf8d76SBarry Smith Mat **B) 1176cddf8d76SBarry Smith { 1177cddf8d76SBarry Smith int ierr,i; 1178cddf8d76SBarry Smith 1179cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 11800452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1181cddf8d76SBarry Smith } 1182cddf8d76SBarry Smith 1183cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1184cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1185cddf8d76SBarry Smith } 1186cddf8d76SBarry Smith return 0; 1187cddf8d76SBarry Smith } 1188cddf8d76SBarry Smith 1189e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 11904dcbc457SBarry Smith { 1191e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 119206763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 11938a047759SSatish Balay int start, end, *ai, *aj; 119406763907SSatish Balay char *table; 11958a047759SSatish Balay shift = a->indexshift; 1196e4d965acSSatish Balay m = a->m; 1197e4d965acSSatish Balay ai = a->i; 11988a047759SSatish Balay aj = a->j+shift; 11998a047759SSatish Balay 1200e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 120106763907SSatish Balay 120206763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 120306763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 120406763907SSatish Balay 1205e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1206e4d965acSSatish Balay /* Initialise the two local arrays */ 1207e4d965acSSatish Balay isz = 0; 120806763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1209e4d965acSSatish Balay 1210e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 12114dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1212e4d965acSSatish Balay ierr = ISGetLocalSize(is[i],&n); CHKERRQ(ierr); 1213e4d965acSSatish Balay 1214e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1215e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 121606763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 12174dcbc457SBarry Smith } 121806763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 121906763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1220e4d965acSSatish Balay 122104a348a9SBarry Smith k = 0; 122204a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 122304a348a9SBarry Smith n = isz; 122406763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1225e4d965acSSatish Balay row = nidx[k]; 1226e4d965acSSatish Balay start = ai[row]; 1227e4d965acSSatish Balay end = ai[row+1]; 122804a348a9SBarry Smith for ( l = start; l<end ; l++){ 12298a047759SSatish Balay val = aj[l] + shift; 123006763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1231e4d965acSSatish Balay } 1232e4d965acSSatish Balay } 1233e4d965acSSatish Balay } 1234e4d965acSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1235e4d965acSSatish Balay } 123604a348a9SBarry Smith PetscFree(table); 123706763907SSatish Balay PetscFree(nidx); 1238e4d965acSSatish Balay return 0; 12394dcbc457SBarry Smith } 124017ab2063SBarry Smith 1241682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1242682d7d0cSBarry Smith { 1243682d7d0cSBarry Smith static int called = 0; 1244682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1245682d7d0cSBarry Smith 1246682d7d0cSBarry Smith if (called) return 0; else called = 1; 12472a7368beSLois Curfman McInnes MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 1248682d7d0cSBarry Smith MPIU_printf(comm," -mat_lu_pivotthreshold <threshold>\n"); 12492a7368beSLois Curfman McInnes MPIU_printf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 1250682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_no_inode - Do not use inodes\n"); 1251682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1252682d7d0cSBarry Smith #if defined(HAVE_ESSL) 1253682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1254682d7d0cSBarry Smith #endif 1255682d7d0cSBarry Smith return 0; 1256682d7d0cSBarry Smith } 12577264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg); 1258682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 125917ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 126017ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1261416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1262416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 126317ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 126417ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 126517ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 126617ab2063SBarry Smith MatRelax_SeqAIJ, 126717ab2063SBarry Smith MatTranspose_SeqAIJ, 12687264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1269f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 127017ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 127117ab2063SBarry Smith MatCompress_SeqAIJ, 127217ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 127317ab2063SBarry Smith MatGetReordering_SeqAIJ, 127417ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 127517ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 127617ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 127717ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1278416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 12793d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1280cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 12817eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1282682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1283f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1284f0b747eeSBarry Smith MatScale_SeqAIJ}; 128517ab2063SBarry Smith 128617ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 128717ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 128817ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 128917ab2063SBarry Smith 129017ab2063SBarry Smith /*@C 1291682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 12920d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 12936e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 12940d15e28bSLois Curfman McInnes (or nzz). By setting these parameters accurately, performance can be 12950d15e28bSLois Curfman McInnes increased by more than a factor of 50. 129617ab2063SBarry Smith 129717ab2063SBarry Smith Input Parameters: 129817ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 129917ab2063SBarry Smith . m - number of rows 130017ab2063SBarry Smith . n - number of columns 130117ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 130217ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 130317ab2063SBarry Smith 130417ab2063SBarry Smith Output Parameter: 1305416022c9SBarry Smith . A - the matrix 130617ab2063SBarry Smith 130717ab2063SBarry Smith Notes: 130817ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 130917ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 13100002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 13110002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 131217ab2063SBarry Smith 131317ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1314a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 13150d15e28bSLois Curfman McInnes allocation. For additional details, see the users manual chapter on 13160d15e28bSLois Curfman McInnes matrices and the file $(PETSC_DIR)/Performance. 131717ab2063SBarry Smith 1318682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1319682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1320682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 13216c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 13226c7ebb05SLois Curfman McInnes 13236c7ebb05SLois Curfman McInnes Options Database Keys: 13246c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 13256e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 13266e62573dSLois Curfman McInnes $ (max limit=5) 13276e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 13286e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 13296e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 133017ab2063SBarry Smith 133117ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 133217ab2063SBarry Smith @*/ 1333416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 133417ab2063SBarry Smith { 1335416022c9SBarry Smith Mat B; 1336416022c9SBarry Smith Mat_SeqAIJ *b; 133769957df2SSatish Balay int i,len,ierr, flg; 1338d5d45c9bSBarry Smith 1339416022c9SBarry Smith *A = 0; 13400452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1341416022c9SBarry Smith PLogObjectCreate(B); 13420452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1343416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1344416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1345416022c9SBarry Smith B->view = MatView_SeqAIJ; 1346416022c9SBarry Smith B->factor = 0; 1347416022c9SBarry Smith B->lupivotthreshold = 1.0; 13487a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 134969957df2SSatish Balay &flg); CHKERRQ(ierr); 13507a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 13517a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 13527a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1353416022c9SBarry Smith b->row = 0; 1354416022c9SBarry Smith b->col = 0; 1355416022c9SBarry Smith b->indexshift = 0; 1356b810aeb4SBarry Smith b->reallocs = 0; 135769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 135869957df2SSatish Balay if (flg) b->indexshift = -1; 135917ab2063SBarry Smith 1360416022c9SBarry Smith b->m = m; 1361416022c9SBarry Smith b->n = n; 13620452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1363b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 13647b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 13657b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1366416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 136717ab2063SBarry Smith nz = nz*m; 136817ab2063SBarry Smith } 136917ab2063SBarry Smith else { 137017ab2063SBarry Smith nz = 0; 1371416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 137217ab2063SBarry Smith } 137317ab2063SBarry Smith 137417ab2063SBarry Smith /* allocate the matrix space */ 1375416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 13760452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1377416022c9SBarry Smith b->j = (int *) (b->a + nz); 1378cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1379416022c9SBarry Smith b->i = b->j + nz; 1380416022c9SBarry Smith b->singlemalloc = 1; 138117ab2063SBarry Smith 1382416022c9SBarry Smith b->i[0] = -b->indexshift; 138317ab2063SBarry Smith for (i=1; i<m+1; i++) { 1384416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 138517ab2063SBarry Smith } 138617ab2063SBarry Smith 1387416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 13880452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1389416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1390416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 139117ab2063SBarry Smith 1392416022c9SBarry Smith b->nz = 0; 1393416022c9SBarry Smith b->maxnz = nz; 1394416022c9SBarry Smith b->sorted = 0; 1395416022c9SBarry Smith b->roworiented = 1; 1396416022c9SBarry Smith b->nonew = 0; 1397416022c9SBarry Smith b->diag = 0; 1398416022c9SBarry Smith b->solve_work = 0; 139971bd300dSLois Curfman McInnes b->spptr = 0; 1400754ec7b1SSatish Balay b->inode.node_count = 0; 1401754ec7b1SSatish Balay b->inode.size = 0; 14026c7ebb05SLois Curfman McInnes b->inode.limit = 5; 14036c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 140417ab2063SBarry Smith 1405416022c9SBarry Smith *A = B; 14064b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 14074b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 140869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 140969957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 14104b14c69eSBarry Smith #endif 141169957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 141269957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 141369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 141469957df2SSatish Balay if (flg) { 1415416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1416416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 141717ab2063SBarry Smith } 141869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 141969957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 142017ab2063SBarry Smith return 0; 142117ab2063SBarry Smith } 142217ab2063SBarry Smith 14233d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 142417ab2063SBarry Smith { 1425416022c9SBarry Smith Mat C; 1426416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 142708480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 142817ab2063SBarry Smith 14294043dd9cSLois Curfman McInnes *B = 0; 14300452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1431416022c9SBarry Smith PLogObjectCreate(C); 14320452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 143341c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1434416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1435416022c9SBarry Smith C->view = MatView_SeqAIJ; 1436416022c9SBarry Smith C->factor = A->factor; 1437416022c9SBarry Smith c->row = 0; 1438416022c9SBarry Smith c->col = 0; 1439416022c9SBarry Smith c->indexshift = shift; 1440c456f294SBarry Smith C->assembled = PETSC_TRUE; 144117ab2063SBarry Smith 1442416022c9SBarry Smith c->m = a->m; 1443416022c9SBarry Smith c->n = a->n; 144417ab2063SBarry Smith 14450452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 14460452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 144717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1448416022c9SBarry Smith c->imax[i] = a->imax[i]; 1449416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 145017ab2063SBarry Smith } 145117ab2063SBarry Smith 145217ab2063SBarry Smith /* allocate the matrix space */ 1453416022c9SBarry Smith c->singlemalloc = 1; 1454416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 14550452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1456416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1457416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1458416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 145917ab2063SBarry Smith if (m > 0) { 1460416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 146108480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1462416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 146317ab2063SBarry Smith } 146408480c60SBarry Smith } 146517ab2063SBarry Smith 1466416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1467416022c9SBarry Smith c->sorted = a->sorted; 1468416022c9SBarry Smith c->roworiented = a->roworiented; 1469416022c9SBarry Smith c->nonew = a->nonew; 14707a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 147117ab2063SBarry Smith 1472416022c9SBarry Smith if (a->diag) { 14730452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1474416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 147517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1476416022c9SBarry Smith c->diag[i] = a->diag[i]; 147717ab2063SBarry Smith } 147817ab2063SBarry Smith } 1479416022c9SBarry Smith else c->diag = 0; 14806c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 14816c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1482754ec7b1SSatish Balay if (a->inode.size){ 1483754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1484754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1485754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1486754ec7b1SSatish Balay } else { 1487754ec7b1SSatish Balay c->inode.size = 0; 1488754ec7b1SSatish Balay c->inode.node_count = 0; 1489754ec7b1SSatish Balay } 1490416022c9SBarry Smith c->nz = a->nz; 1491416022c9SBarry Smith c->maxnz = a->maxnz; 1492416022c9SBarry Smith c->solve_work = 0; 149376dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1494754ec7b1SSatish Balay 1495416022c9SBarry Smith *B = C; 149617ab2063SBarry Smith return 0; 149717ab2063SBarry Smith } 149817ab2063SBarry Smith 1499416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 150017ab2063SBarry Smith { 1501416022c9SBarry Smith Mat_SeqAIJ *a; 1502416022c9SBarry Smith Mat B; 150317699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1504*bcd2baecSBarry Smith MPI_Comm comm; 1505*bcd2baecSBarry Smith ViewerType vtype; 150617ab2063SBarry Smith 1507*bcd2baecSBarry Smith ierr = ViewerGetType(bview,&vtype); CHKERRQ(ierr); 1508*bcd2baecSBarry Smith if (vtype != BINARY_FILE_VIEWER) SETERRQ(1,"MatLoad_SeqAIJ:Binary only"); 1509*bcd2baecSBarry Smith 1510*bcd2baecSBarry Smith PetscObjectGetComm((PetscObject) bview,&comm); 151117699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 151217699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 151317ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1514416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 151548d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 151617ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 151717ab2063SBarry Smith 151817ab2063SBarry Smith /* read in row lengths */ 15190452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1520416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 152117ab2063SBarry Smith 152217ab2063SBarry Smith /* create our matrix */ 1523416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1524416022c9SBarry Smith B = *A; 1525416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1526416022c9SBarry Smith shift = a->indexshift; 152717ab2063SBarry Smith 152817ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1529416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 153017ab2063SBarry Smith if (shift) { 153117ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1532416022c9SBarry Smith a->j[i] += 1; 153317ab2063SBarry Smith } 153417ab2063SBarry Smith } 153517ab2063SBarry Smith 153617ab2063SBarry Smith /* read in nonzero values */ 1537416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 153817ab2063SBarry Smith 153917ab2063SBarry Smith /* set matrix "i" values */ 1540416022c9SBarry Smith a->i[0] = -shift; 154117ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1542416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1543416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 154417ab2063SBarry Smith } 15450452661fSBarry Smith PetscFree(rowlengths); 154617ab2063SBarry Smith 1547416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1548416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 154917ab2063SBarry Smith return 0; 155017ab2063SBarry Smith } 155117ab2063SBarry Smith 15527264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg) 15537264ac53SSatish Balay { 15547264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 15557264ac53SSatish Balay 1556*bcd2baecSBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type"); 15577264ac53SSatish Balay 15587264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 15597264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1560*bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 1561*bcd2baecSBarry Smith *flg = 0 ; return 0; 1562*bcd2baecSBarry Smith } 15637264ac53SSatish Balay 15647264ac53SSatish Balay /* if the a->i are the same */ 1565*bcd2baecSBarry Smith if(PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 15667264ac53SSatish Balay *flg = 0 ; return 0; 15677264ac53SSatish Balay } 15687264ac53SSatish Balay 15697264ac53SSatish Balay /* if a->j are the same */ 1570*bcd2baecSBarry Smith if(PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 1571*bcd2baecSBarry Smith *flg = 0 ; return 0; 1572*bcd2baecSBarry Smith } 1573*bcd2baecSBarry Smith 1574*bcd2baecSBarry Smith /* if a->a are the same */ 1575*bcd2baecSBarry Smith if(PetscMemcmp(a, b->a, (a->nz)*sizeof(Scalar))) { 15767264ac53SSatish Balay *flg = 0 ; return 0; 15777264ac53SSatish Balay } 15787264ac53SSatish Balay *flg =1 ; 15797264ac53SSatish Balay return 0; 15807264ac53SSatish Balay 15817264ac53SSatish Balay } 1582