1*8d195f9aSBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*8d195f9aSBarry Smith static char vcid[] = "$Id: aij.c,v 1.234 1997/08/28 14:23:17 curfman Exp bsmith $"; 417ab2063SBarry Smith #endif 517ab2063SBarry Smith 6d5d45c9bSBarry Smith /* 73369ce9aSBarry Smith Defines the basic matrix operations for the AIJ (compressed row) 8d5d45c9bSBarry Smith matrix storage format. 9d5d45c9bSBarry Smith */ 103369ce9aSBarry Smith 113369ce9aSBarry Smith #include "pinclude/pviewer.h" 123369ce9aSBarry Smith #include "sys.h" 1370f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 14f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 15f5eb4b81SSatish Balay #include "src/inline/spops.h" 16*8d195f9aSBarry Smith #include "src/inline/dot.h" 17f5eb4b81SSatish Balay #include "src/inline/bitarray.h" 1817ab2063SBarry Smith 19a2ce50c7SBarry Smith /* 20a2ce50c7SBarry Smith Basic AIJ format ILU based on drop tolerance 21a2ce50c7SBarry Smith */ 225615d1e5SSatish Balay #undef __FUNC__ 235615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ" 24a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact) 25a2ce50c7SBarry Smith { 26a2ce50c7SBarry Smith /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */ 27a2ce50c7SBarry Smith int ierr = 1; 28a2ce50c7SBarry Smith 29e3372554SBarry Smith SETERRQ(ierr,0,"Not implemented"); 30a2ce50c7SBarry Smith } 31a2ce50c7SBarry Smith 32bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 3317ab2063SBarry Smith 345615d1e5SSatish Balay #undef __FUNC__ 355615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ" 368f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja, 376945ee14SBarry Smith PetscTruth *done) 3817ab2063SBarry Smith { 39416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 406945ee14SBarry Smith int ierr,i,ishift; 4117ab2063SBarry Smith 4231625ec5SSatish Balay *m = A->m; 436945ee14SBarry Smith if (!ia) return 0; 446945ee14SBarry Smith ishift = a->indexshift; 456945ee14SBarry Smith if (symmetric) { 4631625ec5SSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 476945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 4831625ec5SSatish Balay int nz = a->i[a->m]; 493b2fbd54SBarry Smith /* malloc space and subtract 1 from i and j indices */ 5031625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 513b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 523b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1; 5331625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1; 546945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 5531625ec5SSatish Balay int nz = a->i[a->m] + 1; 563b2fbd54SBarry Smith /* malloc space and add 1 to i and j indices */ 5731625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 583b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 593b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1; 6031625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1; 616945ee14SBarry Smith } else { 626945ee14SBarry Smith *ia = a->i; *ja = a->j; 63a2ce50c7SBarry Smith } 64a2ce50c7SBarry Smith 65a2744918SBarry Smith return 0; 66a2744918SBarry Smith } 67a2744918SBarry Smith 685615d1e5SSatish Balay #undef __FUNC__ 695615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ" 708f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 716945ee14SBarry Smith PetscTruth *done) 726945ee14SBarry Smith { 736945ee14SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 743b2fbd54SBarry Smith int ishift = a->indexshift; 756945ee14SBarry Smith 766945ee14SBarry Smith if (!ia) return 0; 773b2fbd54SBarry Smith if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) { 786945ee14SBarry Smith PetscFree(*ia); 796945ee14SBarry Smith PetscFree(*ja); 80bcd2baecSBarry Smith } 8117ab2063SBarry Smith return 0; 8217ab2063SBarry Smith } 8317ab2063SBarry Smith 845615d1e5SSatish Balay #undef __FUNC__ 855615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ" 8643a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 873b2fbd54SBarry Smith PetscTruth *done) 883b2fbd54SBarry Smith { 893b2fbd54SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 90a93ec695SBarry Smith int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m; 91a93ec695SBarry Smith int nz = a->i[m]+ishift,row,*jj,mr,col; 923b2fbd54SBarry Smith 933b2fbd54SBarry Smith *nn = A->n; 943b2fbd54SBarry Smith if (!ia) return 0; 953b2fbd54SBarry Smith if (symmetric) { 96179192dfSSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 973b2fbd54SBarry Smith } else { 9861d2ded1SBarry Smith collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths); 993b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1003b2fbd54SBarry Smith cia = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia); 101a93ec695SBarry Smith cja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja); 1023b2fbd54SBarry Smith jj = a->j; 1033b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) { 1043b2fbd54SBarry Smith collengths[jj[i] + ishift]++; 1053b2fbd54SBarry Smith } 1063b2fbd54SBarry Smith cia[0] = oshift; 1073b2fbd54SBarry Smith for ( i=0; i<n; i++) { 1083b2fbd54SBarry Smith cia[i+1] = cia[i] + collengths[i]; 1093b2fbd54SBarry Smith } 1103b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1113b2fbd54SBarry Smith jj = a->j; 112a93ec695SBarry Smith for ( row=0; row<m; row++ ) { 113a93ec695SBarry Smith mr = a->i[row+1] - a->i[row]; 114a93ec695SBarry Smith for ( i=0; i<mr; i++ ) { 1153b2fbd54SBarry Smith col = *jj++ + ishift; 1163b2fbd54SBarry Smith cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1173b2fbd54SBarry Smith } 1183b2fbd54SBarry Smith } 1193b2fbd54SBarry Smith PetscFree(collengths); 1203b2fbd54SBarry Smith *ia = cia; *ja = cja; 1213b2fbd54SBarry Smith } 1223b2fbd54SBarry Smith 1233b2fbd54SBarry Smith return 0; 1243b2fbd54SBarry Smith } 1253b2fbd54SBarry Smith 1265615d1e5SSatish Balay #undef __FUNC__ 1275615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ" 12843a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia, 1293b2fbd54SBarry Smith int **ja,PetscTruth *done) 1303b2fbd54SBarry Smith { 1313b2fbd54SBarry Smith if (!ia) return 0; 1323b2fbd54SBarry Smith 1333b2fbd54SBarry Smith PetscFree(*ia); 1343b2fbd54SBarry Smith PetscFree(*ja); 1353b2fbd54SBarry Smith 1363b2fbd54SBarry Smith return 0; 1373b2fbd54SBarry Smith } 1383b2fbd54SBarry Smith 139227d817aSBarry Smith #define CHUNKSIZE 15 14017ab2063SBarry Smith 1415615d1e5SSatish Balay #undef __FUNC__ 1425615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ" 14361d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 14417ab2063SBarry Smith { 145416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 146416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 1474b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 148d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 149416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 15017ab2063SBarry Smith 15117ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 152416022c9SBarry Smith row = im[k]; 1533b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 154e3372554SBarry Smith if (row < 0) SETERRQ(1,0,"Negative row"); 155e3372554SBarry Smith if (row >= a->m) SETERRQ(1,0,"Row too large"); 1563b2fbd54SBarry Smith #endif 15717ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 15817ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 159416022c9SBarry Smith low = 0; 16017ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 1613b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 162e3372554SBarry Smith if (in[l] < 0) SETERRQ(1,0,"Negative column"); 163e3372554SBarry Smith if (in[l] >= a->n) SETERRQ(1,0,"Column too large"); 1643b2fbd54SBarry Smith #endif 1654b0e389bSBarry Smith col = in[l] - shift; 1664b0e389bSBarry Smith if (roworiented) { 1674b0e389bSBarry Smith value = *v++; 1684b0e389bSBarry Smith } 1694b0e389bSBarry Smith else { 1704b0e389bSBarry Smith value = v[k + l*m]; 1714b0e389bSBarry Smith } 172416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 173416022c9SBarry Smith while (high-low > 5) { 174416022c9SBarry Smith t = (low+high)/2; 175416022c9SBarry Smith if (rp[t] > col) high = t; 176416022c9SBarry Smith else low = t; 17717ab2063SBarry Smith } 178416022c9SBarry Smith for ( i=low; i<high; i++ ) { 17917ab2063SBarry Smith if (rp[i] > col) break; 18017ab2063SBarry Smith if (rp[i] == col) { 181416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 18217ab2063SBarry Smith else ap[i] = value; 18317ab2063SBarry Smith goto noinsert; 18417ab2063SBarry Smith } 18517ab2063SBarry Smith } 186c2653b3dSLois Curfman McInnes if (nonew == 1) goto noinsert; 18711ebbc71SLois Curfman McInnes else if (nonew == -1) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); 18817ab2063SBarry Smith if (nrow >= rmax) { 18917ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 190416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 19117ab2063SBarry Smith Scalar *new_a; 19217ab2063SBarry Smith 19311ebbc71SLois Curfman McInnes if (nonew == -2) SETERRQ(1,0,"Inserting a new nonzero in the matrix"); 19496854ed6SLois Curfman McInnes 19517ab2063SBarry Smith /* malloc new storage space */ 196416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1970452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 19817ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 19917ab2063SBarry Smith new_i = new_j + new_nz; 20017ab2063SBarry Smith 20117ab2063SBarry Smith /* copy over old data into new slots */ 20217ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 203416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 204416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 205416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 206416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 20717ab2063SBarry Smith len*sizeof(int)); 208416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 209416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 21017ab2063SBarry Smith len*sizeof(Scalar)); 21117ab2063SBarry Smith /* free up old matrix storage */ 2120452661fSBarry Smith PetscFree(a->a); 2130452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 214416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 215416022c9SBarry Smith a->singlemalloc = 1; 21617ab2063SBarry Smith 21717ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 218416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 219416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 220416022c9SBarry Smith a->maxnz += CHUNKSIZE; 221b810aeb4SBarry Smith a->reallocs++; 22217ab2063SBarry Smith } 223416022c9SBarry Smith N = nrow++ - 1; a->nz++; 224416022c9SBarry Smith /* shift up all the later entries in this row */ 225416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 22617ab2063SBarry Smith rp[ii+1] = rp[ii]; 22717ab2063SBarry Smith ap[ii+1] = ap[ii]; 22817ab2063SBarry Smith } 22917ab2063SBarry Smith rp[i] = col; 23017ab2063SBarry Smith ap[i] = value; 23117ab2063SBarry Smith noinsert:; 232416022c9SBarry Smith low = i + 1; 23317ab2063SBarry Smith } 23417ab2063SBarry Smith ailen[row] = nrow; 23517ab2063SBarry Smith } 23617ab2063SBarry Smith return 0; 23717ab2063SBarry Smith } 23817ab2063SBarry Smith 2395615d1e5SSatish Balay #undef __FUNC__ 2405615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ" 2418f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 2427eb43aa7SLois Curfman McInnes { 2437eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 244b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 2457eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 2467eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 2477eb43aa7SLois Curfman McInnes 2487eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 2497eb43aa7SLois Curfman McInnes row = im[k]; 250e3372554SBarry Smith if (row < 0) SETERRQ(1,0,"Negative row"); 251e3372554SBarry Smith if (row >= a->m) SETERRQ(1,0,"Row too large"); 2527eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 2537eb43aa7SLois Curfman McInnes nrow = ailen[row]; 2547eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 255e3372554SBarry Smith if (in[l] < 0) SETERRQ(1,0,"Negative column"); 256e3372554SBarry Smith if (in[l] >= a->n) SETERRQ(1,0,"Column too large"); 2577eb43aa7SLois Curfman McInnes col = in[l] - shift; 2587eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 2597eb43aa7SLois Curfman McInnes while (high-low > 5) { 2607eb43aa7SLois Curfman McInnes t = (low+high)/2; 2617eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 2627eb43aa7SLois Curfman McInnes else low = t; 2637eb43aa7SLois Curfman McInnes } 2647eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 2657eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2667eb43aa7SLois Curfman McInnes if (rp[i] == col) { 267b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2687eb43aa7SLois Curfman McInnes goto finished; 2697eb43aa7SLois Curfman McInnes } 2707eb43aa7SLois Curfman McInnes } 271b49de8d1SLois Curfman McInnes *v++ = zero; 2727eb43aa7SLois Curfman McInnes finished:; 2737eb43aa7SLois Curfman McInnes } 2747eb43aa7SLois Curfman McInnes } 2757eb43aa7SLois Curfman McInnes return 0; 2767eb43aa7SLois Curfman McInnes } 2777eb43aa7SLois Curfman McInnes 27817ab2063SBarry Smith 2795615d1e5SSatish Balay #undef __FUNC__ 2805615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary" 2818f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 28217ab2063SBarry Smith { 283416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 284416022c9SBarry Smith int i, fd, *col_lens, ierr; 28517ab2063SBarry Smith 28690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2870452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 288416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 289416022c9SBarry Smith col_lens[1] = a->m; 290416022c9SBarry Smith col_lens[2] = a->n; 291416022c9SBarry Smith col_lens[3] = a->nz; 292416022c9SBarry Smith 293416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 294416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 295416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 29617ab2063SBarry Smith } 29777c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 2980452661fSBarry Smith PetscFree(col_lens); 299416022c9SBarry Smith 300416022c9SBarry Smith /* store column indices (zero start index) */ 301416022c9SBarry Smith if (a->indexshift) { 302416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 30317ab2063SBarry Smith } 30477c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 305416022c9SBarry Smith if (a->indexshift) { 306416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 30717ab2063SBarry Smith } 308416022c9SBarry Smith 309416022c9SBarry Smith /* store nonzero values */ 31077c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 31117ab2063SBarry Smith return 0; 31217ab2063SBarry Smith } 313416022c9SBarry Smith 3145615d1e5SSatish Balay #undef __FUNC__ 3155615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII" 3168f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 317416022c9SBarry Smith { 318416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 319496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 32017ab2063SBarry Smith FILE *fd; 32117ab2063SBarry Smith char *outputname; 32217ab2063SBarry Smith 32390ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 324416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 32590ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 326a93ec695SBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 32795e01e2fSLois Curfman McInnes return 0; 32895e01e2fSLois Curfman McInnes } 329a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 330496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 331496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 332496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 33395e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 33495e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 33517ab2063SBarry Smith } 336a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 337d00d2cf4SBarry Smith int nofinalvalue = 0; 338d00d2cf4SBarry Smith if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) { 339d00d2cf4SBarry Smith nofinalvalue = 1; 340d00d2cf4SBarry Smith } 341416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3424e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 343d00d2cf4SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue); 34417ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 34517ab2063SBarry Smith 34617ab2063SBarry Smith for (i=0; i<m; i++) { 347416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 34817ab2063SBarry Smith #if defined(PETSC_COMPLEX) 3496945ee14SBarry Smith fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]), 350416022c9SBarry Smith imag(a->a[j])); 35117ab2063SBarry Smith #else 3527a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 35317ab2063SBarry Smith #endif 35417ab2063SBarry Smith } 35517ab2063SBarry Smith } 356d00d2cf4SBarry Smith if (nofinalvalue) { 357d00d2cf4SBarry Smith fprintf(fd,"%d %d %18.16e\n", m, a->n, 0.0); 358d00d2cf4SBarry Smith } 35917ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 36017ab2063SBarry Smith } 361a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 36244cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 36344cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 36444cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 36544cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 366766eeae4SLois Curfman McInnes if (imag(a->a[j]) > 0.0 && real(a->a[j]) != 0.0) 36744cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 368766eeae4SLois Curfman McInnes else if (imag(a->a[j]) < 0.0 && real(a->a[j]) != 0.0) 369766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j])); 37044cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 37144cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 37244cd7ae7SLois Curfman McInnes #else 37344cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 37444cd7ae7SLois Curfman McInnes #endif 37544cd7ae7SLois Curfman McInnes } 37644cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 37744cd7ae7SLois Curfman McInnes } 37844cd7ae7SLois Curfman McInnes } 379496be53dSLois Curfman McInnes else if (format == VIEWER_FORMAT_ASCII_SYMMODU) { 380496be53dSLois Curfman McInnes int nzd=0, fshift=1, *sptr; 3812e44a96cSLois Curfman McInnes sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr); 382496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 383496be53dSLois Curfman McInnes sptr[i] = nzd+1; 384496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 385496be53dSLois Curfman McInnes if (a->j[j] >= i) { 386496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX) 387496be53dSLois Curfman McInnes if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) nzd++; 388496be53dSLois Curfman McInnes #else 389496be53dSLois Curfman McInnes if (a->a[j] != 0.0) nzd++; 390496be53dSLois Curfman McInnes #endif 391496be53dSLois Curfman McInnes } 392496be53dSLois Curfman McInnes } 393496be53dSLois Curfman McInnes } 3942e44a96cSLois Curfman McInnes sptr[m] = nzd+1; 395496be53dSLois Curfman McInnes fprintf(fd," %d %d\n\n",m,nzd); 3962e44a96cSLois Curfman McInnes for ( i=0; i<m+1; i+=6 ) { 3972e44a96cSLois Curfman McInnes if (i+4<m) fprintf(fd," %d %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]); 3982e44a96cSLois Curfman McInnes else if (i+3<m) fprintf(fd," %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]); 3992e44a96cSLois Curfman McInnes else if (i+2<m) fprintf(fd," %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]); 4002e44a96cSLois Curfman McInnes else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]); 4012e44a96cSLois Curfman McInnes else if (i<m) fprintf(fd," %d %d\n",sptr[i],sptr[i+1]); 4027272d637SLois Curfman McInnes else fprintf(fd," %d\n",sptr[i]); 403496be53dSLois Curfman McInnes } 404496be53dSLois Curfman McInnes fprintf(fd,"\n"); 405496be53dSLois Curfman McInnes PetscFree(sptr); 406496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 407496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 408496be53dSLois Curfman McInnes if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift); 409496be53dSLois Curfman McInnes } 410496be53dSLois Curfman McInnes fprintf(fd,"\n"); 411496be53dSLois Curfman McInnes } 412496be53dSLois Curfman McInnes fprintf(fd,"\n"); 413496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 414496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 415496be53dSLois Curfman McInnes if (a->j[j] >= i) { 416496be53dSLois Curfman McInnes #if defined(PETSC_COMPLEX) 417496be53dSLois Curfman McInnes if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) 418496be53dSLois Curfman McInnes fprintf(fd," %18.16e %18.16e ",real(a->a[j]),imag(a->a[j])); 419496be53dSLois Curfman McInnes #else 420496be53dSLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]); 421496be53dSLois Curfman McInnes #endif 422496be53dSLois Curfman McInnes } 423496be53dSLois Curfman McInnes } 424496be53dSLois Curfman McInnes fprintf(fd,"\n"); 425496be53dSLois Curfman McInnes } 426496be53dSLois Curfman McInnes } 42717ab2063SBarry Smith else { 42817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 42917ab2063SBarry Smith fprintf(fd,"row %d:",i); 430416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 43117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 432766eeae4SLois Curfman McInnes if (imag(a->a[j]) > 0.0) { 433416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 434766eeae4SLois Curfman McInnes } else if (imag(a->a[j]) < 0.0) { 435766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j])); 43617ab2063SBarry Smith } 43717ab2063SBarry Smith else { 438416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 43917ab2063SBarry Smith } 44017ab2063SBarry Smith #else 441416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 44217ab2063SBarry Smith #endif 44317ab2063SBarry Smith } 44417ab2063SBarry Smith fprintf(fd,"\n"); 44517ab2063SBarry Smith } 44617ab2063SBarry Smith } 44717ab2063SBarry Smith fflush(fd); 448416022c9SBarry Smith return 0; 449416022c9SBarry Smith } 450416022c9SBarry Smith 4515615d1e5SSatish Balay #undef __FUNC__ 4525615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw" 4538f6be9afSLois Curfman McInnes extern int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 454416022c9SBarry Smith { 455416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 456cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 4570513a670SBarry Smith int format; 45894a9d846SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv = 0.0; 459bcd2baecSBarry Smith Draw draw; 460cddf8d76SBarry Smith DrawButton button; 46119bcc07fSBarry Smith PetscTruth isnull; 462cddf8d76SBarry Smith 4630513a670SBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 4640513a670SBarry Smith ierr = DrawClear(draw); CHKERRQ(ierr); 4650513a670SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 46619bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 46719bcc07fSBarry Smith 468416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 469416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 470416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 471416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 4720513a670SBarry Smith 4730513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 4740513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 475cddf8d76SBarry Smith color = DRAW_BLUE; 476416022c9SBarry Smith for ( i=0; i<m; i++ ) { 477cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 478416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 479cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 480cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 481cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 482cddf8d76SBarry Smith #else 483cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 484cddf8d76SBarry Smith #endif 485cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 486cddf8d76SBarry Smith } 487cddf8d76SBarry Smith } 488cddf8d76SBarry Smith color = DRAW_CYAN; 489cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 490cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 491cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 492cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 493cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 494cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 495cddf8d76SBarry Smith } 496cddf8d76SBarry Smith } 497cddf8d76SBarry Smith color = DRAW_RED; 498cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 499cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 500cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 501cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 502cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 503cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 504cddf8d76SBarry Smith #else 505cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 506cddf8d76SBarry Smith #endif 507cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 508416022c9SBarry Smith } 509416022c9SBarry Smith } 5100513a670SBarry Smith } else { 5110513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 5120513a670SBarry Smith /* first determine max of all nonzero values */ 5130513a670SBarry Smith int nz = a->nz,count; 5140513a670SBarry Smith Draw popup; 5150513a670SBarry Smith 5160513a670SBarry Smith for ( i=0; i<nz; i++ ) { 5170513a670SBarry Smith if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 5180513a670SBarry Smith } 5190513a670SBarry Smith ierr = DrawCreatePopUp(draw,&popup); CHKERRQ(ierr); 5200513a670SBarry Smith ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr); 5210513a670SBarry Smith count = 0; 5220513a670SBarry Smith for ( i=0; i<m; i++ ) { 5230513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 5240513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 5250513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5260513a670SBarry Smith color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv); 5270513a670SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 5280513a670SBarry Smith count++; 5290513a670SBarry Smith } 5300513a670SBarry Smith } 5310513a670SBarry Smith } 532416022c9SBarry Smith DrawFlush(draw); 533cddf8d76SBarry Smith DrawGetPause(draw,&pause); 534cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 535cddf8d76SBarry Smith 536cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 5376945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 538cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 539cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 540cddf8d76SBarry Smith DrawClear(draw); 541cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 542cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 543cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 544cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 545cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 546cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 547cddf8d76SBarry Smith w *= scale; h *= scale; 548cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 5490513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 5500513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 551cddf8d76SBarry Smith color = DRAW_BLUE; 552cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 553cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 554cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 555cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 556cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 557cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 558cddf8d76SBarry Smith #else 559cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 560cddf8d76SBarry Smith #endif 561cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 562cddf8d76SBarry Smith } 563cddf8d76SBarry Smith } 564cddf8d76SBarry Smith color = DRAW_CYAN; 565cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 566cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 567cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 568cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 569cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 570cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 571cddf8d76SBarry Smith } 572cddf8d76SBarry Smith } 573cddf8d76SBarry Smith color = DRAW_RED; 574cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 575cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 576cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 577cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 578cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 579cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 580cddf8d76SBarry Smith #else 581cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 582cddf8d76SBarry Smith #endif 583cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 584cddf8d76SBarry Smith } 585cddf8d76SBarry Smith } 5860513a670SBarry Smith } else { 5870513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 5880513a670SBarry Smith int count = 0; 5890513a670SBarry Smith for ( i=0; i<m; i++ ) { 5900513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 5910513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 5920513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5930513a670SBarry Smith color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv); 5940513a670SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 5950513a670SBarry Smith count++; 5960513a670SBarry Smith } 5970513a670SBarry Smith } 5980513a670SBarry Smith } 5990513a670SBarry Smith 6006945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 601cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 602cddf8d76SBarry Smith } 603416022c9SBarry Smith return 0; 604416022c9SBarry Smith } 605416022c9SBarry Smith 6065615d1e5SSatish Balay #undef __FUNC__ 607d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ" 6088f6be9afSLois Curfman McInnes int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 609416022c9SBarry Smith { 610416022c9SBarry Smith Mat A = (Mat) obj; 611416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 612bcd2baecSBarry Smith ViewerType vtype; 613bcd2baecSBarry Smith int ierr; 614416022c9SBarry Smith 615bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 616bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 617416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 618416022c9SBarry Smith } 619bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 620416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 621416022c9SBarry Smith } 622bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 623416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 624416022c9SBarry Smith } 625bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 626bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 62717ab2063SBarry Smith } 62817ab2063SBarry Smith return 0; 62917ab2063SBarry Smith } 63019bcc07fSBarry Smith 631c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 6325615d1e5SSatish Balay #undef __FUNC__ 6335615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ" 6348f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 63517ab2063SBarry Smith { 636416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 63741c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 63843ee02c3SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0; 639416022c9SBarry Smith Scalar *aa = a->a, *ap; 64017ab2063SBarry Smith 6416d4a8577SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) return 0; 64217ab2063SBarry Smith 64343ee02c3SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 64417ab2063SBarry Smith for ( i=1; i<m; i++ ) { 645416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 64617ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 64794a9d846SBarry Smith rmax = PetscMax(rmax,ailen[i]); 64817ab2063SBarry Smith if (fshift) { 649416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 65017ab2063SBarry Smith N = ailen[i]; 65117ab2063SBarry Smith for ( j=0; j<N; j++ ) { 65217ab2063SBarry Smith ip[j-fshift] = ip[j]; 65317ab2063SBarry Smith ap[j-fshift] = ap[j]; 65417ab2063SBarry Smith } 65517ab2063SBarry Smith } 65617ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 65717ab2063SBarry Smith } 65817ab2063SBarry Smith if (m) { 65917ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 66017ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 66117ab2063SBarry Smith } 66217ab2063SBarry Smith /* reset ilen and imax for each row */ 66317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 66417ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 66517ab2063SBarry Smith } 666416022c9SBarry Smith a->nz = ai[m] + shift; 66717ab2063SBarry Smith 66817ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 669416022c9SBarry Smith if (fshift && a->diag) { 6700452661fSBarry Smith PetscFree(a->diag); 671416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 672416022c9SBarry Smith a->diag = 0; 67317ab2063SBarry Smith } 6744e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 6754e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 6764e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 677b810aeb4SBarry Smith a->reallocs); 67894a9d846SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax); 679dd5f02e7SSatish Balay a->reallocs = 0; 6804e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 6814e220ebcSLois Curfman McInnes 68276dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 68341c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 68417ab2063SBarry Smith return 0; 68517ab2063SBarry Smith } 68617ab2063SBarry Smith 6875615d1e5SSatish Balay #undef __FUNC__ 6885615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ" 6898f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A) 69017ab2063SBarry Smith { 691416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 692cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 69317ab2063SBarry Smith return 0; 69417ab2063SBarry Smith } 695416022c9SBarry Smith 6965615d1e5SSatish Balay #undef __FUNC__ 6975615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ" 69817ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 69917ab2063SBarry Smith { 700416022c9SBarry Smith Mat A = (Mat) obj; 701416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 702d5d45c9bSBarry Smith 70317ab2063SBarry Smith #if defined(PETSC_LOG) 704416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 70517ab2063SBarry Smith #endif 7060452661fSBarry Smith PetscFree(a->a); 7070452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 7080452661fSBarry Smith if (a->diag) PetscFree(a->diag); 7090452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 7100452661fSBarry Smith if (a->imax) PetscFree(a->imax); 7110452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 71276dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 7130452661fSBarry Smith PetscFree(a); 714eed86810SBarry Smith 715f2655603SLois Curfman McInnes PLogObjectDestroy(A); 716f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 71717ab2063SBarry Smith return 0; 71817ab2063SBarry Smith } 71917ab2063SBarry Smith 7205615d1e5SSatish Balay #undef __FUNC__ 7215615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ" 7228f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A) 72317ab2063SBarry Smith { 72417ab2063SBarry Smith return 0; 72517ab2063SBarry Smith } 72617ab2063SBarry Smith 7275615d1e5SSatish Balay #undef __FUNC__ 7285615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ" 7298f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op) 73017ab2063SBarry Smith { 731416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7326d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 7336d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 7346d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 735219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 7366d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 737c2653b3dSLois Curfman McInnes else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 73896854ed6SLois Curfman McInnes else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 7396d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 7406d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 741219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 7426d4a8577SBarry Smith op == MAT_SYMMETRIC || 7436d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 74490f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 7452b362799SSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES) 74694a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 7476d4a8577SBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) 748e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");} 7496d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 7506d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 7516d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 7526d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 7536d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 754e2f28af5SBarry Smith else 755e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"unknown option");} 75617ab2063SBarry Smith return 0; 75717ab2063SBarry Smith } 75817ab2063SBarry Smith 7595615d1e5SSatish Balay #undef __FUNC__ 7605615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ" 7618f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 76217ab2063SBarry Smith { 763416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 764416022c9SBarry Smith int i,j, n,shift = a->indexshift; 76517ab2063SBarry Smith Scalar *x, zero = 0.0; 76617ab2063SBarry Smith 76717ab2063SBarry Smith VecSet(&zero,v); 76890f02eecSBarry Smith VecGetArray_Fast(v,x); VecGetLocalSize(v,&n); 769e3372554SBarry Smith if (n != a->m) SETERRQ(1,0,"Nonconforming matrix and vector"); 770416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 771416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 772416022c9SBarry Smith if (a->j[j]+shift == i) { 773416022c9SBarry Smith x[i] = a->a[j]; 77417ab2063SBarry Smith break; 77517ab2063SBarry Smith } 77617ab2063SBarry Smith } 77717ab2063SBarry Smith } 77817ab2063SBarry Smith return 0; 77917ab2063SBarry Smith } 78017ab2063SBarry Smith 78117ab2063SBarry Smith /* -------------------------------------------------------*/ 78217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 78317ab2063SBarry Smith /* -------------------------------------------------------*/ 7845615d1e5SSatish Balay #undef __FUNC__ 7855615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ" 78644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 78717ab2063SBarry Smith { 788416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 78917ab2063SBarry Smith Scalar *x, *y, *v, alpha; 790416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 79117ab2063SBarry Smith 79290f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 793cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 79417ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 79517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 796416022c9SBarry Smith idx = a->j + a->i[i] + shift; 797416022c9SBarry Smith v = a->a + a->i[i] + shift; 798416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 79917ab2063SBarry Smith alpha = x[i]; 80017ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 80117ab2063SBarry Smith } 802416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 80317ab2063SBarry Smith return 0; 80417ab2063SBarry Smith } 805d5d45c9bSBarry Smith 8065615d1e5SSatish Balay #undef __FUNC__ 8075615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ" 80844cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 80917ab2063SBarry Smith { 810416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 81117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 812416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 81317ab2063SBarry Smith 81490f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 81517ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 81617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 81717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 818416022c9SBarry Smith idx = a->j + a->i[i] + shift; 819416022c9SBarry Smith v = a->a + a->i[i] + shift; 820416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 82117ab2063SBarry Smith alpha = x[i]; 82217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 82317ab2063SBarry Smith } 82490f02eecSBarry Smith PLogFlops(2*a->nz); 82517ab2063SBarry Smith return 0; 82617ab2063SBarry Smith } 82717ab2063SBarry Smith 8285615d1e5SSatish Balay #undef __FUNC__ 8295615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ" 83044cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 83117ab2063SBarry Smith { 832416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 83317ab2063SBarry Smith Scalar *x, *y, *v, sum; 8349ea0dfa2SSatish Balay int m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j; 83517ab2063SBarry Smith 83690f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 83717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 838416022c9SBarry Smith idx = a->j; 8399b88f4a7SBarry Smith v = a->a + shift; /* shift for Fortran start by 1 indexing */ 840416022c9SBarry Smith ii = a->i; 841*8d195f9aSBarry Smith #if defined(USE_FORTRAN_KERNELS) 842*8d195f9aSBarry Smith fortranmultaij_(&m,x,ii,idx,v,y); 843*8d195f9aSBarry Smith #else 84417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8459ea0dfa2SSatish Balay jrow = ii[i]; 8469ea0dfa2SSatish Balay n = ii[i+1] - jrow; 84717ab2063SBarry Smith sum = 0.0; 8489ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8499ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8509ea0dfa2SSatish Balay } 85117ab2063SBarry Smith y[i] = sum; 85217ab2063SBarry Smith } 853*8d195f9aSBarry Smith #endif 854416022c9SBarry Smith PLogFlops(2*a->nz - m); 85517ab2063SBarry Smith return 0; 85617ab2063SBarry Smith } 85717ab2063SBarry Smith 8585615d1e5SSatish Balay #undef __FUNC__ 8595615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ" 86044cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 86117ab2063SBarry Smith { 862416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 86317ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 8649ea0dfa2SSatish Balay int m = a->m, n, i, *idx, shift = a->indexshift,*ii,jrow,j; 8659ea0dfa2SSatish Balay 86617ab2063SBarry Smith 86790f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z); 86817ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 869cddf8d76SBarry Smith idx = a->j; 8709b88f4a7SBarry Smith v = a->a + shift; /* shift for Fortran start by 1 indexing */ 871cddf8d76SBarry Smith ii = a->i; 87217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8739ea0dfa2SSatish Balay jrow = ii[i]; 8749ea0dfa2SSatish Balay n = ii[i+1] - jrow; 87517ab2063SBarry Smith sum = y[i]; 8769ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8779ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8789ea0dfa2SSatish Balay } 87917ab2063SBarry Smith z[i] = sum; 88017ab2063SBarry Smith } 881416022c9SBarry Smith PLogFlops(2*a->nz); 88217ab2063SBarry Smith return 0; 88317ab2063SBarry Smith } 88417ab2063SBarry Smith 88517ab2063SBarry Smith /* 88617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 88717ab2063SBarry Smith */ 88817ab2063SBarry Smith 8895615d1e5SSatish Balay #undef __FUNC__ 8905615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ" 891416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 89217ab2063SBarry Smith { 893416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 894416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 89517ab2063SBarry Smith 8960452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 897416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 898416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 899416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 900416022c9SBarry Smith if (a->j[j]+shift == i) { 90117ab2063SBarry Smith diag[i] = j - shift; 90217ab2063SBarry Smith break; 90317ab2063SBarry Smith } 90417ab2063SBarry Smith } 90517ab2063SBarry Smith } 906416022c9SBarry Smith a->diag = diag; 90717ab2063SBarry Smith return 0; 90817ab2063SBarry Smith } 90917ab2063SBarry Smith 9105615d1e5SSatish Balay #undef __FUNC__ 9115615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ" 91244cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 91317ab2063SBarry Smith double fshift,int its,Vec xx) 91417ab2063SBarry Smith { 915416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 916416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 917d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 91817ab2063SBarry Smith 91990f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b); 920d00d2cf4SBarry Smith if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);} 921416022c9SBarry Smith diag = a->diag; 922416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 92317ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 92417ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 92517ab2063SBarry Smith bs = b + shift; 92617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 927416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 928416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 929416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 930416022c9SBarry Smith v = a->a + diag[i] + (!shift); 93117ab2063SBarry Smith sum = b[i]*d/omega; 93217ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 93317ab2063SBarry Smith x[i] = sum; 93417ab2063SBarry Smith } 93517ab2063SBarry Smith return 0; 93617ab2063SBarry Smith } 93717ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 938e3372554SBarry Smith SETERRQ(1,0,"SOR_APPLY_LOWER is not done"); 93917ab2063SBarry Smith } 940416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 94117ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 94217ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 94317ab2063SBarry Smith 94417ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 94517ab2063SBarry Smith 94617ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 94717ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 94817ab2063SBarry Smith is the relaxation factor. 94917ab2063SBarry Smith */ 9500452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 95117ab2063SBarry Smith scale = (2.0/omega) - 1.0; 95217ab2063SBarry Smith 95317ab2063SBarry Smith /* x = (E + U)^{-1} b */ 95417ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 955416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 956416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 957416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 958416022c9SBarry Smith v = a->a + diag[i] + (!shift); 95917ab2063SBarry Smith sum = b[i]; 96017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 96117ab2063SBarry Smith x[i] = omega*(sum/d); 96217ab2063SBarry Smith } 96317ab2063SBarry Smith 96417ab2063SBarry Smith /* t = b - (2*E - D)x */ 965416022c9SBarry Smith v = a->a; 96617ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 96717ab2063SBarry Smith 96817ab2063SBarry Smith /* t = (E + L)^{-1}t */ 969416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 970416022c9SBarry Smith diag = a->diag; 97117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 972416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 973416022c9SBarry Smith n = diag[i] - a->i[i]; 974416022c9SBarry Smith idx = a->j + a->i[i] + shift; 975416022c9SBarry Smith v = a->a + a->i[i] + shift; 97617ab2063SBarry Smith sum = t[i]; 97717ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 97817ab2063SBarry Smith t[i] = omega*(sum/d); 97917ab2063SBarry Smith } 98017ab2063SBarry Smith 98117ab2063SBarry Smith /* x = x + t */ 98217ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 9830452661fSBarry Smith PetscFree(t); 98417ab2063SBarry Smith return 0; 98517ab2063SBarry Smith } 98617ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 98717ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 98817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 989416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 990416022c9SBarry Smith n = diag[i] - a->i[i]; 991416022c9SBarry Smith idx = a->j + a->i[i] + shift; 992416022c9SBarry Smith v = a->a + a->i[i] + shift; 99317ab2063SBarry Smith sum = b[i]; 99417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 99517ab2063SBarry Smith x[i] = omega*(sum/d); 99617ab2063SBarry Smith } 99717ab2063SBarry Smith xb = x; 99817ab2063SBarry Smith } 99917ab2063SBarry Smith else xb = b; 100017ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 100117ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 100217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1003416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 100417ab2063SBarry Smith } 100517ab2063SBarry Smith } 100617ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 100717ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1008416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1009416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1010416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1011416022c9SBarry Smith v = a->a + diag[i] + (!shift); 101217ab2063SBarry Smith sum = xb[i]; 101317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 101417ab2063SBarry Smith x[i] = omega*(sum/d); 101517ab2063SBarry Smith } 101617ab2063SBarry Smith } 101717ab2063SBarry Smith its--; 101817ab2063SBarry Smith } 101917ab2063SBarry Smith while (its--) { 102017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 102117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1022416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1023416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1024416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1025416022c9SBarry Smith v = a->a + a->i[i] + shift; 102617ab2063SBarry Smith sum = b[i]; 102717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 10287e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 102917ab2063SBarry Smith } 103017ab2063SBarry Smith } 103117ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 103217ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1033416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1034416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1035416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1036416022c9SBarry Smith v = a->a + a->i[i] + shift; 103717ab2063SBarry Smith sum = b[i]; 103817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 10397e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 104017ab2063SBarry Smith } 104117ab2063SBarry Smith } 104217ab2063SBarry Smith } 104317ab2063SBarry Smith return 0; 104417ab2063SBarry Smith } 104517ab2063SBarry Smith 10465615d1e5SSatish Balay #undef __FUNC__ 10475615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ" 10488f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 104917ab2063SBarry Smith { 1050416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 10514e220ebcSLois Curfman McInnes 10524e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 10534e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 10544e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 10554e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 10564e220ebcSLois Curfman McInnes info->block_size = 1.0; 10574e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 10584e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 10594e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 10604e220ebcSLois Curfman McInnes /* if (info->nz_unneeded != A->info.nz_unneeded) 10614e220ebcSLois Curfman McInnes printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 10624e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 10634e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 10644e220ebcSLois Curfman McInnes info->memory = A->mem; 10654e220ebcSLois Curfman McInnes if (A->factor) { 10664e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 10674e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 10684e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 10694e220ebcSLois Curfman McInnes } else { 10704e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 10714e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10724e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10734e220ebcSLois Curfman McInnes } 107417ab2063SBarry Smith return 0; 107517ab2063SBarry Smith } 107617ab2063SBarry Smith 107717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 107817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 107917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 108017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 108117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 108217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 108317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 108417ab2063SBarry Smith 10855615d1e5SSatish Balay #undef __FUNC__ 10865615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ" 10878f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 108817ab2063SBarry Smith { 1089416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1090416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 109117ab2063SBarry Smith 109277c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 109317ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 109417ab2063SBarry Smith if (diag) { 109517ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1096e3372554SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range"); 1097416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 1098416022c9SBarry Smith a->ilen[rows[i]] = 1; 1099416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 1100416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 110117ab2063SBarry Smith } 110217ab2063SBarry Smith else { 110317ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 110417ab2063SBarry Smith CHKERRQ(ierr); 110517ab2063SBarry Smith } 110617ab2063SBarry Smith } 110717ab2063SBarry Smith } 110817ab2063SBarry Smith else { 110917ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1110e3372554SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range"); 1111416022c9SBarry Smith a->ilen[rows[i]] = 0; 111217ab2063SBarry Smith } 111317ab2063SBarry Smith } 111417ab2063SBarry Smith ISRestoreIndices(is,&rows); 111543a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 111617ab2063SBarry Smith return 0; 111717ab2063SBarry Smith } 111817ab2063SBarry Smith 11195615d1e5SSatish Balay #undef __FUNC__ 11205615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ" 11218f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 112217ab2063SBarry Smith { 1123416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1124416022c9SBarry Smith *m = a->m; *n = a->n; 112517ab2063SBarry Smith return 0; 112617ab2063SBarry Smith } 112717ab2063SBarry Smith 11285615d1e5SSatish Balay #undef __FUNC__ 11295615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ" 11308f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 113117ab2063SBarry Smith { 1132416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1133416022c9SBarry Smith *m = 0; *n = a->m; 113417ab2063SBarry Smith return 0; 113517ab2063SBarry Smith } 1136026e39d0SSatish Balay 11375615d1e5SSatish Balay #undef __FUNC__ 11385615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ" 11394e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 114017ab2063SBarry Smith { 1141416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1142c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 114317ab2063SBarry Smith 1144e3372554SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,0,"Row out of range"); 114517ab2063SBarry Smith 1146416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 1147416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 114817ab2063SBarry Smith if (idx) { 1149416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 11504e093b46SBarry Smith if (*nz && shift) { 11510452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 115217ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 11534e093b46SBarry Smith } else if (*nz) { 11544e093b46SBarry Smith *idx = itmp; 115517ab2063SBarry Smith } 115617ab2063SBarry Smith else *idx = 0; 115717ab2063SBarry Smith } 115817ab2063SBarry Smith return 0; 115917ab2063SBarry Smith } 116017ab2063SBarry Smith 11615615d1e5SSatish Balay #undef __FUNC__ 11625615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ" 11634e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 116417ab2063SBarry Smith { 11654e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11664e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 116717ab2063SBarry Smith return 0; 116817ab2063SBarry Smith } 116917ab2063SBarry Smith 11705615d1e5SSatish Balay #undef __FUNC__ 11715615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ" 11728f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 117317ab2063SBarry Smith { 1174416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1175416022c9SBarry Smith Scalar *v = a->a; 117617ab2063SBarry Smith double sum = 0.0; 1177416022c9SBarry Smith int i, j,shift = a->indexshift; 117817ab2063SBarry Smith 117917ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1180416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 118117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 118217ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 118317ab2063SBarry Smith #else 118417ab2063SBarry Smith sum += (*v)*(*v); v++; 118517ab2063SBarry Smith #endif 118617ab2063SBarry Smith } 118717ab2063SBarry Smith *norm = sqrt(sum); 118817ab2063SBarry Smith } 118917ab2063SBarry Smith else if (type == NORM_1) { 119017ab2063SBarry Smith double *tmp; 1191416022c9SBarry Smith int *jj = a->j; 119266963ce1SSatish Balay tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp); 1193cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 119417ab2063SBarry Smith *norm = 0.0; 1195416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1196a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 119717ab2063SBarry Smith } 1198416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 119917ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 120017ab2063SBarry Smith } 12010452661fSBarry Smith PetscFree(tmp); 120217ab2063SBarry Smith } 120317ab2063SBarry Smith else if (type == NORM_INFINITY) { 120417ab2063SBarry Smith *norm = 0.0; 1205416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1206416022c9SBarry Smith v = a->a + a->i[j] + shift; 120717ab2063SBarry Smith sum = 0.0; 1208416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1209cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 121017ab2063SBarry Smith } 121117ab2063SBarry Smith if (sum > *norm) *norm = sum; 121217ab2063SBarry Smith } 121317ab2063SBarry Smith } 121417ab2063SBarry Smith else { 1215e3372554SBarry Smith SETERRQ(1,0,"No support for two norm yet"); 121617ab2063SBarry Smith } 121717ab2063SBarry Smith return 0; 121817ab2063SBarry Smith } 121917ab2063SBarry Smith 12205615d1e5SSatish Balay #undef __FUNC__ 12215615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ" 12228f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B) 122317ab2063SBarry Smith { 1224416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1225416022c9SBarry Smith Mat C; 1226416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1227416022c9SBarry Smith int shift = a->indexshift; 1228d5d45c9bSBarry Smith Scalar *array = a->a; 122917ab2063SBarry Smith 12303638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 1231e3372554SBarry Smith SETERRQ(1,0,"Square matrix only for in-place"); 12320452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1233cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 123417ab2063SBarry Smith if (shift) { 123517ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 123617ab2063SBarry Smith } 123717ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1238416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 12390452661fSBarry Smith PetscFree(col); 124017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 124117ab2063SBarry Smith len = ai[i+1]-ai[i]; 1242416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 124317ab2063SBarry Smith array += len; aj += len; 124417ab2063SBarry Smith } 124517ab2063SBarry Smith if (shift) { 124617ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 124717ab2063SBarry Smith } 124817ab2063SBarry Smith 12496d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 12506d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 125117ab2063SBarry Smith 12523638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1253416022c9SBarry Smith *B = C; 125417ab2063SBarry Smith } else { 1255416022c9SBarry Smith /* This isn't really an in-place transpose */ 12560452661fSBarry Smith PetscFree(a->a); 12570452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 12580452661fSBarry Smith if (a->diag) PetscFree(a->diag); 12590452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 12600452661fSBarry Smith if (a->imax) PetscFree(a->imax); 12610452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 12621073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 12630452661fSBarry Smith PetscFree(a); 1264f09e8eb9SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 12650452661fSBarry Smith PetscHeaderDestroy(C); 126617ab2063SBarry Smith } 126717ab2063SBarry Smith return 0; 126817ab2063SBarry Smith } 126917ab2063SBarry Smith 12705615d1e5SSatish Balay #undef __FUNC__ 12715615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ" 12728f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 127317ab2063SBarry Smith { 1274416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 127517ab2063SBarry Smith Scalar *l,*r,x,*v; 1276d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 127717ab2063SBarry Smith 127817ab2063SBarry Smith if (ll) { 12793ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 12803ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 12819b1297e1SSatish Balay VecGetLocalSize_Fast(ll,m); 1282e3372554SBarry Smith if (m != a->m) SETERRQ(1,0,"Left scaling vector wrong length"); 128390f02eecSBarry Smith VecGetArray_Fast(ll,l); 1284416022c9SBarry Smith v = a->a; 128517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 128617ab2063SBarry Smith x = l[i]; 1287416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 128817ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 128917ab2063SBarry Smith } 129044cd7ae7SLois Curfman McInnes PLogFlops(nz); 129117ab2063SBarry Smith } 129217ab2063SBarry Smith if (rr) { 12939b1297e1SSatish Balay VecGetLocalSize_Fast(rr,n); 1294e3372554SBarry Smith if (n != a->n) SETERRQ(1,0,"Right scaling vector wrong length"); 129590f02eecSBarry Smith VecGetArray_Fast(rr,r); 1296416022c9SBarry Smith v = a->a; jj = a->j; 129717ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 129817ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 129917ab2063SBarry Smith } 130044cd7ae7SLois Curfman McInnes PLogFlops(nz); 130117ab2063SBarry Smith } 130217ab2063SBarry Smith return 0; 130317ab2063SBarry Smith } 130417ab2063SBarry Smith 13055615d1e5SSatish Balay #undef __FUNC__ 13065615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ" 13078f6be9afSLois Curfman McInnes int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 130817ab2063SBarry Smith { 1309db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 131002834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 131199141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1312a2744918SBarry Smith register int sum,lensi; 131399141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 131499141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 131599141d43SSatish Balay Scalar *a_new,*mat_a; 1316416022c9SBarry Smith Mat C; 131717ab2063SBarry Smith 1318b48a1e75SSatish Balay ierr = ISSorted(isrow,(PetscTruth*)&i); 1319e3372554SBarry Smith if (!i) SETERRQ(1,0,"ISrow is not sorted"); 132099141d43SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i); 1321e3372554SBarry Smith if (!i) SETERRQ(1,0,"IScol is not sorted"); 132299141d43SSatish Balay 132317ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 132417ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 132517ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 132617ab2063SBarry Smith 13277264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 132802834360SBarry Smith /* special case of contiguous rows */ 132957faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 133002834360SBarry Smith starts = lens + ncols; 133102834360SBarry Smith /* loop over new rows determining lens and starting points */ 133202834360SBarry Smith for (i=0; i<nrows; i++) { 1333a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1334a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 133502834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1336d8ced48eSBarry Smith if (aj[k]+shift >= first) { 133702834360SBarry Smith starts[i] = k; 133802834360SBarry Smith break; 133902834360SBarry Smith } 134002834360SBarry Smith } 1341a2744918SBarry Smith sum = 0; 134202834360SBarry Smith while (k < kend) { 1343d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1344a2744918SBarry Smith sum++; 134502834360SBarry Smith } 1346a2744918SBarry Smith lens[i] = sum; 134702834360SBarry Smith } 134802834360SBarry Smith /* create submatrix */ 1349cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 135008480c60SBarry Smith int n_cols,n_rows; 135108480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1352e3372554SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,0,""); 1353d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 135408480c60SBarry Smith C = *B; 135508480c60SBarry Smith } 135608480c60SBarry Smith else { 135702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 135808480c60SBarry Smith } 1359db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1360db02288aSLois Curfman McInnes 136102834360SBarry Smith /* loop over rows inserting into submatrix */ 1362db02288aSLois Curfman McInnes a_new = c->a; 1363db02288aSLois Curfman McInnes j_new = c->j; 1364db02288aSLois Curfman McInnes i_new = c->i; 1365db02288aSLois Curfman McInnes i_new[0] = -shift; 136602834360SBarry Smith for (i=0; i<nrows; i++) { 1367a2744918SBarry Smith ii = starts[i]; 1368a2744918SBarry Smith lensi = lens[i]; 1369a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1370a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 137102834360SBarry Smith } 1372a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1373a2744918SBarry Smith a_new += lensi; 1374a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1375a2744918SBarry Smith c->ilen[i] = lensi; 137602834360SBarry Smith } 13770452661fSBarry Smith PetscFree(lens); 137802834360SBarry Smith } 137902834360SBarry Smith else { 138002834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 13810452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 138202834360SBarry Smith ssmap = smap + shift; 138399141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1384cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 138517ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 138602834360SBarry Smith /* determine lens of each row */ 138702834360SBarry Smith for (i=0; i<nrows; i++) { 1388d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 138902834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 139002834360SBarry Smith lens[i] = 0; 139102834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1392d8ced48eSBarry Smith if (ssmap[aj[k]]) { 139302834360SBarry Smith lens[i]++; 139402834360SBarry Smith } 139502834360SBarry Smith } 139602834360SBarry Smith } 139717ab2063SBarry Smith /* Create and fill new matrix */ 1398a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 139999141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 140099141d43SSatish Balay 1401e3372554SBarry Smith if (c->m != nrows || c->n != ncols) SETERRQ(1,0,""); 140299141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1403e3372554SBarry Smith SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros"); 140499141d43SSatish Balay } 140599141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 140608480c60SBarry Smith C = *B; 140708480c60SBarry Smith } 140808480c60SBarry Smith else { 140902834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 141008480c60SBarry Smith } 141199141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 141217ab2063SBarry Smith for (i=0; i<nrows; i++) { 141399141d43SSatish Balay row = irow[i]; 141417ab2063SBarry Smith nznew = 0; 141599141d43SSatish Balay kstart = ai[row]+shift; 141699141d43SSatish Balay kend = kstart + a->ilen[row]; 141799141d43SSatish Balay mat_i = c->i[i]+shift; 141899141d43SSatish Balay mat_j = c->j + mat_i; 141999141d43SSatish Balay mat_a = c->a + mat_i; 142099141d43SSatish Balay mat_ilen = c->ilen + i; 142117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 142299141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 142399141d43SSatish Balay *mat_j++ = tcol - (!shift); 142499141d43SSatish Balay *mat_a++ = a->a[k]; 142599141d43SSatish Balay (*mat_ilen)++; 142699141d43SSatish Balay 142717ab2063SBarry Smith } 142817ab2063SBarry Smith } 142917ab2063SBarry Smith } 143002834360SBarry Smith /* Free work space */ 143102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 143299141d43SSatish Balay PetscFree(smap); PetscFree(lens); 143302834360SBarry Smith } 14346d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 14356d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 143617ab2063SBarry Smith 143717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1438416022c9SBarry Smith *B = C; 143917ab2063SBarry Smith return 0; 144017ab2063SBarry Smith } 144117ab2063SBarry Smith 1442a871dcd8SBarry Smith /* 144363b91edcSBarry Smith note: This can only work for identity for row and col. It would 144463b91edcSBarry Smith be good to check this and otherwise generate an error. 1445a871dcd8SBarry Smith */ 14465615d1e5SSatish Balay #undef __FUNC__ 14475615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ" 14488f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1449a871dcd8SBarry Smith { 145063b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 145108480c60SBarry Smith int ierr; 145263b91edcSBarry Smith Mat outA; 145363b91edcSBarry Smith 1454e3372554SBarry Smith if (fill != 0) SETERRQ(1,0,"Only fill=0 supported"); 1455a871dcd8SBarry Smith 145663b91edcSBarry Smith outA = inA; 145763b91edcSBarry Smith inA->factor = FACTOR_LU; 145863b91edcSBarry Smith a->row = row; 145963b91edcSBarry Smith a->col = col; 146063b91edcSBarry Smith 146194a9d846SBarry Smith if (!a->solve_work) { /* this matrix may have been factored before */ 14620452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 146394a9d846SBarry Smith } 146463b91edcSBarry Smith 146508480c60SBarry Smith if (!a->diag) { 146608480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 146763b91edcSBarry Smith } 146863b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1469a871dcd8SBarry Smith return 0; 1470a871dcd8SBarry Smith } 1471a871dcd8SBarry Smith 1472f0b747eeSBarry Smith #include "pinclude/plapack.h" 14735615d1e5SSatish Balay #undef __FUNC__ 14745615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ" 14758f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1476f0b747eeSBarry Smith { 1477f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1478f0b747eeSBarry Smith int one = 1; 1479f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1480f0b747eeSBarry Smith PLogFlops(a->nz); 1481f0b747eeSBarry Smith return 0; 1482f0b747eeSBarry Smith } 1483f0b747eeSBarry Smith 14845615d1e5SSatish Balay #undef __FUNC__ 14855615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ" 14868f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1487cddf8d76SBarry Smith Mat **B) 1488cddf8d76SBarry Smith { 1489cddf8d76SBarry Smith int ierr,i; 1490cddf8d76SBarry Smith 1491cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 14920452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1493cddf8d76SBarry Smith } 1494cddf8d76SBarry Smith 1495cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1496905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1497cddf8d76SBarry Smith } 1498cddf8d76SBarry Smith return 0; 1499cddf8d76SBarry Smith } 1500cddf8d76SBarry Smith 15015615d1e5SSatish Balay #undef __FUNC__ 15025615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ" 15038f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 15045a838052SSatish Balay { 15055a838052SSatish Balay *bs = 1; 15065a838052SSatish Balay return 0; 15075a838052SSatish Balay } 15085a838052SSatish Balay 15095615d1e5SSatish Balay #undef __FUNC__ 15105615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ" 15118f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 15124dcbc457SBarry Smith { 1513e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 151406763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 15158a047759SSatish Balay int start, end, *ai, *aj; 151606763907SSatish Balay char *table; 15178a047759SSatish Balay shift = a->indexshift; 1518e4d965acSSatish Balay m = a->m; 1519e4d965acSSatish Balay ai = a->i; 15208a047759SSatish Balay aj = a->j+shift; 15218a047759SSatish Balay 1522e3372554SBarry Smith if (ov < 0) SETERRQ(1,0,"illegal overlap value used"); 152306763907SSatish Balay 152406763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 152506763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 152606763907SSatish Balay 1527e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1528b97fc60eSLois Curfman McInnes /* Initialize the two local arrays */ 1529e4d965acSSatish Balay isz = 0; 153006763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1531e4d965acSSatish Balay 1532e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 15334dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 153477c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1535e4d965acSSatish Balay 1536dd097bc3SLois Curfman McInnes /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 1537e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 153806763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 15394dcbc457SBarry Smith } 154006763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 154106763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1542e4d965acSSatish Balay 154304a348a9SBarry Smith k = 0; 154404a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap */ 154504a348a9SBarry Smith n = isz; 154606763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1547e4d965acSSatish Balay row = nidx[k]; 1548e4d965acSSatish Balay start = ai[row]; 1549e4d965acSSatish Balay end = ai[row+1]; 155004a348a9SBarry Smith for ( l = start; l<end ; l++){ 15518a047759SSatish Balay val = aj[l] + shift; 155206763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1553e4d965acSSatish Balay } 1554e4d965acSSatish Balay } 1555e4d965acSSatish Balay } 1556029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1557e4d965acSSatish Balay } 155804a348a9SBarry Smith PetscFree(table); 155906763907SSatish Balay PetscFree(nidx); 1560e4d965acSSatish Balay return 0; 15614dcbc457SBarry Smith } 156217ab2063SBarry Smith 15630513a670SBarry Smith /* -------------------------------------------------------------- */ 15640513a670SBarry Smith #undef __FUNC__ 15650513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ" 15660513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 15670513a670SBarry Smith { 15680513a670SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 15690513a670SBarry Smith Scalar *vwork; 15700513a670SBarry Smith int i, ierr, nz, m = a->m, n = a->n, *cwork; 15710513a670SBarry Smith int *row,*col,*cnew,j,*lens; 157256cd22aeSBarry Smith IS icolp,irowp; 15730513a670SBarry Smith 157456cd22aeSBarry Smith ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr); 157556cd22aeSBarry Smith ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr); 157656cd22aeSBarry Smith ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr); 157756cd22aeSBarry Smith ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr); 15780513a670SBarry Smith 15790513a670SBarry Smith /* determine lengths of permuted rows */ 15800513a670SBarry Smith lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens); 15810513a670SBarry Smith for (i=0; i<m; i++ ) { 15820513a670SBarry Smith lens[row[i]] = a->i[i+1] - a->i[i]; 15830513a670SBarry Smith } 15840513a670SBarry Smith ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr); 15850513a670SBarry Smith PetscFree(lens); 15860513a670SBarry Smith 15870513a670SBarry Smith cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew); 15880513a670SBarry Smith for (i=0; i<m; i++) { 15890513a670SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 15900513a670SBarry Smith for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];} 15910513a670SBarry Smith ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr); 15920513a670SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 15930513a670SBarry Smith } 15940513a670SBarry Smith PetscFree(cnew); 15950513a670SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 15960513a670SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 159756cd22aeSBarry Smith ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr); 159856cd22aeSBarry Smith ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr); 159956cd22aeSBarry Smith ierr = ISDestroy(irowp); CHKERRQ(ierr); 160056cd22aeSBarry Smith ierr = ISDestroy(icolp); CHKERRQ(ierr); 16010513a670SBarry Smith return 0; 16020513a670SBarry Smith } 16030513a670SBarry Smith 16045615d1e5SSatish Balay #undef __FUNC__ 16055615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ" 1606682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1607682d7d0cSBarry Smith { 1608682d7d0cSBarry Smith static int called = 0; 1609682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1610682d7d0cSBarry Smith 1611682d7d0cSBarry Smith if (called) return 0; else called = 1; 161277c4ece6SBarry Smith PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 16130f665d81SLois Curfman McInnes PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n"); 16140f665d81SLois Curfman McInnes PetscPrintf(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n"); 16150f665d81SLois Curfman McInnes PetscPrintf(comm," -mat_aij_no_inode: Do not use inodes\n"); 16160f665d81SLois Curfman McInnes PetscPrintf(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n"); 1617682d7d0cSBarry Smith #if defined(HAVE_ESSL) 16180f665d81SLois Curfman McInnes PetscPrintf(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.\n"); 1619682d7d0cSBarry Smith #endif 1620682d7d0cSBarry Smith return 0; 1621682d7d0cSBarry Smith } 16228f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1623a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1624a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1625a93ec695SBarry Smith 1626682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 162717ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 162817ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1629416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1630416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 163117ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 163217ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 163317ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 163417ab2063SBarry Smith MatRelax_SeqAIJ, 163517ab2063SBarry Smith MatTranspose_SeqAIJ, 16367264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1637f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 163817ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 163917ab2063SBarry Smith MatCompress_SeqAIJ, 164017ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 164117ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 164217ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 164317ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 164494a9d846SBarry Smith 0,0, 16453d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1646cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 16477eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1648682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1649f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 16505a838052SSatish Balay MatScale_SeqAIJ,0,0, 16516945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 16526945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 16533b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 16543b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 16553b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1656a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1657a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 16580513a670SBarry Smith MatColoringPatch_SeqAIJ, 16590513a670SBarry Smith 0, 16600513a670SBarry Smith MatPermute_SeqAIJ}; 166117ab2063SBarry Smith 166217ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 166317ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 166417ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 166517ab2063SBarry Smith 16665615d1e5SSatish Balay #undef __FUNC__ 16675615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ" 166817ab2063SBarry Smith /*@C 1669682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 16700d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 16716e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 16722bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 16732bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 167417ab2063SBarry Smith 167517ab2063SBarry Smith Input Parameters: 1676029af93fSBarry Smith . comm - MPI communicator, set to PETSC_COMM_SELF 167717ab2063SBarry Smith . m - number of rows 167817ab2063SBarry Smith . n - number of columns 167917ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 16802bd5e0b2SLois Curfman McInnes . nzz - array containing the number of nonzeros in the various rows 16812bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 168217ab2063SBarry Smith 168317ab2063SBarry Smith Output Parameter: 1684416022c9SBarry Smith . A - the matrix 168517ab2063SBarry Smith 168617ab2063SBarry Smith Notes: 168717ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 168817ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 16890002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 169044cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 169117ab2063SBarry Smith 169217ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1693a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 16943d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 16956da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 169617ab2063SBarry Smith 1697682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 16984fca80b9SLois Curfman McInnes improve numerical efficiency of matrix-vector products and solves. We 1699682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 17006c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 17016c7ebb05SLois Curfman McInnes 17026c7ebb05SLois Curfman McInnes Options Database Keys: 17036c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 17046e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 17056e62573dSLois Curfman McInnes $ (max limit=5) 17066e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 17076e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 17086e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 170917ab2063SBarry Smith 171017ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 171117ab2063SBarry Smith @*/ 1712416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 171317ab2063SBarry Smith { 1714416022c9SBarry Smith Mat B; 1715416022c9SBarry Smith Mat_SeqAIJ *b; 17166945ee14SBarry Smith int i, len, ierr, flg,size; 17176945ee14SBarry Smith 17186945ee14SBarry Smith MPI_Comm_size(comm,&size); 1719e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"Comm must be of size 1"); 1720d5d45c9bSBarry Smith 1721416022c9SBarry Smith *A = 0; 1722d4bb536fSBarry Smith PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView); 1723416022c9SBarry Smith PLogObjectCreate(B); 17240452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 172544cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1726416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1727416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1728416022c9SBarry Smith B->view = MatView_SeqAIJ; 1729416022c9SBarry Smith B->factor = 0; 1730416022c9SBarry Smith B->lupivotthreshold = 1.0; 173190f02eecSBarry Smith B->mapping = 0; 17327a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 173369957df2SSatish Balay &flg); CHKERRQ(ierr); 17347a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 17357a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 17367a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1737416022c9SBarry Smith b->row = 0; 1738416022c9SBarry Smith b->col = 0; 1739416022c9SBarry Smith b->indexshift = 0; 1740b810aeb4SBarry Smith b->reallocs = 0; 174169957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 174269957df2SSatish Balay if (flg) b->indexshift = -1; 174317ab2063SBarry Smith 174444cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 174544cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 17460452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1747b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 17487b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 17497b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1750416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 175117ab2063SBarry Smith nz = nz*m; 175217ab2063SBarry Smith } 175317ab2063SBarry Smith else { 175417ab2063SBarry Smith nz = 0; 1755416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 175617ab2063SBarry Smith } 175717ab2063SBarry Smith 175817ab2063SBarry Smith /* allocate the matrix space */ 1759416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 17600452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1761416022c9SBarry Smith b->j = (int *) (b->a + nz); 1762cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1763416022c9SBarry Smith b->i = b->j + nz; 1764416022c9SBarry Smith b->singlemalloc = 1; 176517ab2063SBarry Smith 1766416022c9SBarry Smith b->i[0] = -b->indexshift; 176717ab2063SBarry Smith for (i=1; i<m+1; i++) { 1768416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 176917ab2063SBarry Smith } 177017ab2063SBarry Smith 1771416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 17720452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1773f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 1774416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 177517ab2063SBarry Smith 1776416022c9SBarry Smith b->nz = 0; 1777416022c9SBarry Smith b->maxnz = nz; 1778416022c9SBarry Smith b->sorted = 0; 1779416022c9SBarry Smith b->roworiented = 1; 1780416022c9SBarry Smith b->nonew = 0; 1781416022c9SBarry Smith b->diag = 0; 1782416022c9SBarry Smith b->solve_work = 0; 178371bd300dSLois Curfman McInnes b->spptr = 0; 1784754ec7b1SSatish Balay b->inode.node_count = 0; 1785754ec7b1SSatish Balay b->inode.size = 0; 17866c7ebb05SLois Curfman McInnes b->inode.limit = 5; 17876c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 17884e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 178917ab2063SBarry Smith 1790416022c9SBarry Smith *A = B; 17914e220ebcSLois Curfman McInnes 17924b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 17934b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 179469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 179569957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 17964b14c69eSBarry Smith #endif 179769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 179869957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 179969957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 180069957df2SSatish Balay if (flg) { 1801e3372554SBarry Smith if (!b->indexshift) SETERRQ(1,0,"need -mat_aij_oneindex with -mat_aij_dxml"); 1802416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 180317ab2063SBarry Smith } 180469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 180569957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 180617ab2063SBarry Smith return 0; 180717ab2063SBarry Smith } 180817ab2063SBarry Smith 18095615d1e5SSatish Balay #undef __FUNC__ 18105615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ" 18113d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 181217ab2063SBarry Smith { 1813416022c9SBarry Smith Mat C; 1814416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 181508480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 181617ab2063SBarry Smith 18174043dd9cSLois Curfman McInnes *B = 0; 1818d4bb536fSBarry Smith PetscHeaderCreate(C,_p_Mat,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView); 1819416022c9SBarry Smith PLogObjectCreate(C); 18200452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 182141c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1822416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1823416022c9SBarry Smith C->view = MatView_SeqAIJ; 1824416022c9SBarry Smith C->factor = A->factor; 1825416022c9SBarry Smith c->row = 0; 1826416022c9SBarry Smith c->col = 0; 1827416022c9SBarry Smith c->indexshift = shift; 1828c456f294SBarry Smith C->assembled = PETSC_TRUE; 182917ab2063SBarry Smith 183044cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 183144cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 183244cd7ae7SLois Curfman McInnes C->M = a->m; 183344cd7ae7SLois Curfman McInnes C->N = a->n; 183417ab2063SBarry Smith 18350452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 18360452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 183717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1838416022c9SBarry Smith c->imax[i] = a->imax[i]; 1839416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 184017ab2063SBarry Smith } 184117ab2063SBarry Smith 184217ab2063SBarry Smith /* allocate the matrix space */ 1843416022c9SBarry Smith c->singlemalloc = 1; 1844416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 18450452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1846416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1847416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1848416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 184917ab2063SBarry Smith if (m > 0) { 1850416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 185108480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1852416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 185317ab2063SBarry Smith } 185408480c60SBarry Smith } 185517ab2063SBarry Smith 1856f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 1857416022c9SBarry Smith c->sorted = a->sorted; 1858416022c9SBarry Smith c->roworiented = a->roworiented; 1859416022c9SBarry Smith c->nonew = a->nonew; 18607a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 186117ab2063SBarry Smith 1862416022c9SBarry Smith if (a->diag) { 18630452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1864416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 186517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1866416022c9SBarry Smith c->diag[i] = a->diag[i]; 186717ab2063SBarry Smith } 186817ab2063SBarry Smith } 1869416022c9SBarry Smith else c->diag = 0; 18706c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 18716c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1872754ec7b1SSatish Balay if (a->inode.size){ 1873daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 1874754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1875daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 1876754ec7b1SSatish Balay } else { 1877754ec7b1SSatish Balay c->inode.size = 0; 1878754ec7b1SSatish Balay c->inode.node_count = 0; 1879754ec7b1SSatish Balay } 1880416022c9SBarry Smith c->nz = a->nz; 1881416022c9SBarry Smith c->maxnz = a->maxnz; 1882416022c9SBarry Smith c->solve_work = 0; 188376dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1884754ec7b1SSatish Balay 1885416022c9SBarry Smith *B = C; 188617ab2063SBarry Smith return 0; 188717ab2063SBarry Smith } 188817ab2063SBarry Smith 18895615d1e5SSatish Balay #undef __FUNC__ 18905615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ" 189119bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 189217ab2063SBarry Smith { 1893416022c9SBarry Smith Mat_SeqAIJ *a; 1894416022c9SBarry Smith Mat B; 189517699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1896bcd2baecSBarry Smith MPI_Comm comm; 189717ab2063SBarry Smith 189819bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 189917699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 1900e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"view must have one processor"); 190190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 190277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 1903e3372554SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object in file"); 190417ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 190517ab2063SBarry Smith 190617ab2063SBarry Smith /* read in row lengths */ 19070452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 190877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 190917ab2063SBarry Smith 191017ab2063SBarry Smith /* create our matrix */ 1911416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1912416022c9SBarry Smith B = *A; 1913416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1914416022c9SBarry Smith shift = a->indexshift; 191517ab2063SBarry Smith 191617ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 191777c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 191817ab2063SBarry Smith if (shift) { 191917ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1920416022c9SBarry Smith a->j[i] += 1; 192117ab2063SBarry Smith } 192217ab2063SBarry Smith } 192317ab2063SBarry Smith 192417ab2063SBarry Smith /* read in nonzero values */ 192577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 192617ab2063SBarry Smith 192717ab2063SBarry Smith /* set matrix "i" values */ 1928416022c9SBarry Smith a->i[0] = -shift; 192917ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1930416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1931416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 193217ab2063SBarry Smith } 19330452661fSBarry Smith PetscFree(rowlengths); 193417ab2063SBarry Smith 19356d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 19366d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 193717ab2063SBarry Smith return 0; 193817ab2063SBarry Smith } 193917ab2063SBarry Smith 19405615d1e5SSatish Balay #undef __FUNC__ 19415615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ" 19428f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 19437264ac53SSatish Balay { 19447264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 19457264ac53SSatish Balay 1946e3372554SBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,0,"Matrices must be same type"); 19477264ac53SSatish Balay 19487264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 19497264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1950bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 195177c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1952bcd2baecSBarry Smith } 19537264ac53SSatish Balay 19547264ac53SSatish Balay /* if the a->i are the same */ 19558108c231SLois Curfman McInnes if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) { 195677c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 19577264ac53SSatish Balay } 19587264ac53SSatish Balay 19597264ac53SSatish Balay /* if a->j are the same */ 1960bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 196177c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1962bcd2baecSBarry Smith } 1963bcd2baecSBarry Smith 1964bcd2baecSBarry Smith /* if a->a are the same */ 196519bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 196677c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 19677264ac53SSatish Balay } 196877c4ece6SBarry Smith *flg = PETSC_TRUE; 19697264ac53SSatish Balay return 0; 19707264ac53SSatish Balay 19717264ac53SSatish Balay } 1972