16d84be18SBarry Smith 26945ee14SBarry Smith 317ab2063SBarry Smith #ifndef lint 4*dd5f02e7SSatish Balay static char vcid[] = "$Id: aij.c,v 1.198 1996/11/29 22:26:43 curfman Exp balay $"; 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 3031625ec5SSatish Balay static int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,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 3631625ec5SSatish Balay *m = A->m; 376945ee14SBarry Smith if (!ia) return 0; 386945ee14SBarry Smith ishift = a->indexshift; 396945ee14SBarry Smith if (symmetric) { 4031625ec5SSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 416945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 4231625ec5SSatish Balay int nz = a->i[a->m]; 433b2fbd54SBarry Smith /* malloc space and subtract 1 from i and j indices */ 4431625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+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; 4731625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1; 486945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 4931625ec5SSatish Balay int nz = a->i[a->m] + 1; 503b2fbd54SBarry Smith /* malloc space and add 1 to i and j indices */ 5131625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+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; 5431625ec5SSatish Balay for ( i=0; i<a->m+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 7643a90d84SBarry Smith 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; 80a93ec695SBarry Smith int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m; 81a93ec695SBarry Smith int nz = a->i[m]+ishift,row,*jj,mr,col; 823b2fbd54SBarry Smith 833b2fbd54SBarry Smith *nn = A->n; 843b2fbd54SBarry Smith if (!ia) return 0; 853b2fbd54SBarry Smith if (symmetric) { 86179192dfSSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 873b2fbd54SBarry Smith } else { 8861d2ded1SBarry Smith collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths); 893b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 903b2fbd54SBarry Smith cia = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia); 91a93ec695SBarry Smith cja = (int *) PetscMalloc( (nz+1)*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; 102a93ec695SBarry Smith for ( row=0; row<m; row++ ) { 103a93ec695SBarry Smith mr = a->i[row+1] - a->i[row]; 104a93ec695SBarry Smith for ( i=0; i<mr; i++ ) { 1053b2fbd54SBarry Smith col = *jj++ + ishift; 1063b2fbd54SBarry Smith cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1073b2fbd54SBarry Smith } 1083b2fbd54SBarry Smith } 1093b2fbd54SBarry Smith PetscFree(collengths); 1103b2fbd54SBarry Smith *ia = cia; *ja = cja; 1113b2fbd54SBarry Smith } 1123b2fbd54SBarry Smith 1133b2fbd54SBarry Smith return 0; 1143b2fbd54SBarry Smith } 1153b2fbd54SBarry Smith 11643a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia, 1173b2fbd54SBarry Smith int **ja,PetscTruth *done) 1183b2fbd54SBarry Smith { 1193b2fbd54SBarry Smith if (!ia) return 0; 1203b2fbd54SBarry Smith 1213b2fbd54SBarry Smith PetscFree(*ia); 1223b2fbd54SBarry Smith PetscFree(*ja); 1233b2fbd54SBarry Smith 1243b2fbd54SBarry Smith return 0; 1253b2fbd54SBarry Smith } 1263b2fbd54SBarry Smith 127227d817aSBarry Smith #define CHUNKSIZE 15 12817ab2063SBarry Smith 12961d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 13017ab2063SBarry Smith { 131416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 132416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 1334b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 134d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 135416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 13617ab2063SBarry Smith 13717ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 138416022c9SBarry Smith row = im[k]; 1393b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 14017ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 141416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 1423b2fbd54SBarry Smith #endif 14317ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 14417ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 145416022c9SBarry Smith low = 0; 14617ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 1473b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 148416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 149416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 1503b2fbd54SBarry Smith #endif 1514b0e389bSBarry Smith col = in[l] - shift; 1524b0e389bSBarry Smith if (roworiented) { 1534b0e389bSBarry Smith value = *v++; 1544b0e389bSBarry Smith } 1554b0e389bSBarry Smith else { 1564b0e389bSBarry Smith value = v[k + l*m]; 1574b0e389bSBarry Smith } 158416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 159416022c9SBarry Smith while (high-low > 5) { 160416022c9SBarry Smith t = (low+high)/2; 161416022c9SBarry Smith if (rp[t] > col) high = t; 162416022c9SBarry Smith else low = t; 16317ab2063SBarry Smith } 164416022c9SBarry Smith for ( i=low; i<high; i++ ) { 16517ab2063SBarry Smith if (rp[i] > col) break; 16617ab2063SBarry Smith if (rp[i] == col) { 167416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 16817ab2063SBarry Smith else ap[i] = value; 16917ab2063SBarry Smith goto noinsert; 17017ab2063SBarry Smith } 17117ab2063SBarry Smith } 17217ab2063SBarry Smith if (nonew) goto noinsert; 17317ab2063SBarry Smith if (nrow >= rmax) { 17417ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 175416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 17617ab2063SBarry Smith Scalar *new_a; 17717ab2063SBarry Smith 17817ab2063SBarry Smith /* malloc new storage space */ 179416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1800452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 18117ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 18217ab2063SBarry Smith new_i = new_j + new_nz; 18317ab2063SBarry Smith 18417ab2063SBarry Smith /* copy over old data into new slots */ 18517ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 186416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 187416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 188416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 189416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 19017ab2063SBarry Smith len*sizeof(int)); 191416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 192416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 19317ab2063SBarry Smith len*sizeof(Scalar)); 19417ab2063SBarry Smith /* free up old matrix storage */ 1950452661fSBarry Smith PetscFree(a->a); 1960452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 197416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 198416022c9SBarry Smith a->singlemalloc = 1; 19917ab2063SBarry Smith 20017ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 201416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 202416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 203416022c9SBarry Smith a->maxnz += CHUNKSIZE; 204b810aeb4SBarry Smith a->reallocs++; 20517ab2063SBarry Smith } 206416022c9SBarry Smith N = nrow++ - 1; a->nz++; 207416022c9SBarry Smith /* shift up all the later entries in this row */ 208416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 20917ab2063SBarry Smith rp[ii+1] = rp[ii]; 21017ab2063SBarry Smith ap[ii+1] = ap[ii]; 21117ab2063SBarry Smith } 21217ab2063SBarry Smith rp[i] = col; 21317ab2063SBarry Smith ap[i] = value; 21417ab2063SBarry Smith noinsert:; 215416022c9SBarry Smith low = i + 1; 21617ab2063SBarry Smith } 21717ab2063SBarry Smith ailen[row] = nrow; 21817ab2063SBarry Smith } 21917ab2063SBarry Smith return 0; 22017ab2063SBarry Smith } 22117ab2063SBarry Smith 2227eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 2237eb43aa7SLois Curfman McInnes { 2247eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 225b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 2267eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 2277eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 2287eb43aa7SLois Curfman McInnes 2297eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 2307eb43aa7SLois Curfman McInnes row = im[k]; 2317eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 2327eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 2337eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 2347eb43aa7SLois Curfman McInnes nrow = ailen[row]; 2357eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 2367eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 2377eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 2387eb43aa7SLois Curfman McInnes col = in[l] - shift; 2397eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 2407eb43aa7SLois Curfman McInnes while (high-low > 5) { 2417eb43aa7SLois Curfman McInnes t = (low+high)/2; 2427eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 2437eb43aa7SLois Curfman McInnes else low = t; 2447eb43aa7SLois Curfman McInnes } 2457eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 2467eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2477eb43aa7SLois Curfman McInnes if (rp[i] == col) { 248b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2497eb43aa7SLois Curfman McInnes goto finished; 2507eb43aa7SLois Curfman McInnes } 2517eb43aa7SLois Curfman McInnes } 252b49de8d1SLois Curfman McInnes *v++ = zero; 2537eb43aa7SLois Curfman McInnes finished:; 2547eb43aa7SLois Curfman McInnes } 2557eb43aa7SLois Curfman McInnes } 2567eb43aa7SLois Curfman McInnes return 0; 2577eb43aa7SLois Curfman McInnes } 2587eb43aa7SLois Curfman McInnes 25917ab2063SBarry Smith #include "draw.h" 26017ab2063SBarry Smith #include "pinclude/pviewer.h" 26177c4ece6SBarry Smith #include "sys.h" 26217ab2063SBarry Smith 263416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 26417ab2063SBarry Smith { 265416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 266416022c9SBarry Smith int i, fd, *col_lens, ierr; 26717ab2063SBarry Smith 26890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2690452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 270416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 271416022c9SBarry Smith col_lens[1] = a->m; 272416022c9SBarry Smith col_lens[2] = a->n; 273416022c9SBarry Smith col_lens[3] = a->nz; 274416022c9SBarry Smith 275416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 276416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 277416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 27817ab2063SBarry Smith } 27977c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 2800452661fSBarry Smith PetscFree(col_lens); 281416022c9SBarry Smith 282416022c9SBarry Smith /* store column indices (zero start index) */ 283416022c9SBarry Smith if (a->indexshift) { 284416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 28517ab2063SBarry Smith } 28677c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 287416022c9SBarry Smith if (a->indexshift) { 288416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 28917ab2063SBarry Smith } 290416022c9SBarry Smith 291416022c9SBarry Smith /* store nonzero values */ 29277c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 29317ab2063SBarry Smith return 0; 29417ab2063SBarry Smith } 295416022c9SBarry Smith 296416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 297416022c9SBarry Smith { 298416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 299496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 30017ab2063SBarry Smith FILE *fd; 30117ab2063SBarry Smith char *outputname; 30217ab2063SBarry Smith 30390ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 304416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 30590ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 306a93ec695SBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 30795e01e2fSLois Curfman McInnes return 0; 30895e01e2fSLois Curfman McInnes } 309a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 310496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 311496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 312496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 31395e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 31495e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 31517ab2063SBarry Smith } 316a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 317416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3184e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 3194e220ebcSLois Curfman McInnes fprintf(fd,"zzz = zeros(%d,3);\n",a->nz); 32017ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 32117ab2063SBarry Smith 32217ab2063SBarry Smith for (i=0; i<m; i++) { 323416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 32417ab2063SBarry Smith #if defined(PETSC_COMPLEX) 3256945ee14SBarry Smith fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]), 326416022c9SBarry Smith imag(a->a[j])); 32717ab2063SBarry Smith #else 3287a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 32917ab2063SBarry Smith #endif 33017ab2063SBarry Smith } 33117ab2063SBarry Smith } 33217ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 33317ab2063SBarry Smith } 334a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 33544cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 33644cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 33744cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 33844cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 33944cd7ae7SLois Curfman McInnes if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0) 34044cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 34144cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 34244cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 34344cd7ae7SLois Curfman McInnes #else 34444cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 34544cd7ae7SLois Curfman McInnes #endif 34644cd7ae7SLois Curfman McInnes } 34744cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 34844cd7ae7SLois Curfman McInnes } 34944cd7ae7SLois Curfman McInnes } 35017ab2063SBarry Smith else { 35117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 35217ab2063SBarry Smith fprintf(fd,"row %d:",i); 353416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 35417ab2063SBarry Smith #if defined(PETSC_COMPLEX) 355416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 356416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 35717ab2063SBarry Smith } 35817ab2063SBarry Smith else { 359416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 36017ab2063SBarry Smith } 36117ab2063SBarry Smith #else 362416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 36317ab2063SBarry Smith #endif 36417ab2063SBarry Smith } 36517ab2063SBarry Smith fprintf(fd,"\n"); 36617ab2063SBarry Smith } 36717ab2063SBarry Smith } 36817ab2063SBarry Smith fflush(fd); 369416022c9SBarry Smith return 0; 370416022c9SBarry Smith } 371416022c9SBarry Smith 372416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 373416022c9SBarry Smith { 374416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 375cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 376cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 377bcd2baecSBarry Smith Draw draw; 378cddf8d76SBarry Smith DrawButton button; 37919bcc07fSBarry Smith PetscTruth isnull; 380cddf8d76SBarry Smith 381bcd2baecSBarry Smith ViewerDrawGetDraw(viewer,&draw); 38219bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 38319bcc07fSBarry Smith 384416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 385416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 386416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 387416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 388cddf8d76SBarry Smith color = DRAW_BLUE; 389416022c9SBarry Smith for ( i=0; i<m; i++ ) { 390cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 391416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 392cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 393cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 394cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 395cddf8d76SBarry Smith #else 396cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 397cddf8d76SBarry Smith #endif 398cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 399cddf8d76SBarry Smith } 400cddf8d76SBarry Smith } 401cddf8d76SBarry Smith color = DRAW_CYAN; 402cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 403cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 404cddf8d76SBarry 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 (a->a[j] != 0.) continue; 407cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 408cddf8d76SBarry Smith } 409cddf8d76SBarry Smith } 410cddf8d76SBarry Smith color = DRAW_RED; 411cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 412cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 413cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 414cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 415cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 416cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 417cddf8d76SBarry Smith #else 418cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 419cddf8d76SBarry Smith #endif 420cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 421416022c9SBarry Smith } 422416022c9SBarry Smith } 423416022c9SBarry Smith DrawFlush(draw); 424cddf8d76SBarry Smith DrawGetPause(draw,&pause); 425cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 426cddf8d76SBarry Smith 427cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 4286945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 429cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 430cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 431cddf8d76SBarry Smith DrawClear(draw); 432cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 433cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 434cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 435cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 436cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 437cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 438cddf8d76SBarry Smith w *= scale; h *= scale; 439cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 440cddf8d76SBarry Smith color = DRAW_BLUE; 441cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 442cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 443cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 444cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 445cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 446cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 447cddf8d76SBarry Smith #else 448cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 449cddf8d76SBarry Smith #endif 450cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 451cddf8d76SBarry Smith } 452cddf8d76SBarry Smith } 453cddf8d76SBarry Smith color = DRAW_CYAN; 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 (a->a[j] != 0.) continue; 459cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 460cddf8d76SBarry Smith } 461cddf8d76SBarry Smith } 462cddf8d76SBarry Smith color = DRAW_RED; 463cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 464cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 465cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 466cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 467cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 468cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 469cddf8d76SBarry Smith #else 470cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 471cddf8d76SBarry Smith #endif 472cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 473cddf8d76SBarry Smith } 474cddf8d76SBarry Smith } 4756945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 476cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 477cddf8d76SBarry Smith } 478416022c9SBarry Smith return 0; 479416022c9SBarry Smith } 480416022c9SBarry Smith 481416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 482416022c9SBarry Smith { 483416022c9SBarry Smith Mat A = (Mat) obj; 484416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 485bcd2baecSBarry Smith ViewerType vtype; 486bcd2baecSBarry Smith int ierr; 487416022c9SBarry Smith 488bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 489bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 490416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 491416022c9SBarry Smith } 492bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 493416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 494416022c9SBarry Smith } 495bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 496416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 497416022c9SBarry Smith } 498bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 499bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 50017ab2063SBarry Smith } 50117ab2063SBarry Smith return 0; 50217ab2063SBarry Smith } 50319bcc07fSBarry Smith 504c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 505416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 50617ab2063SBarry Smith { 507416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 50841c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 509416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 510416022c9SBarry Smith Scalar *aa = a->a, *ap; 51117ab2063SBarry Smith 5126d4a8577SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) return 0; 51317ab2063SBarry Smith 51417ab2063SBarry Smith for ( i=1; i<m; i++ ) { 515416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 51617ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 51717ab2063SBarry Smith if (fshift) { 518416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 51917ab2063SBarry Smith N = ailen[i]; 52017ab2063SBarry Smith for ( j=0; j<N; j++ ) { 52117ab2063SBarry Smith ip[j-fshift] = ip[j]; 52217ab2063SBarry Smith ap[j-fshift] = ap[j]; 52317ab2063SBarry Smith } 52417ab2063SBarry Smith } 52517ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 52617ab2063SBarry Smith } 52717ab2063SBarry Smith if (m) { 52817ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 52917ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 53017ab2063SBarry Smith } 53117ab2063SBarry Smith /* reset ilen and imax for each row */ 53217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 53317ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 53417ab2063SBarry Smith } 535416022c9SBarry Smith a->nz = ai[m] + shift; 53617ab2063SBarry Smith 53717ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 538416022c9SBarry Smith if (fshift && a->diag) { 5390452661fSBarry Smith PetscFree(a->diag); 540416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 541416022c9SBarry Smith a->diag = 0; 54217ab2063SBarry Smith } 5434e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 5444e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 5454e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 546b810aeb4SBarry Smith a->reallocs); 547*dd5f02e7SSatish Balay a->reallocs = 0; 5484e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 5494e220ebcSLois Curfman McInnes 55076dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 55141c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 55217ab2063SBarry Smith return 0; 55317ab2063SBarry Smith } 55417ab2063SBarry Smith 555416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 55617ab2063SBarry Smith { 557416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 558cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 55917ab2063SBarry Smith return 0; 56017ab2063SBarry Smith } 561416022c9SBarry Smith 56217ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 56317ab2063SBarry Smith { 564416022c9SBarry Smith Mat A = (Mat) obj; 565416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 56690f02eecSBarry Smith int ierr; 567d5d45c9bSBarry Smith 56817ab2063SBarry Smith #if defined(PETSC_LOG) 569416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 57017ab2063SBarry Smith #endif 5710452661fSBarry Smith PetscFree(a->a); 5720452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 5730452661fSBarry Smith if (a->diag) PetscFree(a->diag); 5740452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 5750452661fSBarry Smith if (a->imax) PetscFree(a->imax); 5760452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 57776dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 5780452661fSBarry Smith PetscFree(a); 57990f02eecSBarry Smith if (A->mapping) { 58090f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 58190f02eecSBarry Smith } 582f2655603SLois Curfman McInnes PLogObjectDestroy(A); 583f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 58417ab2063SBarry Smith return 0; 58517ab2063SBarry Smith } 58617ab2063SBarry Smith 587416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 58817ab2063SBarry Smith { 58917ab2063SBarry Smith return 0; 59017ab2063SBarry Smith } 59117ab2063SBarry Smith 592416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 59317ab2063SBarry Smith { 594416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 5956d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 5966d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 5976d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 598219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 5996d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 6006d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 6016d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 602219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 6036d4a8577SBarry Smith op == MAT_SYMMETRIC || 6046d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 60590f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 60690f02eecSBarry Smith op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES) 60794a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 6086d4a8577SBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) 6096d4a8577SBarry Smith {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:MAT_NO_NEW_DIAGONALS");} 6106d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 6116d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 6126d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 6136d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 6146d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 615e2f28af5SBarry Smith else 6164d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 61717ab2063SBarry Smith return 0; 61817ab2063SBarry Smith } 61917ab2063SBarry Smith 620416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 62117ab2063SBarry Smith { 622416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 623416022c9SBarry Smith int i,j, n,shift = a->indexshift; 62417ab2063SBarry Smith Scalar *x, zero = 0.0; 62517ab2063SBarry Smith 62617ab2063SBarry Smith VecSet(&zero,v); 62790f02eecSBarry Smith VecGetArray_Fast(v,x); VecGetLocalSize(v,&n); 628416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 629416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 630416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 631416022c9SBarry Smith if (a->j[j]+shift == i) { 632416022c9SBarry Smith x[i] = a->a[j]; 63317ab2063SBarry Smith break; 63417ab2063SBarry Smith } 63517ab2063SBarry Smith } 63617ab2063SBarry Smith } 63717ab2063SBarry Smith return 0; 63817ab2063SBarry Smith } 63917ab2063SBarry Smith 64017ab2063SBarry Smith /* -------------------------------------------------------*/ 64117ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 64217ab2063SBarry Smith /* -------------------------------------------------------*/ 64344cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 64417ab2063SBarry Smith { 645416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 64617ab2063SBarry Smith Scalar *x, *y, *v, alpha; 647416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 64817ab2063SBarry Smith 64990f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 650cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 65117ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 65217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 653416022c9SBarry Smith idx = a->j + a->i[i] + shift; 654416022c9SBarry Smith v = a->a + a->i[i] + shift; 655416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 65617ab2063SBarry Smith alpha = x[i]; 65717ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 65817ab2063SBarry Smith } 659416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 66017ab2063SBarry Smith return 0; 66117ab2063SBarry Smith } 662d5d45c9bSBarry Smith 66344cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 66417ab2063SBarry Smith { 665416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 66617ab2063SBarry Smith Scalar *x, *y, *v, alpha; 667416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 66817ab2063SBarry Smith 66990f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 67017ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 67117ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 67217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 673416022c9SBarry Smith idx = a->j + a->i[i] + shift; 674416022c9SBarry Smith v = a->a + a->i[i] + shift; 675416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 67617ab2063SBarry Smith alpha = x[i]; 67717ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 67817ab2063SBarry Smith } 67990f02eecSBarry Smith PLogFlops(2*a->nz); 68017ab2063SBarry Smith return 0; 68117ab2063SBarry Smith } 68217ab2063SBarry Smith 68344cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 68417ab2063SBarry Smith { 685416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 68617ab2063SBarry Smith Scalar *x, *y, *v, sum; 687416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 68817ab2063SBarry Smith 68990f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 69017ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 691416022c9SBarry Smith idx = a->j; 692416022c9SBarry Smith v = a->a; 693416022c9SBarry Smith ii = a->i; 69417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 695416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 69617ab2063SBarry Smith sum = 0.0; 69717ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 69817ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 699416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 70017ab2063SBarry Smith y[i] = sum; 70117ab2063SBarry Smith } 702416022c9SBarry Smith PLogFlops(2*a->nz - m); 70317ab2063SBarry Smith return 0; 70417ab2063SBarry Smith } 70517ab2063SBarry Smith 70644cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 70717ab2063SBarry Smith { 708416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 70917ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 710cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 71117ab2063SBarry Smith 71290f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z); 71317ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 714cddf8d76SBarry Smith idx = a->j; 715cddf8d76SBarry Smith v = a->a; 716cddf8d76SBarry Smith ii = a->i; 71717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 718cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 71917ab2063SBarry Smith sum = y[i]; 720cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 721cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 72217ab2063SBarry Smith z[i] = sum; 72317ab2063SBarry Smith } 724416022c9SBarry Smith PLogFlops(2*a->nz); 72517ab2063SBarry Smith return 0; 72617ab2063SBarry Smith } 72717ab2063SBarry Smith 72817ab2063SBarry Smith /* 72917ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 73017ab2063SBarry Smith */ 73117ab2063SBarry Smith 732416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 73317ab2063SBarry Smith { 734416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 735416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 73617ab2063SBarry Smith 7370452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 738416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 739416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 740416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 741416022c9SBarry Smith if (a->j[j]+shift == i) { 74217ab2063SBarry Smith diag[i] = j - shift; 74317ab2063SBarry Smith break; 74417ab2063SBarry Smith } 74517ab2063SBarry Smith } 74617ab2063SBarry Smith } 747416022c9SBarry Smith a->diag = diag; 74817ab2063SBarry Smith return 0; 74917ab2063SBarry Smith } 75017ab2063SBarry Smith 75144cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 75217ab2063SBarry Smith double fshift,int its,Vec xx) 75317ab2063SBarry Smith { 754416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 755416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 756d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 75717ab2063SBarry Smith 75890f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b); 759416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 760416022c9SBarry Smith diag = a->diag; 761416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 76217ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 76317ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 76417ab2063SBarry Smith bs = b + shift; 76517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 766416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 767416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 768416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 769416022c9SBarry Smith v = a->a + diag[i] + (!shift); 77017ab2063SBarry Smith sum = b[i]*d/omega; 77117ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 77217ab2063SBarry Smith x[i] = sum; 77317ab2063SBarry Smith } 77417ab2063SBarry Smith return 0; 77517ab2063SBarry Smith } 77617ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 77717ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 77817ab2063SBarry Smith } 779416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 78017ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 78117ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 78217ab2063SBarry Smith 78317ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 78417ab2063SBarry Smith 78517ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 78617ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 78717ab2063SBarry Smith is the relaxation factor. 78817ab2063SBarry Smith */ 7890452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 79017ab2063SBarry Smith scale = (2.0/omega) - 1.0; 79117ab2063SBarry Smith 79217ab2063SBarry Smith /* x = (E + U)^{-1} b */ 79317ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 794416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 795416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 796416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 797416022c9SBarry Smith v = a->a + diag[i] + (!shift); 79817ab2063SBarry Smith sum = b[i]; 79917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 80017ab2063SBarry Smith x[i] = omega*(sum/d); 80117ab2063SBarry Smith } 80217ab2063SBarry Smith 80317ab2063SBarry Smith /* t = b - (2*E - D)x */ 804416022c9SBarry Smith v = a->a; 80517ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 80617ab2063SBarry Smith 80717ab2063SBarry Smith /* t = (E + L)^{-1}t */ 808416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 809416022c9SBarry Smith diag = a->diag; 81017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 811416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 812416022c9SBarry Smith n = diag[i] - a->i[i]; 813416022c9SBarry Smith idx = a->j + a->i[i] + shift; 814416022c9SBarry Smith v = a->a + a->i[i] + shift; 81517ab2063SBarry Smith sum = t[i]; 81617ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 81717ab2063SBarry Smith t[i] = omega*(sum/d); 81817ab2063SBarry Smith } 81917ab2063SBarry Smith 82017ab2063SBarry Smith /* x = x + t */ 82117ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 8220452661fSBarry Smith PetscFree(t); 82317ab2063SBarry Smith return 0; 82417ab2063SBarry Smith } 82517ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 82617ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 82717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 828416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 829416022c9SBarry Smith n = diag[i] - a->i[i]; 830416022c9SBarry Smith idx = a->j + a->i[i] + shift; 831416022c9SBarry Smith v = a->a + a->i[i] + shift; 83217ab2063SBarry Smith sum = b[i]; 83317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 83417ab2063SBarry Smith x[i] = omega*(sum/d); 83517ab2063SBarry Smith } 83617ab2063SBarry Smith xb = x; 83717ab2063SBarry Smith } 83817ab2063SBarry Smith else xb = b; 83917ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 84017ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 84117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 842416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 84317ab2063SBarry Smith } 84417ab2063SBarry Smith } 84517ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 84617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 847416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 848416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 849416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 850416022c9SBarry Smith v = a->a + diag[i] + (!shift); 85117ab2063SBarry Smith sum = xb[i]; 85217ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 85317ab2063SBarry Smith x[i] = omega*(sum/d); 85417ab2063SBarry Smith } 85517ab2063SBarry Smith } 85617ab2063SBarry Smith its--; 85717ab2063SBarry Smith } 85817ab2063SBarry Smith while (its--) { 85917ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 86017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 861416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 862416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 863416022c9SBarry Smith idx = a->j + a->i[i] + shift; 864416022c9SBarry Smith v = a->a + a->i[i] + shift; 86517ab2063SBarry Smith sum = b[i]; 86617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 8677e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 86817ab2063SBarry Smith } 86917ab2063SBarry Smith } 87017ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 87117ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 872416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 873416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 874416022c9SBarry Smith idx = a->j + a->i[i] + shift; 875416022c9SBarry Smith v = a->a + a->i[i] + shift; 87617ab2063SBarry Smith sum = b[i]; 87717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 8787e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 87917ab2063SBarry Smith } 88017ab2063SBarry Smith } 88117ab2063SBarry Smith } 88217ab2063SBarry Smith return 0; 88317ab2063SBarry Smith } 88417ab2063SBarry Smith 8854e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 88617ab2063SBarry Smith { 887416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 8884e220ebcSLois Curfman McInnes 8894e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 8904e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 8914e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 8924e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 8934e220ebcSLois Curfman McInnes info->block_size = 1.0; 8944e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 8954e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 8964e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 8974e220ebcSLois Curfman McInnes /* if (info->nz_unneeded != A->info.nz_unneeded) 8984e220ebcSLois Curfman McInnes printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 8994e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 9004e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 9014e220ebcSLois Curfman McInnes info->memory = A->mem; 9024e220ebcSLois Curfman McInnes if (A->factor) { 9034e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 9044e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 9054e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 9064e220ebcSLois Curfman McInnes } else { 9074e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 9084e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 9094e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 9104e220ebcSLois Curfman McInnes } 91117ab2063SBarry Smith return 0; 91217ab2063SBarry Smith } 91317ab2063SBarry Smith 91417ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 91517ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 91617ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 91717ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 91817ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 91917ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 92017ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 92117ab2063SBarry Smith 92217ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 92317ab2063SBarry Smith { 924416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 925416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 92617ab2063SBarry Smith 92777c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 92817ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 92917ab2063SBarry Smith if (diag) { 93017ab2063SBarry Smith for ( i=0; i<N; i++ ) { 931416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 932416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 933416022c9SBarry Smith a->ilen[rows[i]] = 1; 934416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 935416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 93617ab2063SBarry Smith } 93717ab2063SBarry Smith else { 93817ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 93917ab2063SBarry Smith CHKERRQ(ierr); 94017ab2063SBarry Smith } 94117ab2063SBarry Smith } 94217ab2063SBarry Smith } 94317ab2063SBarry Smith else { 94417ab2063SBarry Smith for ( i=0; i<N; i++ ) { 945416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 946416022c9SBarry Smith a->ilen[rows[i]] = 0; 94717ab2063SBarry Smith } 94817ab2063SBarry Smith } 94917ab2063SBarry Smith ISRestoreIndices(is,&rows); 95043a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 95117ab2063SBarry Smith return 0; 95217ab2063SBarry Smith } 95317ab2063SBarry Smith 954416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 95517ab2063SBarry Smith { 956416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 957416022c9SBarry Smith *m = a->m; *n = a->n; 95817ab2063SBarry Smith return 0; 95917ab2063SBarry Smith } 96017ab2063SBarry Smith 961416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 96217ab2063SBarry Smith { 963416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 964416022c9SBarry Smith *m = 0; *n = a->m; 96517ab2063SBarry Smith return 0; 96617ab2063SBarry Smith } 9674e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 96817ab2063SBarry Smith { 969416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 970c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 97117ab2063SBarry Smith 972416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 97317ab2063SBarry Smith 974416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 975416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 97617ab2063SBarry Smith if (idx) { 977416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 9784e093b46SBarry Smith if (*nz && shift) { 9790452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 98017ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 9814e093b46SBarry Smith } else if (*nz) { 9824e093b46SBarry Smith *idx = itmp; 98317ab2063SBarry Smith } 98417ab2063SBarry Smith else *idx = 0; 98517ab2063SBarry Smith } 98617ab2063SBarry Smith return 0; 98717ab2063SBarry Smith } 98817ab2063SBarry Smith 9894e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 99017ab2063SBarry Smith { 9914e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 9924e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 99317ab2063SBarry Smith return 0; 99417ab2063SBarry Smith } 99517ab2063SBarry Smith 996cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 99717ab2063SBarry Smith { 998416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 999416022c9SBarry Smith Scalar *v = a->a; 100017ab2063SBarry Smith double sum = 0.0; 1001416022c9SBarry Smith int i, j,shift = a->indexshift; 100217ab2063SBarry Smith 100317ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1004416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 100517ab2063SBarry Smith #if defined(PETSC_COMPLEX) 100617ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 100717ab2063SBarry Smith #else 100817ab2063SBarry Smith sum += (*v)*(*v); v++; 100917ab2063SBarry Smith #endif 101017ab2063SBarry Smith } 101117ab2063SBarry Smith *norm = sqrt(sum); 101217ab2063SBarry Smith } 101317ab2063SBarry Smith else if (type == NORM_1) { 101417ab2063SBarry Smith double *tmp; 1015416022c9SBarry Smith int *jj = a->j; 10160452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 1017cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 101817ab2063SBarry Smith *norm = 0.0; 1019416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1020a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 102117ab2063SBarry Smith } 1022416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 102317ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 102417ab2063SBarry Smith } 10250452661fSBarry Smith PetscFree(tmp); 102617ab2063SBarry Smith } 102717ab2063SBarry Smith else if (type == NORM_INFINITY) { 102817ab2063SBarry Smith *norm = 0.0; 1029416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1030416022c9SBarry Smith v = a->a + a->i[j] + shift; 103117ab2063SBarry Smith sum = 0.0; 1032416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1033cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 103417ab2063SBarry Smith } 103517ab2063SBarry Smith if (sum > *norm) *norm = sum; 103617ab2063SBarry Smith } 103717ab2063SBarry Smith } 103817ab2063SBarry Smith else { 103948d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 104017ab2063SBarry Smith } 104117ab2063SBarry Smith return 0; 104217ab2063SBarry Smith } 104317ab2063SBarry Smith 1044416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 104517ab2063SBarry Smith { 1046416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1047416022c9SBarry Smith Mat C; 1048416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1049416022c9SBarry Smith int shift = a->indexshift; 1050d5d45c9bSBarry Smith Scalar *array = a->a; 105117ab2063SBarry Smith 10523638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 1053dfa27b74SSatish Balay SETERRQ(1,"MatTranspose_SeqAIJ:Square matrix only for in-place"); 10540452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1055cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 105617ab2063SBarry Smith if (shift) { 105717ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 105817ab2063SBarry Smith } 105917ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1060416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 10610452661fSBarry Smith PetscFree(col); 106217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 106317ab2063SBarry Smith len = ai[i+1]-ai[i]; 1064416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 106517ab2063SBarry Smith array += len; aj += len; 106617ab2063SBarry Smith } 106717ab2063SBarry Smith if (shift) { 106817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 106917ab2063SBarry Smith } 107017ab2063SBarry Smith 10716d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10726d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 107317ab2063SBarry Smith 10743638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1075416022c9SBarry Smith *B = C; 107617ab2063SBarry Smith } else { 1077416022c9SBarry Smith /* This isn't really an in-place transpose */ 10780452661fSBarry Smith PetscFree(a->a); 10790452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 10800452661fSBarry Smith if (a->diag) PetscFree(a->diag); 10810452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 10820452661fSBarry Smith if (a->imax) PetscFree(a->imax); 10830452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 10841073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 10850452661fSBarry Smith PetscFree(a); 1086416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 10870452661fSBarry Smith PetscHeaderDestroy(C); 108817ab2063SBarry Smith } 108917ab2063SBarry Smith return 0; 109017ab2063SBarry Smith } 109117ab2063SBarry Smith 1092f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 109317ab2063SBarry Smith { 1094416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 109517ab2063SBarry Smith Scalar *l,*r,x,*v; 1096d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 109717ab2063SBarry Smith 109817ab2063SBarry Smith if (ll) { 10993ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 11003ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 11019b1297e1SSatish Balay VecGetLocalSize_Fast(ll,m); 1102f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 110390f02eecSBarry Smith VecGetArray_Fast(ll,l); 1104416022c9SBarry Smith v = a->a; 110517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 110617ab2063SBarry Smith x = l[i]; 1107416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 110817ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 110917ab2063SBarry Smith } 111044cd7ae7SLois Curfman McInnes PLogFlops(nz); 111117ab2063SBarry Smith } 111217ab2063SBarry Smith if (rr) { 11139b1297e1SSatish Balay VecGetLocalSize_Fast(rr,n); 1114f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 111590f02eecSBarry Smith VecGetArray_Fast(rr,r); 1116416022c9SBarry Smith v = a->a; jj = a->j; 111717ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 111817ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 111917ab2063SBarry Smith } 112044cd7ae7SLois Curfman McInnes PLogFlops(nz); 112117ab2063SBarry Smith } 112217ab2063SBarry Smith return 0; 112317ab2063SBarry Smith } 112417ab2063SBarry Smith 1125cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 112617ab2063SBarry Smith { 1127db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 112802834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 112999141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1130a2744918SBarry Smith register int sum,lensi; 113199141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 113299141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 113399141d43SSatish Balay Scalar *a_new,*mat_a; 1134416022c9SBarry Smith Mat C; 113517ab2063SBarry Smith 1136b48a1e75SSatish Balay ierr = ISSorted(isrow,(PetscTruth*)&i); 1137b48a1e75SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:ISrow is not sorted"); 113899141d43SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i); 1139b48a1e75SSatish Balay if (!i) SETERRQ(1,"MatGetSubmatrices_SeqAIJ:IScol is not sorted"); 114099141d43SSatish Balay 114117ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 114217ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 114317ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 114417ab2063SBarry Smith 11457264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 114602834360SBarry Smith /* special case of contiguous rows */ 114757faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 114802834360SBarry Smith starts = lens + ncols; 114902834360SBarry Smith /* loop over new rows determining lens and starting points */ 115002834360SBarry Smith for (i=0; i<nrows; i++) { 1151a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1152a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 115302834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1154d8ced48eSBarry Smith if (aj[k]+shift >= first) { 115502834360SBarry Smith starts[i] = k; 115602834360SBarry Smith break; 115702834360SBarry Smith } 115802834360SBarry Smith } 1159a2744918SBarry Smith sum = 0; 116002834360SBarry Smith while (k < kend) { 1161d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1162a2744918SBarry Smith sum++; 116302834360SBarry Smith } 1164a2744918SBarry Smith lens[i] = sum; 116502834360SBarry Smith } 116602834360SBarry Smith /* create submatrix */ 1167cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 116808480c60SBarry Smith int n_cols,n_rows; 116908480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 117008480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1171d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 117208480c60SBarry Smith C = *B; 117308480c60SBarry Smith } 117408480c60SBarry Smith else { 117502834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 117608480c60SBarry Smith } 1177db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1178db02288aSLois Curfman McInnes 117902834360SBarry Smith /* loop over rows inserting into submatrix */ 1180db02288aSLois Curfman McInnes a_new = c->a; 1181db02288aSLois Curfman McInnes j_new = c->j; 1182db02288aSLois Curfman McInnes i_new = c->i; 1183db02288aSLois Curfman McInnes i_new[0] = -shift; 118402834360SBarry Smith for (i=0; i<nrows; i++) { 1185a2744918SBarry Smith ii = starts[i]; 1186a2744918SBarry Smith lensi = lens[i]; 1187a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1188a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 118902834360SBarry Smith } 1190a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1191a2744918SBarry Smith a_new += lensi; 1192a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1193a2744918SBarry Smith c->ilen[i] = lensi; 119402834360SBarry Smith } 11950452661fSBarry Smith PetscFree(lens); 119602834360SBarry Smith } 119702834360SBarry Smith else { 119802834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 11990452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 120002834360SBarry Smith ssmap = smap + shift; 120199141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1202cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 120317ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 120402834360SBarry Smith /* determine lens of each row */ 120502834360SBarry Smith for (i=0; i<nrows; i++) { 1206d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 120702834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 120802834360SBarry Smith lens[i] = 0; 120902834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1210d8ced48eSBarry Smith if (ssmap[aj[k]]) { 121102834360SBarry Smith lens[i]++; 121202834360SBarry Smith } 121302834360SBarry Smith } 121402834360SBarry Smith } 121517ab2063SBarry Smith /* Create and fill new matrix */ 1216a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 121799141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 121899141d43SSatish Balay 121999141d43SSatish Balay if (c->m != nrows || c->n != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 122099141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 122199141d43SSatish Balay SETERRQ(1,"MatGetSubmatrices_SeqAIJ:Cannot reuse matrix. wrong no of nonzeros"); 122299141d43SSatish Balay } 122399141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 122408480c60SBarry Smith C = *B; 122508480c60SBarry Smith } 122608480c60SBarry Smith else { 122702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 122808480c60SBarry Smith } 122999141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 123017ab2063SBarry Smith for (i=0; i<nrows; i++) { 123199141d43SSatish Balay row = irow[i]; 123217ab2063SBarry Smith nznew = 0; 123399141d43SSatish Balay kstart = ai[row]+shift; 123499141d43SSatish Balay kend = kstart + a->ilen[row]; 123599141d43SSatish Balay mat_i = c->i[i]+shift; 123699141d43SSatish Balay mat_j = c->j + mat_i; 123799141d43SSatish Balay mat_a = c->a + mat_i; 123899141d43SSatish Balay mat_ilen = c->ilen + i; 123917ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 124099141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 124199141d43SSatish Balay *mat_j++ = tcol - (!shift); 124299141d43SSatish Balay *mat_a++ = a->a[k]; 124399141d43SSatish Balay (*mat_ilen)++; 124499141d43SSatish Balay 124517ab2063SBarry Smith } 124617ab2063SBarry Smith } 124717ab2063SBarry Smith } 124802834360SBarry Smith /* Free work space */ 124902834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 125099141d43SSatish Balay PetscFree(smap); PetscFree(lens); 125102834360SBarry Smith } 12526d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 12536d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 125417ab2063SBarry Smith 125517ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1256416022c9SBarry Smith *B = C; 125717ab2063SBarry Smith return 0; 125817ab2063SBarry Smith } 125917ab2063SBarry Smith 1260a871dcd8SBarry Smith /* 126163b91edcSBarry Smith note: This can only work for identity for row and col. It would 126263b91edcSBarry Smith be good to check this and otherwise generate an error. 1263a871dcd8SBarry Smith */ 126463b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1265a871dcd8SBarry Smith { 126663b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 126708480c60SBarry Smith int ierr; 126863b91edcSBarry Smith Mat outA; 126963b91edcSBarry Smith 1270a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1271a871dcd8SBarry Smith 127263b91edcSBarry Smith outA = inA; 127363b91edcSBarry Smith inA->factor = FACTOR_LU; 127463b91edcSBarry Smith a->row = row; 127563b91edcSBarry Smith a->col = col; 127663b91edcSBarry Smith 12770452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 127863b91edcSBarry Smith 127908480c60SBarry Smith if (!a->diag) { 128008480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 128163b91edcSBarry Smith } 128263b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1283a871dcd8SBarry Smith return 0; 1284a871dcd8SBarry Smith } 1285a871dcd8SBarry Smith 1286f0b747eeSBarry Smith #include "pinclude/plapack.h" 1287f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1288f0b747eeSBarry Smith { 1289f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1290f0b747eeSBarry Smith int one = 1; 1291f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1292f0b747eeSBarry Smith PLogFlops(a->nz); 1293f0b747eeSBarry Smith return 0; 1294f0b747eeSBarry Smith } 1295f0b747eeSBarry Smith 1296cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1297cddf8d76SBarry Smith Mat **B) 1298cddf8d76SBarry Smith { 1299cddf8d76SBarry Smith int ierr,i; 1300cddf8d76SBarry Smith 1301cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 13020452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1303cddf8d76SBarry Smith } 1304cddf8d76SBarry Smith 1305cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1306905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1307cddf8d76SBarry Smith } 1308cddf8d76SBarry Smith return 0; 1309cddf8d76SBarry Smith } 1310cddf8d76SBarry Smith 13115a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 13125a838052SSatish Balay { 13135a838052SSatish Balay *bs = 1; 13145a838052SSatish Balay return 0; 13155a838052SSatish Balay } 13165a838052SSatish Balay 1317e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 13184dcbc457SBarry Smith { 1319e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 132006763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 13218a047759SSatish Balay int start, end, *ai, *aj; 132206763907SSatish Balay char *table; 13238a047759SSatish Balay shift = a->indexshift; 1324e4d965acSSatish Balay m = a->m; 1325e4d965acSSatish Balay ai = a->i; 13268a047759SSatish Balay aj = a->j+shift; 13278a047759SSatish Balay 1328e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 132906763907SSatish Balay 133006763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 133106763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 133206763907SSatish Balay 1333e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1334e4d965acSSatish Balay /* Initialise the two local arrays */ 1335e4d965acSSatish Balay isz = 0; 133606763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1337e4d965acSSatish Balay 1338e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 13394dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 134077c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1341e4d965acSSatish Balay 1342e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1343e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 134406763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 13454dcbc457SBarry Smith } 134606763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 134706763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1348e4d965acSSatish Balay 134904a348a9SBarry Smith k = 0; 135004a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 135104a348a9SBarry Smith n = isz; 135206763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1353e4d965acSSatish Balay row = nidx[k]; 1354e4d965acSSatish Balay start = ai[row]; 1355e4d965acSSatish Balay end = ai[row+1]; 135604a348a9SBarry Smith for ( l = start; l<end ; l++){ 13578a047759SSatish Balay val = aj[l] + shift; 135806763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1359e4d965acSSatish Balay } 1360e4d965acSSatish Balay } 1361e4d965acSSatish Balay } 1362537820f0SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1363e4d965acSSatish Balay } 136404a348a9SBarry Smith PetscFree(table); 136506763907SSatish Balay PetscFree(nidx); 1366e4d965acSSatish Balay return 0; 13674dcbc457SBarry Smith } 136817ab2063SBarry Smith 1369682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1370682d7d0cSBarry Smith { 1371682d7d0cSBarry Smith static int called = 0; 1372682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1373682d7d0cSBarry Smith 1374682d7d0cSBarry Smith if (called) return 0; else called = 1; 137577c4ece6SBarry Smith PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 137677c4ece6SBarry Smith PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>\n"); 137777c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 137877c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_no_inode - Do not use inodes\n"); 137977c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1380682d7d0cSBarry Smith #if defined(HAVE_ESSL) 138177c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1382682d7d0cSBarry Smith #endif 1383682d7d0cSBarry Smith return 0; 1384682d7d0cSBarry Smith } 1385205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1386a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1387a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1388a93ec695SBarry Smith 1389682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 139017ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 139117ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1392416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1393416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 139417ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 139517ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 139617ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 139717ab2063SBarry Smith MatRelax_SeqAIJ, 139817ab2063SBarry Smith MatTranspose_SeqAIJ, 13997264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1400f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 140117ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 140217ab2063SBarry Smith MatCompress_SeqAIJ, 140317ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 140417ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 140517ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 140617ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 140717ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 14083d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1409cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 14107eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1411682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1412f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 14135a838052SSatish Balay MatScale_SeqAIJ,0,0, 14146945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 14156945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 14163b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 14173b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 14183b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1419a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1420a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 1421a93ec695SBarry Smith MatColoringPatch_SeqAIJ}; 142217ab2063SBarry Smith 142317ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 142417ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 142517ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 142617ab2063SBarry Smith 142717ab2063SBarry Smith /*@C 1428682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 14290d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 14306e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 14312bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 14322bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 143317ab2063SBarry Smith 143417ab2063SBarry Smith Input Parameters: 143517ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 143617ab2063SBarry Smith . m - number of rows 143717ab2063SBarry Smith . n - number of columns 143817ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 14392bd5e0b2SLois Curfman McInnes . nzz - array containing the number of nonzeros in the various rows 14402bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 144117ab2063SBarry Smith 144217ab2063SBarry Smith Output Parameter: 1443416022c9SBarry Smith . A - the matrix 144417ab2063SBarry Smith 144517ab2063SBarry Smith Notes: 144617ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 144717ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 14480002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 144944cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 145017ab2063SBarry Smith 145117ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1452a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 14533d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 14546da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 145517ab2063SBarry Smith 1456682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1457682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1458682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 14596c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 14606c7ebb05SLois Curfman McInnes 14616c7ebb05SLois Curfman McInnes Options Database Keys: 14626c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 14636e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 14646e62573dSLois Curfman McInnes $ (max limit=5) 14656e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 14666e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 14676e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 146817ab2063SBarry Smith 146917ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 147017ab2063SBarry Smith @*/ 1471416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 147217ab2063SBarry Smith { 1473416022c9SBarry Smith Mat B; 1474416022c9SBarry Smith Mat_SeqAIJ *b; 14756945ee14SBarry Smith int i, len, ierr, flg,size; 14766945ee14SBarry Smith 14776945ee14SBarry Smith MPI_Comm_size(comm,&size); 14786945ee14SBarry Smith if (size > 1) SETERRQ(1,"MatCreateSeqAIJ:Comm must be of size 1"); 1479d5d45c9bSBarry Smith 1480416022c9SBarry Smith *A = 0; 14810452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1482416022c9SBarry Smith PLogObjectCreate(B); 14830452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 148444cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1485416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1486416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1487416022c9SBarry Smith B->view = MatView_SeqAIJ; 1488416022c9SBarry Smith B->factor = 0; 1489416022c9SBarry Smith B->lupivotthreshold = 1.0; 149090f02eecSBarry Smith B->mapping = 0; 14917a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 149269957df2SSatish Balay &flg); CHKERRQ(ierr); 14937a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 14947a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 14957a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1496416022c9SBarry Smith b->row = 0; 1497416022c9SBarry Smith b->col = 0; 1498416022c9SBarry Smith b->indexshift = 0; 1499b810aeb4SBarry Smith b->reallocs = 0; 150069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 150169957df2SSatish Balay if (flg) b->indexshift = -1; 150217ab2063SBarry Smith 150344cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 150444cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 15050452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1506b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 15077b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 15087b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1509416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 151017ab2063SBarry Smith nz = nz*m; 151117ab2063SBarry Smith } 151217ab2063SBarry Smith else { 151317ab2063SBarry Smith nz = 0; 1514416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 151517ab2063SBarry Smith } 151617ab2063SBarry Smith 151717ab2063SBarry Smith /* allocate the matrix space */ 1518416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 15190452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1520416022c9SBarry Smith b->j = (int *) (b->a + nz); 1521cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1522416022c9SBarry Smith b->i = b->j + nz; 1523416022c9SBarry Smith b->singlemalloc = 1; 152417ab2063SBarry Smith 1525416022c9SBarry Smith b->i[0] = -b->indexshift; 152617ab2063SBarry Smith for (i=1; i<m+1; i++) { 1527416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 152817ab2063SBarry Smith } 152917ab2063SBarry Smith 1530416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 15310452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1532416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1533416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 153417ab2063SBarry Smith 1535416022c9SBarry Smith b->nz = 0; 1536416022c9SBarry Smith b->maxnz = nz; 1537416022c9SBarry Smith b->sorted = 0; 1538416022c9SBarry Smith b->roworiented = 1; 1539416022c9SBarry Smith b->nonew = 0; 1540416022c9SBarry Smith b->diag = 0; 1541416022c9SBarry Smith b->solve_work = 0; 154271bd300dSLois Curfman McInnes b->spptr = 0; 1543754ec7b1SSatish Balay b->inode.node_count = 0; 1544754ec7b1SSatish Balay b->inode.size = 0; 15456c7ebb05SLois Curfman McInnes b->inode.limit = 5; 15466c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 15474e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 154817ab2063SBarry Smith 1549416022c9SBarry Smith *A = B; 15504e220ebcSLois Curfman McInnes 15514b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 15524b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 155369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 155469957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 15554b14c69eSBarry Smith #endif 155669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 155769957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 155869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 155969957df2SSatish Balay if (flg) { 1560416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1561416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 156217ab2063SBarry Smith } 156369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 156469957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 156517ab2063SBarry Smith return 0; 156617ab2063SBarry Smith } 156717ab2063SBarry Smith 15683d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 156917ab2063SBarry Smith { 1570416022c9SBarry Smith Mat C; 1571416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 157208480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 157317ab2063SBarry Smith 15744043dd9cSLois Curfman McInnes *B = 0; 15750452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1576416022c9SBarry Smith PLogObjectCreate(C); 15770452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 157841c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1579416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1580416022c9SBarry Smith C->view = MatView_SeqAIJ; 1581416022c9SBarry Smith C->factor = A->factor; 1582416022c9SBarry Smith c->row = 0; 1583416022c9SBarry Smith c->col = 0; 1584416022c9SBarry Smith c->indexshift = shift; 1585c456f294SBarry Smith C->assembled = PETSC_TRUE; 158617ab2063SBarry Smith 158744cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 158844cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 158944cd7ae7SLois Curfman McInnes C->M = a->m; 159044cd7ae7SLois Curfman McInnes C->N = a->n; 159117ab2063SBarry Smith 15920452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 15930452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 159417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1595416022c9SBarry Smith c->imax[i] = a->imax[i]; 1596416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 159717ab2063SBarry Smith } 159817ab2063SBarry Smith 159917ab2063SBarry Smith /* allocate the matrix space */ 1600416022c9SBarry Smith c->singlemalloc = 1; 1601416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 16020452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1603416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1604416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1605416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 160617ab2063SBarry Smith if (m > 0) { 1607416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 160808480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1609416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 161017ab2063SBarry Smith } 161108480c60SBarry Smith } 161217ab2063SBarry Smith 1613416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1614416022c9SBarry Smith c->sorted = a->sorted; 1615416022c9SBarry Smith c->roworiented = a->roworiented; 1616416022c9SBarry Smith c->nonew = a->nonew; 16177a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 161817ab2063SBarry Smith 1619416022c9SBarry Smith if (a->diag) { 16200452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1621416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 162217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1623416022c9SBarry Smith c->diag[i] = a->diag[i]; 162417ab2063SBarry Smith } 162517ab2063SBarry Smith } 1626416022c9SBarry Smith else c->diag = 0; 16276c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 16286c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1629754ec7b1SSatish Balay if (a->inode.size){ 1630daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 1631754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1632daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 1633754ec7b1SSatish Balay } else { 1634754ec7b1SSatish Balay c->inode.size = 0; 1635754ec7b1SSatish Balay c->inode.node_count = 0; 1636754ec7b1SSatish Balay } 1637416022c9SBarry Smith c->nz = a->nz; 1638416022c9SBarry Smith c->maxnz = a->maxnz; 1639416022c9SBarry Smith c->solve_work = 0; 164076dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1641754ec7b1SSatish Balay 1642416022c9SBarry Smith *B = C; 164317ab2063SBarry Smith return 0; 164417ab2063SBarry Smith } 164517ab2063SBarry Smith 164619bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 164717ab2063SBarry Smith { 1648416022c9SBarry Smith Mat_SeqAIJ *a; 1649416022c9SBarry Smith Mat B; 165017699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1651bcd2baecSBarry Smith MPI_Comm comm; 165217ab2063SBarry Smith 165319bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 165417699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 165517699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 165690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 165777c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 165848d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 165917ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 166017ab2063SBarry Smith 166117ab2063SBarry Smith /* read in row lengths */ 16620452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 166377c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 166417ab2063SBarry Smith 166517ab2063SBarry Smith /* create our matrix */ 1666416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1667416022c9SBarry Smith B = *A; 1668416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1669416022c9SBarry Smith shift = a->indexshift; 167017ab2063SBarry Smith 167117ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 167277c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 167317ab2063SBarry Smith if (shift) { 167417ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1675416022c9SBarry Smith a->j[i] += 1; 167617ab2063SBarry Smith } 167717ab2063SBarry Smith } 167817ab2063SBarry Smith 167917ab2063SBarry Smith /* read in nonzero values */ 168077c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 168117ab2063SBarry Smith 168217ab2063SBarry Smith /* set matrix "i" values */ 1683416022c9SBarry Smith a->i[0] = -shift; 168417ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1685416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1686416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 168717ab2063SBarry Smith } 16880452661fSBarry Smith PetscFree(rowlengths); 168917ab2063SBarry Smith 16906d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 16916d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 169217ab2063SBarry Smith return 0; 169317ab2063SBarry Smith } 169417ab2063SBarry Smith 1695686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 16967264ac53SSatish Balay { 16977264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 16987264ac53SSatish Balay 1699bcd2baecSBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,"MatEqual_SeqAIJ:Matrices must be same type"); 17007264ac53SSatish Balay 17017264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 17027264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1703bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 170477c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1705bcd2baecSBarry Smith } 17067264ac53SSatish Balay 17077264ac53SSatish Balay /* if the a->i are the same */ 1708bcd2baecSBarry Smith if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 170977c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 17107264ac53SSatish Balay } 17117264ac53SSatish Balay 17127264ac53SSatish Balay /* if a->j are the same */ 1713bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 171477c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1715bcd2baecSBarry Smith } 1716bcd2baecSBarry Smith 1717bcd2baecSBarry Smith /* if a->a are the same */ 171819bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 171977c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 17207264ac53SSatish Balay } 172177c4ece6SBarry Smith *flg = PETSC_TRUE; 17227264ac53SSatish Balay return 0; 17237264ac53SSatish Balay 17247264ac53SSatish Balay } 1725