16d84be18SBarry Smith 26945ee14SBarry Smith 317ab2063SBarry Smith #ifndef lint 4*7e33a6baSBarry Smith static char vcid[] = "$Id: aij.c,v 1.186 1996/09/20 03:56:04 balay Exp bsmith $"; 517ab2063SBarry Smith #endif 617ab2063SBarry Smith 7d5d45c9bSBarry Smith /* 85a838052SSatish Balay B Defines the basic matrix operations for the AIJ (compressed row) 9d5d45c9bSBarry Smith matrix storage format. 10d5d45c9bSBarry Smith */ 1170f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 12f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 13f5eb4b81SSatish Balay #include "src/inline/spops.h" 14e4d965acSSatish Balay #include "petsc.h" 15f5eb4b81SSatish Balay #include "src/inline/bitarray.h" 1617ab2063SBarry Smith 17a2ce50c7SBarry Smith /* 18a2ce50c7SBarry Smith Basic AIJ format ILU based on drop tolerance 19a2ce50c7SBarry Smith */ 20a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact) 21a2ce50c7SBarry Smith { 22a2ce50c7SBarry Smith /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */ 23a2ce50c7SBarry Smith int ierr = 1; 24a2ce50c7SBarry Smith 25a2ce50c7SBarry Smith SETERRQ(ierr,"MatILUDTFactor_SeqAIJ:Not implemented"); 26a2ce50c7SBarry Smith } 27a2ce50c7SBarry Smith 28bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 2917ab2063SBarry Smith 303b2fbd54SBarry Smith static int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 316945ee14SBarry Smith PetscTruth *done) 3217ab2063SBarry Smith { 33416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 346945ee14SBarry Smith int ierr,i,ishift; 3517ab2063SBarry Smith 366945ee14SBarry Smith *n = A->n; 376945ee14SBarry Smith if (!ia) return 0; 386945ee14SBarry Smith ishift = a->indexshift; 396945ee14SBarry Smith if (symmetric) { 406945ee14SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 416945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 426945ee14SBarry Smith int nz = a->i[a->n]; 433b2fbd54SBarry Smith /* malloc space and subtract 1 from i and j indices */ 443b2fbd54SBarry Smith *ia = (int *) PetscMalloc( (a->n+1)*sizeof(int) ); CHKPTRQ(*ia); 453b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 463b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1; 473b2fbd54SBarry Smith for ( i=0; i<a->n+1; i++ ) (*ia)[i] = a->i[i] - 1; 486945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 496945ee14SBarry Smith int nz = a->i[a->n] + 1; 503b2fbd54SBarry Smith /* malloc space and add 1 to i and j indices */ 513b2fbd54SBarry Smith *ia = (int *) PetscMalloc( (a->n+1)*sizeof(int) ); CHKPTRQ(*ia); 523b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 533b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1; 543b2fbd54SBarry Smith for ( i=0; i<a->n+1; i++ ) (*ia)[i] = a->i[i] + 1; 556945ee14SBarry Smith } else { 566945ee14SBarry Smith *ia = a->i; *ja = a->j; 57a2ce50c7SBarry Smith } 58a2ce50c7SBarry Smith 59a2744918SBarry Smith return 0; 60a2744918SBarry Smith } 61a2744918SBarry Smith 623b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 636945ee14SBarry Smith PetscTruth *done) 646945ee14SBarry Smith { 656945ee14SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 663b2fbd54SBarry Smith int ishift = a->indexshift; 676945ee14SBarry Smith 686945ee14SBarry Smith if (!ia) return 0; 693b2fbd54SBarry Smith if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) { 706945ee14SBarry Smith PetscFree(*ia); 716945ee14SBarry Smith PetscFree(*ja); 72bcd2baecSBarry Smith } 7317ab2063SBarry Smith return 0; 7417ab2063SBarry Smith } 7517ab2063SBarry Smith 763b2fbd54SBarry Smith static int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 773b2fbd54SBarry Smith PetscTruth *done) 783b2fbd54SBarry Smith { 793b2fbd54SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 803b2fbd54SBarry Smith int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n; 813b2fbd54SBarry Smith int nz = a->i[n]+ishift,row,*jj,m,col; 823b2fbd54SBarry Smith 833b2fbd54SBarry Smith *nn = A->n; 843b2fbd54SBarry Smith if (!ia) return 0; 853b2fbd54SBarry Smith if (symmetric) { 863b2fbd54SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 873b2fbd54SBarry Smith } else { 883b2fbd54SBarry Smith collengths = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(collengths); 893b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 903b2fbd54SBarry Smith cia = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia); 913b2fbd54SBarry Smith cja = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cja); 923b2fbd54SBarry Smith jj = a->j; 933b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) { 943b2fbd54SBarry Smith collengths[jj[i] + ishift]++; 953b2fbd54SBarry Smith } 963b2fbd54SBarry Smith cia[0] = oshift; 973b2fbd54SBarry Smith for ( i=0; i<n; i++) { 983b2fbd54SBarry Smith cia[i+1] = cia[i] + collengths[i]; 993b2fbd54SBarry Smith } 1003b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1013b2fbd54SBarry Smith jj = a->j; 1023b2fbd54SBarry Smith for ( row=0; row<n; row++ ) { 1033b2fbd54SBarry Smith m = a->i[row+1] - a->i[row]; 1043b2fbd54SBarry Smith for ( i=0; i<m; i++ ) { 1053b2fbd54SBarry Smith 1063b2fbd54SBarry Smith col = *jj++ + ishift; 1073b2fbd54SBarry Smith printf("row col %d %d\n",row,col); 1083b2fbd54SBarry Smith cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1093b2fbd54SBarry Smith } 1103b2fbd54SBarry Smith } 1113b2fbd54SBarry Smith 1123b2fbd54SBarry Smith jj = cja; 1133b2fbd54SBarry Smith for ( i=0; i<n; i++ ) { 1143b2fbd54SBarry Smith m = cia[i+1] - cia[i]; 1153b2fbd54SBarry Smith printf("Column %d length %d\n",i,m); 1163b2fbd54SBarry Smith for ( row = 0; row < m; row++ ) { 1173b2fbd54SBarry Smith printf("%d ",*jj++); 1183b2fbd54SBarry Smith } 1193b2fbd54SBarry Smith printf("\n"); 1203b2fbd54SBarry Smith } 1213b2fbd54SBarry Smith PetscFree(collengths); 1223b2fbd54SBarry Smith *ia = cia; *ja = cja; 1233b2fbd54SBarry Smith } 1243b2fbd54SBarry Smith 1253b2fbd54SBarry Smith return 0; 1263b2fbd54SBarry Smith } 1273b2fbd54SBarry Smith 1283b2fbd54SBarry Smith static 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 14117ab2063SBarry Smith /* This version has row oriented v */ 142416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 14317ab2063SBarry Smith { 144416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 145416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 1464b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 147d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 148416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 14917ab2063SBarry Smith 15017ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 151416022c9SBarry Smith row = im[k]; 1523b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 15317ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 154416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 1553b2fbd54SBarry Smith #endif 15617ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 15717ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 158416022c9SBarry Smith low = 0; 15917ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 1603b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 161416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 162416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 1633b2fbd54SBarry Smith #endif 1644b0e389bSBarry Smith col = in[l] - shift; 1654b0e389bSBarry Smith if (roworiented) { 1664b0e389bSBarry Smith value = *v++; 1674b0e389bSBarry Smith } 1684b0e389bSBarry Smith else { 1694b0e389bSBarry Smith value = v[k + l*m]; 1704b0e389bSBarry Smith } 171416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 172416022c9SBarry Smith while (high-low > 5) { 173416022c9SBarry Smith t = (low+high)/2; 174416022c9SBarry Smith if (rp[t] > col) high = t; 175416022c9SBarry Smith else low = t; 17617ab2063SBarry Smith } 177416022c9SBarry Smith for ( i=low; i<high; i++ ) { 17817ab2063SBarry Smith if (rp[i] > col) break; 17917ab2063SBarry Smith if (rp[i] == col) { 180416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 18117ab2063SBarry Smith else ap[i] = value; 18217ab2063SBarry Smith goto noinsert; 18317ab2063SBarry Smith } 18417ab2063SBarry Smith } 18517ab2063SBarry Smith if (nonew) goto noinsert; 18617ab2063SBarry Smith if (nrow >= rmax) { 18717ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 188416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 18917ab2063SBarry Smith Scalar *new_a; 19017ab2063SBarry Smith 19117ab2063SBarry Smith /* malloc new storage space */ 192416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1930452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 19417ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 19517ab2063SBarry Smith new_i = new_j + new_nz; 19617ab2063SBarry Smith 19717ab2063SBarry Smith /* copy over old data into new slots */ 19817ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 199416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 200416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 201416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 202416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 20317ab2063SBarry Smith len*sizeof(int)); 204416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 205416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 20617ab2063SBarry Smith len*sizeof(Scalar)); 20717ab2063SBarry Smith /* free up old matrix storage */ 2080452661fSBarry Smith PetscFree(a->a); 2090452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 210416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 211416022c9SBarry Smith a->singlemalloc = 1; 21217ab2063SBarry Smith 21317ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 214416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 215416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 216416022c9SBarry Smith a->maxnz += CHUNKSIZE; 217b810aeb4SBarry Smith a->reallocs++; 21817ab2063SBarry Smith } 219416022c9SBarry Smith N = nrow++ - 1; a->nz++; 220416022c9SBarry Smith /* shift up all the later entries in this row */ 221416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 22217ab2063SBarry Smith rp[ii+1] = rp[ii]; 22317ab2063SBarry Smith ap[ii+1] = ap[ii]; 22417ab2063SBarry Smith } 22517ab2063SBarry Smith rp[i] = col; 22617ab2063SBarry Smith ap[i] = value; 22717ab2063SBarry Smith noinsert:; 228416022c9SBarry Smith low = i + 1; 22917ab2063SBarry Smith } 23017ab2063SBarry Smith ailen[row] = nrow; 23117ab2063SBarry Smith } 23217ab2063SBarry Smith return 0; 23317ab2063SBarry Smith } 23417ab2063SBarry Smith 2357eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 2367eb43aa7SLois Curfman McInnes { 2377eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 238b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 2397eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 2407eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 2417eb43aa7SLois Curfman McInnes 2427eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 2437eb43aa7SLois Curfman McInnes row = im[k]; 2447eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 2457eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 2467eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 2477eb43aa7SLois Curfman McInnes nrow = ailen[row]; 2487eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 2497eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 2507eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 2517eb43aa7SLois Curfman McInnes col = in[l] - shift; 2527eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 2537eb43aa7SLois Curfman McInnes while (high-low > 5) { 2547eb43aa7SLois Curfman McInnes t = (low+high)/2; 2557eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 2567eb43aa7SLois Curfman McInnes else low = t; 2577eb43aa7SLois Curfman McInnes } 2587eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 2597eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2607eb43aa7SLois Curfman McInnes if (rp[i] == col) { 261b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2627eb43aa7SLois Curfman McInnes goto finished; 2637eb43aa7SLois Curfman McInnes } 2647eb43aa7SLois Curfman McInnes } 265b49de8d1SLois Curfman McInnes *v++ = zero; 2667eb43aa7SLois Curfman McInnes finished:; 2677eb43aa7SLois Curfman McInnes } 2687eb43aa7SLois Curfman McInnes } 2697eb43aa7SLois Curfman McInnes return 0; 2707eb43aa7SLois Curfman McInnes } 2717eb43aa7SLois Curfman McInnes 27217ab2063SBarry Smith #include "draw.h" 27317ab2063SBarry Smith #include "pinclude/pviewer.h" 27477c4ece6SBarry Smith #include "sys.h" 27517ab2063SBarry Smith 276416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 27717ab2063SBarry Smith { 278416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 279416022c9SBarry Smith int i, fd, *col_lens, ierr; 28017ab2063SBarry Smith 28190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2820452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 283416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 284416022c9SBarry Smith col_lens[1] = a->m; 285416022c9SBarry Smith col_lens[2] = a->n; 286416022c9SBarry Smith col_lens[3] = a->nz; 287416022c9SBarry Smith 288416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 289416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 290416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 29117ab2063SBarry Smith } 29277c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 2930452661fSBarry Smith PetscFree(col_lens); 294416022c9SBarry Smith 295416022c9SBarry Smith /* store column indices (zero start index) */ 296416022c9SBarry Smith if (a->indexshift) { 297416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 29817ab2063SBarry Smith } 29977c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 300416022c9SBarry Smith if (a->indexshift) { 301416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 30217ab2063SBarry Smith } 303416022c9SBarry Smith 304416022c9SBarry Smith /* store nonzero values */ 30577c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 30617ab2063SBarry Smith return 0; 30717ab2063SBarry Smith } 308416022c9SBarry Smith 309416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 310416022c9SBarry Smith { 311416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 312496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 31317ab2063SBarry Smith FILE *fd; 31417ab2063SBarry Smith char *outputname; 31517ab2063SBarry Smith 31690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 317416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 31890ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 31990ace30eSBarry Smith if (format == ASCII_FORMAT_INFO) { 32095e01e2fSLois Curfman McInnes return 0; 32195e01e2fSLois Curfman McInnes } 32290ace30eSBarry Smith else if (format == ASCII_FORMAT_INFO_DETAILED) { 323496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 324496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 325496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 32695e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 32795e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 32817ab2063SBarry Smith } 32990ace30eSBarry Smith else if (format == ASCII_FORMAT_MATLAB) { 330416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3314e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 3324e220ebcSLois Curfman McInnes fprintf(fd,"zzz = zeros(%d,3);\n",a->nz); 33317ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 33417ab2063SBarry Smith 33517ab2063SBarry Smith for (i=0; i<m; i++) { 336416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 33717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 3386945ee14SBarry Smith fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]), 339416022c9SBarry Smith imag(a->a[j])); 34017ab2063SBarry Smith #else 3417a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 34217ab2063SBarry Smith #endif 34317ab2063SBarry Smith } 34417ab2063SBarry Smith } 34517ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 34617ab2063SBarry Smith } 34744cd7ae7SLois Curfman McInnes else if (format == ASCII_FORMAT_COMMON) { 34844cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 34944cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 35044cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 35144cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 35244cd7ae7SLois Curfman McInnes if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0) 35344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 35444cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 35544cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 35644cd7ae7SLois Curfman McInnes #else 35744cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 35844cd7ae7SLois Curfman McInnes #endif 35944cd7ae7SLois Curfman McInnes } 36044cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 36144cd7ae7SLois Curfman McInnes } 36244cd7ae7SLois Curfman McInnes } 36317ab2063SBarry Smith else { 36417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 36517ab2063SBarry Smith fprintf(fd,"row %d:",i); 366416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 36717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 368416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 369416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 37017ab2063SBarry Smith } 37117ab2063SBarry Smith else { 372416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 37317ab2063SBarry Smith } 37417ab2063SBarry Smith #else 375416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 37617ab2063SBarry Smith #endif 37717ab2063SBarry Smith } 37817ab2063SBarry Smith fprintf(fd,"\n"); 37917ab2063SBarry Smith } 38017ab2063SBarry Smith } 38117ab2063SBarry Smith fflush(fd); 382416022c9SBarry Smith return 0; 383416022c9SBarry Smith } 384416022c9SBarry Smith 385416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 386416022c9SBarry Smith { 387416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 388cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 389cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 390bcd2baecSBarry Smith Draw draw; 391cddf8d76SBarry Smith DrawButton button; 39219bcc07fSBarry Smith PetscTruth isnull; 393cddf8d76SBarry Smith 394bcd2baecSBarry Smith ViewerDrawGetDraw(viewer,&draw); 39519bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 39619bcc07fSBarry Smith 397416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 398416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 399416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 400416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 401cddf8d76SBarry Smith color = DRAW_BLUE; 402416022c9SBarry Smith for ( i=0; i<m; i++ ) { 403cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 404416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 405cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 406cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 407cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 408cddf8d76SBarry Smith #else 409cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 410cddf8d76SBarry Smith #endif 411cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 412cddf8d76SBarry Smith } 413cddf8d76SBarry Smith } 414cddf8d76SBarry Smith color = DRAW_CYAN; 415cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 416cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 417cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 418cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 419cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 420cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 421cddf8d76SBarry Smith } 422cddf8d76SBarry Smith } 423cddf8d76SBarry Smith color = DRAW_RED; 424cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 425cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 426cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 427cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 428cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 429cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 430cddf8d76SBarry Smith #else 431cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 432cddf8d76SBarry Smith #endif 433cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 434416022c9SBarry Smith } 435416022c9SBarry Smith } 436416022c9SBarry Smith DrawFlush(draw); 437cddf8d76SBarry Smith DrawGetPause(draw,&pause); 438cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 439cddf8d76SBarry Smith 440cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 4416945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 442cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 443cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 444cddf8d76SBarry Smith DrawClear(draw); 445cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 446cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 447cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 448cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 449cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 450cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 451cddf8d76SBarry Smith w *= scale; h *= scale; 452cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 453cddf8d76SBarry Smith color = DRAW_BLUE; 454cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 455cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 456cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 457cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 458cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 459cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 460cddf8d76SBarry Smith #else 461cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 462cddf8d76SBarry Smith #endif 463cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 464cddf8d76SBarry Smith } 465cddf8d76SBarry Smith } 466cddf8d76SBarry Smith color = DRAW_CYAN; 467cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 468cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 469cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 470cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 471cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 472cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 473cddf8d76SBarry Smith } 474cddf8d76SBarry Smith } 475cddf8d76SBarry Smith color = DRAW_RED; 476cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 477cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 478cddf8d76SBarry 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 } 4886945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 489cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 490cddf8d76SBarry Smith } 491416022c9SBarry Smith return 0; 492416022c9SBarry Smith } 493416022c9SBarry Smith 494416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 495416022c9SBarry Smith { 496416022c9SBarry Smith Mat A = (Mat) obj; 497416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 498bcd2baecSBarry Smith ViewerType vtype; 499bcd2baecSBarry Smith int ierr; 500416022c9SBarry Smith 501bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 502bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 503416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 504416022c9SBarry Smith } 505bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 506416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 507416022c9SBarry Smith } 508bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 509416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 510416022c9SBarry Smith } 511bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 512bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 51317ab2063SBarry Smith } 51417ab2063SBarry Smith return 0; 51517ab2063SBarry Smith } 51619bcc07fSBarry Smith 517c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 518416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 51917ab2063SBarry Smith { 520416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 52141c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 522416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 523416022c9SBarry Smith Scalar *aa = a->a, *ap; 52417ab2063SBarry Smith 5256d4a8577SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) return 0; 52617ab2063SBarry Smith 52717ab2063SBarry Smith for ( i=1; i<m; i++ ) { 528416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 52917ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 53017ab2063SBarry Smith if (fshift) { 531416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 53217ab2063SBarry Smith N = ailen[i]; 53317ab2063SBarry Smith for ( j=0; j<N; j++ ) { 53417ab2063SBarry Smith ip[j-fshift] = ip[j]; 53517ab2063SBarry Smith ap[j-fshift] = ap[j]; 53617ab2063SBarry Smith } 53717ab2063SBarry Smith } 53817ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 53917ab2063SBarry Smith } 54017ab2063SBarry Smith if (m) { 54117ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 54217ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 54317ab2063SBarry Smith } 54417ab2063SBarry Smith /* reset ilen and imax for each row */ 54517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 54617ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 54717ab2063SBarry Smith } 548416022c9SBarry Smith a->nz = ai[m] + shift; 54917ab2063SBarry Smith 55017ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 551416022c9SBarry Smith if (fshift && a->diag) { 5520452661fSBarry Smith PetscFree(a->diag); 553416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 554416022c9SBarry Smith a->diag = 0; 55517ab2063SBarry Smith } 5564e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 5574e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 5584e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 559b810aeb4SBarry Smith a->reallocs); 5604e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 5614e220ebcSLois Curfman McInnes 56276dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 56341c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 56417ab2063SBarry Smith return 0; 56517ab2063SBarry Smith } 56617ab2063SBarry Smith 567416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 56817ab2063SBarry Smith { 569416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 570cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 57117ab2063SBarry Smith return 0; 57217ab2063SBarry Smith } 573416022c9SBarry Smith 57417ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 57517ab2063SBarry Smith { 576416022c9SBarry Smith Mat A = (Mat) obj; 577416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 578d5d45c9bSBarry Smith 57917ab2063SBarry Smith #if defined(PETSC_LOG) 580416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 58117ab2063SBarry Smith #endif 5820452661fSBarry Smith PetscFree(a->a); 5830452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 5840452661fSBarry Smith if (a->diag) PetscFree(a->diag); 5850452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 5860452661fSBarry Smith if (a->imax) PetscFree(a->imax); 5870452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 58876dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 5890452661fSBarry Smith PetscFree(a); 590f2655603SLois Curfman McInnes PLogObjectDestroy(A); 591f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 59217ab2063SBarry Smith return 0; 59317ab2063SBarry Smith } 59417ab2063SBarry Smith 595416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 59617ab2063SBarry Smith { 59717ab2063SBarry Smith return 0; 59817ab2063SBarry Smith } 59917ab2063SBarry Smith 600416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 60117ab2063SBarry Smith { 602416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6036d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 6046d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 6056d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 6066d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 6076d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 6086d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 6096d4a8577SBarry Smith op == MAT_SYMMETRIC || 6106d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 6116d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 61294a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 6136d4a8577SBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) 6146d4a8577SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:MAT_NO_NEW_DIAGONALS");} 6156d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 6166d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 6176d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 6186d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 6196d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 620e2f28af5SBarry Smith else 6214d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 62217ab2063SBarry Smith return 0; 62317ab2063SBarry Smith } 62417ab2063SBarry Smith 625416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 62617ab2063SBarry Smith { 627416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 628416022c9SBarry Smith int i,j, n,shift = a->indexshift; 62917ab2063SBarry Smith Scalar *x, zero = 0.0; 63017ab2063SBarry Smith 63117ab2063SBarry Smith VecSet(&zero,v); 63217ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 633416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 634416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 635416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 636416022c9SBarry Smith if (a->j[j]+shift == i) { 637416022c9SBarry Smith x[i] = a->a[j]; 63817ab2063SBarry Smith break; 63917ab2063SBarry Smith } 64017ab2063SBarry Smith } 64117ab2063SBarry Smith } 64217ab2063SBarry Smith return 0; 64317ab2063SBarry Smith } 64417ab2063SBarry Smith 64517ab2063SBarry Smith /* -------------------------------------------------------*/ 64617ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 64717ab2063SBarry Smith /* -------------------------------------------------------*/ 64844cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 64917ab2063SBarry Smith { 650416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 65117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 652416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 65317ab2063SBarry Smith 65417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 655cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 65617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 65717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 658416022c9SBarry Smith idx = a->j + a->i[i] + shift; 659416022c9SBarry Smith v = a->a + a->i[i] + shift; 660416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 66117ab2063SBarry Smith alpha = x[i]; 66217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 66317ab2063SBarry Smith } 664416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 66517ab2063SBarry Smith return 0; 66617ab2063SBarry Smith } 667d5d45c9bSBarry Smith 66844cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 66917ab2063SBarry Smith { 670416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 67117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 672416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 67317ab2063SBarry Smith 67417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 67517ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 67617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 67717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 678416022c9SBarry Smith idx = a->j + a->i[i] + shift; 679416022c9SBarry Smith v = a->a + a->i[i] + shift; 680416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 68117ab2063SBarry Smith alpha = x[i]; 68217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 68317ab2063SBarry Smith } 68417ab2063SBarry Smith return 0; 68517ab2063SBarry Smith } 68617ab2063SBarry Smith 68744cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 68817ab2063SBarry Smith { 689416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 69017ab2063SBarry Smith Scalar *x, *y, *v, sum; 691416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 69217ab2063SBarry Smith 69317ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 69417ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 695416022c9SBarry Smith idx = a->j; 696416022c9SBarry Smith v = a->a; 697416022c9SBarry Smith ii = a->i; 69817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 699416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 70017ab2063SBarry Smith sum = 0.0; 70117ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 70217ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 703416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 70417ab2063SBarry Smith y[i] = sum; 70517ab2063SBarry Smith } 706416022c9SBarry Smith PLogFlops(2*a->nz - m); 70717ab2063SBarry Smith return 0; 70817ab2063SBarry Smith } 70917ab2063SBarry Smith 71044cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 71117ab2063SBarry Smith { 712416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 71317ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 714cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 71517ab2063SBarry Smith 71617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 71717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 718cddf8d76SBarry Smith idx = a->j; 719cddf8d76SBarry Smith v = a->a; 720cddf8d76SBarry Smith ii = a->i; 72117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 722cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 72317ab2063SBarry Smith sum = y[i]; 724cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 725cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 72617ab2063SBarry Smith z[i] = sum; 72717ab2063SBarry Smith } 728416022c9SBarry Smith PLogFlops(2*a->nz); 72917ab2063SBarry Smith return 0; 73017ab2063SBarry Smith } 73117ab2063SBarry Smith 73217ab2063SBarry Smith /* 73317ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 73417ab2063SBarry Smith */ 73517ab2063SBarry Smith 736416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 73717ab2063SBarry Smith { 738416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 739416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 74017ab2063SBarry Smith 7410452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 742416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 743416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 744416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 745416022c9SBarry Smith if (a->j[j]+shift == i) { 74617ab2063SBarry Smith diag[i] = j - shift; 74717ab2063SBarry Smith break; 74817ab2063SBarry Smith } 74917ab2063SBarry Smith } 75017ab2063SBarry Smith } 751416022c9SBarry Smith a->diag = diag; 75217ab2063SBarry Smith return 0; 75317ab2063SBarry Smith } 75417ab2063SBarry Smith 75544cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 75617ab2063SBarry Smith double fshift,int its,Vec xx) 75717ab2063SBarry Smith { 758416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 759416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 760d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 76117ab2063SBarry Smith 76217ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 763416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 764416022c9SBarry Smith diag = a->diag; 765416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 76617ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 76717ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 76817ab2063SBarry Smith bs = b + shift; 76917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 770416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 771416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 772416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 773416022c9SBarry Smith v = a->a + diag[i] + (!shift); 77417ab2063SBarry Smith sum = b[i]*d/omega; 77517ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 77617ab2063SBarry Smith x[i] = sum; 77717ab2063SBarry Smith } 77817ab2063SBarry Smith return 0; 77917ab2063SBarry Smith } 78017ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 78117ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 78217ab2063SBarry Smith } 783416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 78417ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 78517ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 78617ab2063SBarry Smith 78717ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 78817ab2063SBarry Smith 78917ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 79017ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 79117ab2063SBarry Smith is the relaxation factor. 79217ab2063SBarry Smith */ 7930452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 79417ab2063SBarry Smith scale = (2.0/omega) - 1.0; 79517ab2063SBarry Smith 79617ab2063SBarry Smith /* x = (E + U)^{-1} b */ 79717ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 798416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 799416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 800416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 801416022c9SBarry Smith v = a->a + diag[i] + (!shift); 80217ab2063SBarry Smith sum = b[i]; 80317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 80417ab2063SBarry Smith x[i] = omega*(sum/d); 80517ab2063SBarry Smith } 80617ab2063SBarry Smith 80717ab2063SBarry Smith /* t = b - (2*E - D)x */ 808416022c9SBarry Smith v = a->a; 80917ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 81017ab2063SBarry Smith 81117ab2063SBarry Smith /* t = (E + L)^{-1}t */ 812416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 813416022c9SBarry Smith diag = a->diag; 81417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 815416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 816416022c9SBarry Smith n = diag[i] - a->i[i]; 817416022c9SBarry Smith idx = a->j + a->i[i] + shift; 818416022c9SBarry Smith v = a->a + a->i[i] + shift; 81917ab2063SBarry Smith sum = t[i]; 82017ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 82117ab2063SBarry Smith t[i] = omega*(sum/d); 82217ab2063SBarry Smith } 82317ab2063SBarry Smith 82417ab2063SBarry Smith /* x = x + t */ 82517ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 8260452661fSBarry Smith PetscFree(t); 82717ab2063SBarry Smith return 0; 82817ab2063SBarry Smith } 82917ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 83017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 83117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 832416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 833416022c9SBarry Smith n = diag[i] - a->i[i]; 834416022c9SBarry Smith idx = a->j + a->i[i] + shift; 835416022c9SBarry Smith v = a->a + a->i[i] + shift; 83617ab2063SBarry Smith sum = b[i]; 83717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 83817ab2063SBarry Smith x[i] = omega*(sum/d); 83917ab2063SBarry Smith } 84017ab2063SBarry Smith xb = x; 84117ab2063SBarry Smith } 84217ab2063SBarry Smith else xb = b; 84317ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 84417ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 84517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 846416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 84717ab2063SBarry Smith } 84817ab2063SBarry Smith } 84917ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 85017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 851416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 852416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 853416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 854416022c9SBarry Smith v = a->a + diag[i] + (!shift); 85517ab2063SBarry Smith sum = xb[i]; 85617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 85717ab2063SBarry Smith x[i] = omega*(sum/d); 85817ab2063SBarry Smith } 85917ab2063SBarry Smith } 86017ab2063SBarry Smith its--; 86117ab2063SBarry Smith } 86217ab2063SBarry Smith while (its--) { 86317ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 86417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 865416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 866416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 867416022c9SBarry Smith idx = a->j + a->i[i] + shift; 868416022c9SBarry Smith v = a->a + a->i[i] + shift; 86917ab2063SBarry Smith sum = b[i]; 87017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 871*7e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 87217ab2063SBarry Smith } 87317ab2063SBarry Smith } 87417ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 87517ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 876416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 877416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 878416022c9SBarry Smith idx = a->j + a->i[i] + shift; 879416022c9SBarry Smith v = a->a + a->i[i] + shift; 88017ab2063SBarry Smith sum = b[i]; 88117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 882*7e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 88317ab2063SBarry Smith } 88417ab2063SBarry Smith } 88517ab2063SBarry Smith } 88617ab2063SBarry Smith return 0; 88717ab2063SBarry Smith } 88817ab2063SBarry Smith 8894e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 89017ab2063SBarry Smith { 891416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 8924e220ebcSLois Curfman McInnes 8934e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 8944e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 8954e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 8964e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 8974e220ebcSLois Curfman McInnes info->block_size = 1.0; 8984e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 8994e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 9004e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 9014e220ebcSLois Curfman McInnes /* if (info->nz_unneeded != A->info.nz_unneeded) 9024e220ebcSLois Curfman McInnes printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 9034e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 9044e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 9054e220ebcSLois Curfman McInnes info->memory = A->mem; 9064e220ebcSLois Curfman McInnes if (A->factor) { 9074e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 9084e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 9094e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 9104e220ebcSLois Curfman McInnes } else { 9114e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 9124e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 9134e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 9144e220ebcSLois Curfman McInnes } 91517ab2063SBarry Smith return 0; 91617ab2063SBarry Smith } 91717ab2063SBarry Smith 91817ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 91917ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 92017ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 92117ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 92217ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 92317ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 92417ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 92517ab2063SBarry Smith 92617ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 92717ab2063SBarry Smith { 928416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 929416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 93017ab2063SBarry Smith 93177c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 93217ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 93317ab2063SBarry Smith if (diag) { 93417ab2063SBarry Smith for ( i=0; i<N; i++ ) { 935416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 936416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 937416022c9SBarry Smith a->ilen[rows[i]] = 1; 938416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 939416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 94017ab2063SBarry Smith } 94117ab2063SBarry Smith else { 94217ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 94317ab2063SBarry Smith CHKERRQ(ierr); 94417ab2063SBarry Smith } 94517ab2063SBarry Smith } 94617ab2063SBarry Smith } 94717ab2063SBarry Smith else { 94817ab2063SBarry Smith for ( i=0; i<N; i++ ) { 949416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 950416022c9SBarry Smith a->ilen[rows[i]] = 0; 95117ab2063SBarry Smith } 95217ab2063SBarry Smith } 95317ab2063SBarry Smith ISRestoreIndices(is,&rows); 9546d4a8577SBarry Smith ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 9556d4a8577SBarry Smith ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 95617ab2063SBarry Smith return 0; 95717ab2063SBarry Smith } 95817ab2063SBarry Smith 959416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 96017ab2063SBarry Smith { 961416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 962416022c9SBarry Smith *m = a->m; *n = a->n; 96317ab2063SBarry Smith return 0; 96417ab2063SBarry Smith } 96517ab2063SBarry Smith 966416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 96717ab2063SBarry Smith { 968416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 969416022c9SBarry Smith *m = 0; *n = a->m; 97017ab2063SBarry Smith return 0; 97117ab2063SBarry Smith } 9724e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 97317ab2063SBarry Smith { 974416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 975c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 97617ab2063SBarry Smith 977416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 97817ab2063SBarry Smith 979416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 980416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 98117ab2063SBarry Smith if (idx) { 982416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 9834e093b46SBarry Smith if (*nz && shift) { 9840452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 98517ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 9864e093b46SBarry Smith } else if (*nz) { 9874e093b46SBarry Smith *idx = itmp; 98817ab2063SBarry Smith } 98917ab2063SBarry Smith else *idx = 0; 99017ab2063SBarry Smith } 99117ab2063SBarry Smith return 0; 99217ab2063SBarry Smith } 99317ab2063SBarry Smith 9944e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 99517ab2063SBarry Smith { 9964e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 9974e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 99817ab2063SBarry Smith return 0; 99917ab2063SBarry Smith } 100017ab2063SBarry Smith 1001cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 100217ab2063SBarry Smith { 1003416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1004416022c9SBarry Smith Scalar *v = a->a; 100517ab2063SBarry Smith double sum = 0.0; 1006416022c9SBarry Smith int i, j,shift = a->indexshift; 100717ab2063SBarry Smith 100817ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1009416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 101017ab2063SBarry Smith #if defined(PETSC_COMPLEX) 101117ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 101217ab2063SBarry Smith #else 101317ab2063SBarry Smith sum += (*v)*(*v); v++; 101417ab2063SBarry Smith #endif 101517ab2063SBarry Smith } 101617ab2063SBarry Smith *norm = sqrt(sum); 101717ab2063SBarry Smith } 101817ab2063SBarry Smith else if (type == NORM_1) { 101917ab2063SBarry Smith double *tmp; 1020416022c9SBarry Smith int *jj = a->j; 10210452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 1022cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 102317ab2063SBarry Smith *norm = 0.0; 1024416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1025a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 102617ab2063SBarry Smith } 1027416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 102817ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 102917ab2063SBarry Smith } 10300452661fSBarry Smith PetscFree(tmp); 103117ab2063SBarry Smith } 103217ab2063SBarry Smith else if (type == NORM_INFINITY) { 103317ab2063SBarry Smith *norm = 0.0; 1034416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1035416022c9SBarry Smith v = a->a + a->i[j] + shift; 103617ab2063SBarry Smith sum = 0.0; 1037416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1038cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 103917ab2063SBarry Smith } 104017ab2063SBarry Smith if (sum > *norm) *norm = sum; 104117ab2063SBarry Smith } 104217ab2063SBarry Smith } 104317ab2063SBarry Smith else { 104448d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 104517ab2063SBarry Smith } 104617ab2063SBarry Smith return 0; 104717ab2063SBarry Smith } 104817ab2063SBarry Smith 1049416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 105017ab2063SBarry Smith { 1051416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1052416022c9SBarry Smith Mat C; 1053416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1054416022c9SBarry Smith int shift = a->indexshift; 1055d5d45c9bSBarry Smith Scalar *array = a->a; 105617ab2063SBarry Smith 10573638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 1058dfa27b74SSatish Balay SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place"); 10590452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1060cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 106117ab2063SBarry Smith if (shift) { 106217ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 106317ab2063SBarry Smith } 106417ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1065416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 10660452661fSBarry Smith PetscFree(col); 106717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 106817ab2063SBarry Smith len = ai[i+1]-ai[i]; 1069416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 107017ab2063SBarry Smith array += len; aj += len; 107117ab2063SBarry Smith } 107217ab2063SBarry Smith if (shift) { 107317ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 107417ab2063SBarry Smith } 107517ab2063SBarry Smith 10766d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10776d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 107817ab2063SBarry Smith 10793638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1080416022c9SBarry Smith *B = C; 108117ab2063SBarry Smith } else { 1082416022c9SBarry Smith /* This isn't really an in-place transpose */ 10830452661fSBarry Smith PetscFree(a->a); 10840452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 10850452661fSBarry Smith if (a->diag) PetscFree(a->diag); 10860452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 10870452661fSBarry Smith if (a->imax) PetscFree(a->imax); 10880452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 10891073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 10900452661fSBarry Smith PetscFree(a); 1091416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 10920452661fSBarry Smith PetscHeaderDestroy(C); 109317ab2063SBarry Smith } 109417ab2063SBarry Smith return 0; 109517ab2063SBarry Smith } 109617ab2063SBarry Smith 1097f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 109817ab2063SBarry Smith { 1099416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 110017ab2063SBarry Smith Scalar *l,*r,x,*v; 1101d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 110217ab2063SBarry Smith 110317ab2063SBarry Smith if (ll) { 110417ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 1105f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 1106416022c9SBarry Smith v = a->a; 110717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 110817ab2063SBarry Smith x = l[i]; 1109416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 111017ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 111117ab2063SBarry Smith } 111244cd7ae7SLois Curfman McInnes PLogFlops(nz); 111317ab2063SBarry Smith } 111417ab2063SBarry Smith if (rr) { 111517ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 1116f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 1117416022c9SBarry Smith v = a->a; jj = a->j; 111817ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 111917ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 112017ab2063SBarry Smith } 112144cd7ae7SLois Curfman McInnes PLogFlops(nz); 112217ab2063SBarry Smith } 112317ab2063SBarry Smith return 0; 112417ab2063SBarry Smith } 112517ab2063SBarry Smith 1126cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 112717ab2063SBarry Smith { 1128db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 112902834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 113099141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1131a2744918SBarry Smith register int sum,lensi; 113299141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 113399141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 113499141d43SSatish Balay Scalar *a_new,*mat_a; 1135416022c9SBarry Smith Mat C; 113617ab2063SBarry Smith 1137b48a1e75SSatish Balay ierr = ISSorted(isrow,(PetscTruth*)&i); 1138b48a1e75SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:ISrow is not sorted"); 113999141d43SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i); 1140b48a1e75SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IScol is not sorted"); 114199141d43SSatish Balay 114217ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 114317ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 114417ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 114517ab2063SBarry Smith 11467264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 114702834360SBarry Smith /* special case of contiguous rows */ 114857faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 114902834360SBarry Smith starts = lens + ncols; 115002834360SBarry Smith /* loop over new rows determining lens and starting points */ 115102834360SBarry Smith for (i=0; i<nrows; i++) { 1152a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1153a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 115402834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1155d8ced48eSBarry Smith if (aj[k]+shift >= first) { 115602834360SBarry Smith starts[i] = k; 115702834360SBarry Smith break; 115802834360SBarry Smith } 115902834360SBarry Smith } 1160a2744918SBarry Smith sum = 0; 116102834360SBarry Smith while (k < kend) { 1162d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1163a2744918SBarry Smith sum++; 116402834360SBarry Smith } 1165a2744918SBarry Smith lens[i] = sum; 116602834360SBarry Smith } 116702834360SBarry Smith /* create submatrix */ 1168cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 116908480c60SBarry Smith int n_cols,n_rows; 117008480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 117108480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1172d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 117308480c60SBarry Smith C = *B; 117408480c60SBarry Smith } 117508480c60SBarry Smith else { 117602834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 117708480c60SBarry Smith } 1178db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1179db02288aSLois Curfman McInnes 118002834360SBarry Smith /* loop over rows inserting into submatrix */ 1181db02288aSLois Curfman McInnes a_new = c->a; 1182db02288aSLois Curfman McInnes j_new = c->j; 1183db02288aSLois Curfman McInnes i_new = c->i; 1184db02288aSLois Curfman McInnes i_new[0] = -shift; 118502834360SBarry Smith for (i=0; i<nrows; i++) { 1186a2744918SBarry Smith ii = starts[i]; 1187a2744918SBarry Smith lensi = lens[i]; 1188a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1189a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 119002834360SBarry Smith } 1191a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1192a2744918SBarry Smith a_new += lensi; 1193a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1194a2744918SBarry Smith c->ilen[i] = lensi; 119502834360SBarry Smith } 11960452661fSBarry Smith PetscFree(lens); 119702834360SBarry Smith } 119802834360SBarry Smith else { 119902834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 12000452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 120102834360SBarry Smith ssmap = smap + shift; 120299141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1203cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 120417ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 120502834360SBarry Smith /* determine lens of each row */ 120602834360SBarry Smith for (i=0; i<nrows; i++) { 1207d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 120802834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 120902834360SBarry Smith lens[i] = 0; 121002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1211d8ced48eSBarry Smith if (ssmap[aj[k]]) { 121202834360SBarry Smith lens[i]++; 121302834360SBarry Smith } 121402834360SBarry Smith } 121502834360SBarry Smith } 121617ab2063SBarry Smith /* Create and fill new matrix */ 1217a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 121899141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 121999141d43SSatish Balay 122099141d43SSatish Balay if (c->m != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 122199141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 122299141d43SSatish Balay SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros"); 122399141d43SSatish Balay } 122499141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 122508480c60SBarry Smith C = *B; 122608480c60SBarry Smith } 122708480c60SBarry Smith else { 122802834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 122908480c60SBarry Smith } 123099141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 123117ab2063SBarry Smith for (i=0; i<nrows; i++) { 123299141d43SSatish Balay row = irow[i]; 123317ab2063SBarry Smith nznew = 0; 123499141d43SSatish Balay kstart = ai[row]+shift; 123599141d43SSatish Balay kend = kstart + a->ilen[row]; 123699141d43SSatish Balay mat_i = c->i[i]+shift; 123799141d43SSatish Balay mat_j = c->j + mat_i; 123899141d43SSatish Balay mat_a = c->a + mat_i; 123999141d43SSatish Balay mat_ilen = c->ilen + i; 124017ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 124199141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 124299141d43SSatish Balay *mat_j++ = tcol - (!shift); 124399141d43SSatish Balay *mat_a++ = a->a[k]; 124499141d43SSatish Balay (*mat_ilen)++; 124599141d43SSatish Balay 124617ab2063SBarry Smith } 124717ab2063SBarry Smith } 124817ab2063SBarry Smith } 124902834360SBarry Smith /* Free work space */ 125002834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 125199141d43SSatish Balay PetscFree(smap); PetscFree(lens); 125202834360SBarry Smith } 12536d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 12546d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 125517ab2063SBarry Smith 125617ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1257416022c9SBarry Smith *B = C; 125817ab2063SBarry Smith return 0; 125917ab2063SBarry Smith } 126017ab2063SBarry Smith 1261a871dcd8SBarry Smith /* 126263b91edcSBarry Smith note: This can only work for identity for row and col. It would 126363b91edcSBarry Smith be good to check this and otherwise generate an error. 1264a871dcd8SBarry Smith */ 126563b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1266a871dcd8SBarry Smith { 126763b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 126808480c60SBarry Smith int ierr; 126963b91edcSBarry Smith Mat outA; 127063b91edcSBarry Smith 1271a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1272a871dcd8SBarry Smith 127363b91edcSBarry Smith outA = inA; 127463b91edcSBarry Smith inA->factor = FACTOR_LU; 127563b91edcSBarry Smith a->row = row; 127663b91edcSBarry Smith a->col = col; 127763b91edcSBarry Smith 12780452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 127963b91edcSBarry Smith 128008480c60SBarry Smith if (!a->diag) { 128108480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 128263b91edcSBarry Smith } 128363b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1284a871dcd8SBarry Smith return 0; 1285a871dcd8SBarry Smith } 1286a871dcd8SBarry Smith 1287f0b747eeSBarry Smith #include "pinclude/plapack.h" 1288f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1289f0b747eeSBarry Smith { 1290f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1291f0b747eeSBarry Smith int one = 1; 1292f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1293f0b747eeSBarry Smith PLogFlops(a->nz); 1294f0b747eeSBarry Smith return 0; 1295f0b747eeSBarry Smith } 1296f0b747eeSBarry Smith 1297cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1298cddf8d76SBarry Smith Mat **B) 1299cddf8d76SBarry Smith { 1300cddf8d76SBarry Smith int ierr,i; 1301cddf8d76SBarry Smith 1302cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 13030452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1304cddf8d76SBarry Smith } 1305cddf8d76SBarry Smith 1306cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1307905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1308cddf8d76SBarry Smith } 1309cddf8d76SBarry Smith return 0; 1310cddf8d76SBarry Smith } 1311cddf8d76SBarry Smith 13125a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 13135a838052SSatish Balay { 13145a838052SSatish Balay *bs = 1; 13155a838052SSatish Balay return 0; 13165a838052SSatish Balay } 13175a838052SSatish Balay 1318e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 13194dcbc457SBarry Smith { 1320e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 132106763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 13228a047759SSatish Balay int start, end, *ai, *aj; 132306763907SSatish Balay char *table; 13248a047759SSatish Balay shift = a->indexshift; 1325e4d965acSSatish Balay m = a->m; 1326e4d965acSSatish Balay ai = a->i; 13278a047759SSatish Balay aj = a->j+shift; 13288a047759SSatish Balay 1329e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 133006763907SSatish Balay 133106763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 133206763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 133306763907SSatish Balay 1334e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1335e4d965acSSatish Balay /* Initialise the two local arrays */ 1336e4d965acSSatish Balay isz = 0; 133706763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1338e4d965acSSatish Balay 1339e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 13404dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 134177c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1342e4d965acSSatish Balay 1343e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1344e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 134506763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 13464dcbc457SBarry Smith } 134706763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 134806763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1349e4d965acSSatish Balay 135004a348a9SBarry Smith k = 0; 135104a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 135204a348a9SBarry Smith n = isz; 135306763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1354e4d965acSSatish Balay row = nidx[k]; 1355e4d965acSSatish Balay start = ai[row]; 1356e4d965acSSatish Balay end = ai[row+1]; 135704a348a9SBarry Smith for ( l = start; l<end ; l++){ 13588a047759SSatish Balay val = aj[l] + shift; 135906763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1360e4d965acSSatish Balay } 1361e4d965acSSatish Balay } 1362e4d965acSSatish Balay } 1363537820f0SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1364e4d965acSSatish Balay } 136504a348a9SBarry Smith PetscFree(table); 136606763907SSatish Balay PetscFree(nidx); 1367e4d965acSSatish Balay return 0; 13684dcbc457SBarry Smith } 136917ab2063SBarry Smith 1370682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1371682d7d0cSBarry Smith { 1372682d7d0cSBarry Smith static int called = 0; 1373682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1374682d7d0cSBarry Smith 1375682d7d0cSBarry Smith if (called) return 0; else called = 1; 137677c4ece6SBarry Smith PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 137777c4ece6SBarry Smith PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>\n"); 137877c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 137977c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_no_inode - Do not use inodes\n"); 138077c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1381682d7d0cSBarry Smith #if defined(HAVE_ESSL) 138277c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1383682d7d0cSBarry Smith #endif 1384682d7d0cSBarry Smith return 0; 1385682d7d0cSBarry Smith } 1386205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1387682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 138817ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 138917ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1390416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1391416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 139217ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 139317ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 139417ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 139517ab2063SBarry Smith MatRelax_SeqAIJ, 139617ab2063SBarry Smith MatTranspose_SeqAIJ, 13977264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1398f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 139917ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 140017ab2063SBarry Smith MatCompress_SeqAIJ, 140117ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 140217ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 140317ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 140417ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 140517ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 14063d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1407cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 14087eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1409682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1410f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 14115a838052SSatish Balay MatScale_SeqAIJ,0,0, 14126945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 14136945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 14143b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 14153b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 14163b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 14173b2fbd54SBarry Smith MatRestoreColumnIJ_SeqAIJ}; 141817ab2063SBarry Smith 141917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 142017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 142117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 142217ab2063SBarry Smith 142317ab2063SBarry Smith /*@C 1424682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 14250d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 14266e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 14272bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 14282bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 142917ab2063SBarry Smith 143017ab2063SBarry Smith Input Parameters: 143117ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 143217ab2063SBarry Smith . m - number of rows 143317ab2063SBarry Smith . n - number of columns 143417ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 14352bd5e0b2SLois Curfman McInnes . nzz - array containing the number of nonzeros in the various rows 14362bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 143717ab2063SBarry Smith 143817ab2063SBarry Smith Output Parameter: 1439416022c9SBarry Smith . A - the matrix 144017ab2063SBarry Smith 144117ab2063SBarry Smith Notes: 144217ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 144317ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 14440002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 144544cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 144617ab2063SBarry Smith 144717ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1448a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14493d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 14503d323bbdSBarry Smith will get TERRIBLE performance, see the users' manual chapter on 14510d15e28bSLois Curfman McInnes matrices and the file $(PETSC_DIR)/Performance. 145217ab2063SBarry Smith 1453682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1454682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1455682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 14566c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 14576c7ebb05SLois Curfman McInnes 14586c7ebb05SLois Curfman McInnes Options Database Keys: 14596c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 14606e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 14616e62573dSLois Curfman McInnes $ (max limit=5) 14626e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 14636e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 14646e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 146517ab2063SBarry Smith 146617ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 146717ab2063SBarry Smith @*/ 1468416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 146917ab2063SBarry Smith { 1470416022c9SBarry Smith Mat B; 1471416022c9SBarry Smith Mat_SeqAIJ *b; 14726945ee14SBarry Smith int i, len, ierr, flg,size; 14736945ee14SBarry Smith 14746945ee14SBarry Smith MPI_Comm_size(comm,&size); 14756945ee14SBarry Smith if (size > 1) SETERRQ(1,"MatCreateSeqAIJ:Comm must be of size 1"); 1476d5d45c9bSBarry Smith 1477416022c9SBarry Smith *A = 0; 14780452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1479416022c9SBarry Smith PLogObjectCreate(B); 14800452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 148144cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1482416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1483416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1484416022c9SBarry Smith B->view = MatView_SeqAIJ; 1485416022c9SBarry Smith B->factor = 0; 1486416022c9SBarry Smith B->lupivotthreshold = 1.0; 14877a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 148869957df2SSatish Balay &flg); CHKERRQ(ierr); 14897a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 14907a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 14917a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1492416022c9SBarry Smith b->row = 0; 1493416022c9SBarry Smith b->col = 0; 1494416022c9SBarry Smith b->indexshift = 0; 1495b810aeb4SBarry Smith b->reallocs = 0; 149669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 149769957df2SSatish Balay if (flg) b->indexshift = -1; 149817ab2063SBarry Smith 149944cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 150044cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 15010452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1502b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 15037b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 15047b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1505416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 150617ab2063SBarry Smith nz = nz*m; 150717ab2063SBarry Smith } 150817ab2063SBarry Smith else { 150917ab2063SBarry Smith nz = 0; 1510416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 151117ab2063SBarry Smith } 151217ab2063SBarry Smith 151317ab2063SBarry Smith /* allocate the matrix space */ 1514416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 15150452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1516416022c9SBarry Smith b->j = (int *) (b->a + nz); 1517cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1518416022c9SBarry Smith b->i = b->j + nz; 1519416022c9SBarry Smith b->singlemalloc = 1; 152017ab2063SBarry Smith 1521416022c9SBarry Smith b->i[0] = -b->indexshift; 152217ab2063SBarry Smith for (i=1; i<m+1; i++) { 1523416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 152417ab2063SBarry Smith } 152517ab2063SBarry Smith 1526416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 15270452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1528416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1529416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 153017ab2063SBarry Smith 1531416022c9SBarry Smith b->nz = 0; 1532416022c9SBarry Smith b->maxnz = nz; 1533416022c9SBarry Smith b->sorted = 0; 1534416022c9SBarry Smith b->roworiented = 1; 1535416022c9SBarry Smith b->nonew = 0; 1536416022c9SBarry Smith b->diag = 0; 1537416022c9SBarry Smith b->solve_work = 0; 153871bd300dSLois Curfman McInnes b->spptr = 0; 1539754ec7b1SSatish Balay b->inode.node_count = 0; 1540754ec7b1SSatish Balay b->inode.size = 0; 15416c7ebb05SLois Curfman McInnes b->inode.limit = 5; 15426c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 15434e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 154417ab2063SBarry Smith 1545416022c9SBarry Smith *A = B; 15464e220ebcSLois Curfman McInnes 15474b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 15484b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 154969957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 155069957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 15514b14c69eSBarry Smith #endif 155269957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 155369957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 155469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 155569957df2SSatish Balay if (flg) { 1556416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1557416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 155817ab2063SBarry Smith } 155969957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 156069957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 156117ab2063SBarry Smith return 0; 156217ab2063SBarry Smith } 156317ab2063SBarry Smith 15643d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 156517ab2063SBarry Smith { 1566416022c9SBarry Smith Mat C; 1567416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 156808480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 156917ab2063SBarry Smith 15704043dd9cSLois Curfman McInnes *B = 0; 15710452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1572416022c9SBarry Smith PLogObjectCreate(C); 15730452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 157441c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1575416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1576416022c9SBarry Smith C->view = MatView_SeqAIJ; 1577416022c9SBarry Smith C->factor = A->factor; 1578416022c9SBarry Smith c->row = 0; 1579416022c9SBarry Smith c->col = 0; 1580416022c9SBarry Smith c->indexshift = shift; 1581c456f294SBarry Smith C->assembled = PETSC_TRUE; 158217ab2063SBarry Smith 158344cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 158444cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 158544cd7ae7SLois Curfman McInnes C->M = a->m; 158644cd7ae7SLois Curfman McInnes C->N = a->n; 158717ab2063SBarry Smith 15880452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 15890452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 159017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1591416022c9SBarry Smith c->imax[i] = a->imax[i]; 1592416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 159317ab2063SBarry Smith } 159417ab2063SBarry Smith 159517ab2063SBarry Smith /* allocate the matrix space */ 1596416022c9SBarry Smith c->singlemalloc = 1; 1597416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 15980452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1599416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1600416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1601416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 160217ab2063SBarry Smith if (m > 0) { 1603416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 160408480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1605416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 160617ab2063SBarry Smith } 160708480c60SBarry Smith } 160817ab2063SBarry Smith 1609416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1610416022c9SBarry Smith c->sorted = a->sorted; 1611416022c9SBarry Smith c->roworiented = a->roworiented; 1612416022c9SBarry Smith c->nonew = a->nonew; 16137a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 161417ab2063SBarry Smith 1615416022c9SBarry Smith if (a->diag) { 16160452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1617416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 161817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1619416022c9SBarry Smith c->diag[i] = a->diag[i]; 162017ab2063SBarry Smith } 162117ab2063SBarry Smith } 1622416022c9SBarry Smith else c->diag = 0; 16236c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 16246c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1625754ec7b1SSatish Balay if (a->inode.size){ 1626daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 1627754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1628daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 1629754ec7b1SSatish Balay } else { 1630754ec7b1SSatish Balay c->inode.size = 0; 1631754ec7b1SSatish Balay c->inode.node_count = 0; 1632754ec7b1SSatish Balay } 1633416022c9SBarry Smith c->nz = a->nz; 1634416022c9SBarry Smith c->maxnz = a->maxnz; 1635416022c9SBarry Smith c->solve_work = 0; 163676dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1637754ec7b1SSatish Balay 1638416022c9SBarry Smith *B = C; 163917ab2063SBarry Smith return 0; 164017ab2063SBarry Smith } 164117ab2063SBarry Smith 164219bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 164317ab2063SBarry Smith { 1644416022c9SBarry Smith Mat_SeqAIJ *a; 1645416022c9SBarry Smith Mat B; 164617699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1647bcd2baecSBarry Smith MPI_Comm comm; 164817ab2063SBarry Smith 164919bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 165017699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 165117699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 165290ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 165377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 165448d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 165517ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 165617ab2063SBarry Smith 165717ab2063SBarry Smith /* read in row lengths */ 16580452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 165977c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 166017ab2063SBarry Smith 166117ab2063SBarry Smith /* create our matrix */ 1662416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1663416022c9SBarry Smith B = *A; 1664416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1665416022c9SBarry Smith shift = a->indexshift; 166617ab2063SBarry Smith 166717ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 166877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 166917ab2063SBarry Smith if (shift) { 167017ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1671416022c9SBarry Smith a->j[i] += 1; 167217ab2063SBarry Smith } 167317ab2063SBarry Smith } 167417ab2063SBarry Smith 167517ab2063SBarry Smith /* read in nonzero values */ 167677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 167717ab2063SBarry Smith 167817ab2063SBarry Smith /* set matrix "i" values */ 1679416022c9SBarry Smith a->i[0] = -shift; 168017ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1681416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1682416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 168317ab2063SBarry Smith } 16840452661fSBarry Smith PetscFree(rowlengths); 168517ab2063SBarry Smith 16866d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 16876d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 168817ab2063SBarry Smith return 0; 168917ab2063SBarry Smith } 169017ab2063SBarry Smith 1691686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 16927264ac53SSatish Balay { 16937264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 16947264ac53SSatish Balay 1695bcd2baecSBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type"); 16967264ac53SSatish Balay 16977264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 16987264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1699bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 170077c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1701bcd2baecSBarry Smith } 17027264ac53SSatish Balay 17037264ac53SSatish Balay /* if the a->i are the same */ 1704bcd2baecSBarry Smith if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 170577c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 17067264ac53SSatish Balay } 17077264ac53SSatish Balay 17087264ac53SSatish Balay /* if a->j are the same */ 1709bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 171077c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1711bcd2baecSBarry Smith } 1712bcd2baecSBarry Smith 1713bcd2baecSBarry Smith /* if a->a are the same */ 171419bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 171577c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 17167264ac53SSatish Balay } 171777c4ece6SBarry Smith *flg = PETSC_TRUE; 17187264ac53SSatish Balay return 0; 17197264ac53SSatish Balay 17207264ac53SSatish Balay } 1721