1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*480ef9eaSBarry Smith static char vcid[] = "$Id: aij.c,v 1.265 1998/05/13 14:17:23 bsmith Exp bsmith $"; 317ab2063SBarry Smith #endif 417ab2063SBarry Smith 5d5d45c9bSBarry Smith /* 63369ce9aSBarry Smith Defines the basic matrix operations for the AIJ (compressed row) 7d5d45c9bSBarry Smith matrix storage format. 8d5d45c9bSBarry Smith */ 93369ce9aSBarry Smith 103369ce9aSBarry Smith #include "pinclude/pviewer.h" 113369ce9aSBarry Smith #include "sys.h" 1270f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 13f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 14f5eb4b81SSatish Balay #include "src/inline/spops.h" 158d195f9aSBarry Smith #include "src/inline/dot.h" 16f5eb4b81SSatish Balay #include "src/inline/bitarray.h" 1717ab2063SBarry Smith 18a2ce50c7SBarry Smith /* 19a2ce50c7SBarry Smith Basic AIJ format ILU based on drop tolerance 20a2ce50c7SBarry Smith */ 215615d1e5SSatish Balay #undef __FUNC__ 225615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ" 23a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact) 24a2ce50c7SBarry Smith { 25a2ce50c7SBarry Smith /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */ 26a2ce50c7SBarry Smith int ierr = 1; 27a2ce50c7SBarry Smith 283a40ed3dSBarry Smith PetscFunctionBegin; 29e3372554SBarry Smith SETERRQ(ierr,0,"Not implemented"); 30d64ed03dSBarry Smith #if !defined(USE_PETSC_DEBUG) 31d64ed03dSBarry Smith PetscFunctionReturn(0); 32d64ed03dSBarry Smith #endif 33a2ce50c7SBarry Smith } 34a2ce50c7SBarry Smith 35bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 3617ab2063SBarry Smith 375615d1e5SSatish Balay #undef __FUNC__ 385615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ" 398f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja, 406945ee14SBarry Smith PetscTruth *done) 4117ab2063SBarry Smith { 42416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 436945ee14SBarry Smith int ierr,i,ishift; 4417ab2063SBarry Smith 453a40ed3dSBarry Smith PetscFunctionBegin; 4631625ec5SSatish Balay *m = A->m; 473a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 486945ee14SBarry Smith ishift = a->indexshift; 496945ee14SBarry Smith if (symmetric) { 5031625ec5SSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 516945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 5231625ec5SSatish Balay int nz = a->i[a->m]; 533b2fbd54SBarry Smith /* malloc space and subtract 1 from i and j indices */ 5431625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 553b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 563b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1; 5731625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1; 586945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 5931625ec5SSatish Balay int nz = a->i[a->m] + 1; 603b2fbd54SBarry Smith /* malloc space and add 1 to i and j indices */ 6131625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 623b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 633b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1; 6431625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1; 656945ee14SBarry Smith } else { 666945ee14SBarry Smith *ia = a->i; *ja = a->j; 67a2ce50c7SBarry Smith } 68a2ce50c7SBarry Smith 693a40ed3dSBarry Smith PetscFunctionReturn(0); 70a2744918SBarry Smith } 71a2744918SBarry Smith 725615d1e5SSatish Balay #undef __FUNC__ 735615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ" 748f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 756945ee14SBarry Smith PetscTruth *done) 766945ee14SBarry Smith { 776945ee14SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 783b2fbd54SBarry Smith int ishift = a->indexshift; 796945ee14SBarry Smith 803a40ed3dSBarry Smith PetscFunctionBegin; 813a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 823b2fbd54SBarry Smith if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) { 836945ee14SBarry Smith PetscFree(*ia); 846945ee14SBarry Smith PetscFree(*ja); 85bcd2baecSBarry Smith } 863a40ed3dSBarry Smith PetscFunctionReturn(0); 8717ab2063SBarry Smith } 8817ab2063SBarry Smith 895615d1e5SSatish Balay #undef __FUNC__ 905615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ" 9143a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 923b2fbd54SBarry Smith PetscTruth *done) 933b2fbd54SBarry Smith { 943b2fbd54SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 95a93ec695SBarry Smith int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m; 96a93ec695SBarry Smith int nz = a->i[m]+ishift,row,*jj,mr,col; 973b2fbd54SBarry Smith 983a40ed3dSBarry Smith PetscFunctionBegin; 993b2fbd54SBarry Smith *nn = A->n; 1003a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 1013b2fbd54SBarry Smith if (symmetric) { 102179192dfSSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 1033b2fbd54SBarry Smith } else { 10461d2ded1SBarry Smith collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths); 1053b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1063b2fbd54SBarry Smith cia = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia); 107a93ec695SBarry Smith cja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja); 1083b2fbd54SBarry Smith jj = a->j; 1093b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) { 1103b2fbd54SBarry Smith collengths[jj[i] + ishift]++; 1113b2fbd54SBarry Smith } 1123b2fbd54SBarry Smith cia[0] = oshift; 1133b2fbd54SBarry Smith for ( i=0; i<n; i++) { 1143b2fbd54SBarry Smith cia[i+1] = cia[i] + collengths[i]; 1153b2fbd54SBarry Smith } 1163b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1173b2fbd54SBarry Smith jj = a->j; 118a93ec695SBarry Smith for ( row=0; row<m; row++ ) { 119a93ec695SBarry Smith mr = a->i[row+1] - a->i[row]; 120a93ec695SBarry Smith for ( i=0; i<mr; i++ ) { 1213b2fbd54SBarry Smith col = *jj++ + ishift; 1223b2fbd54SBarry Smith cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1233b2fbd54SBarry Smith } 1243b2fbd54SBarry Smith } 1253b2fbd54SBarry Smith PetscFree(collengths); 1263b2fbd54SBarry Smith *ia = cia; *ja = cja; 1273b2fbd54SBarry Smith } 1283b2fbd54SBarry Smith 1293a40ed3dSBarry Smith PetscFunctionReturn(0); 1303b2fbd54SBarry Smith } 1313b2fbd54SBarry Smith 1325615d1e5SSatish Balay #undef __FUNC__ 1335615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ" 13443a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia, 1353b2fbd54SBarry Smith int **ja,PetscTruth *done) 1363b2fbd54SBarry Smith { 1373a40ed3dSBarry Smith PetscFunctionBegin; 1383a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 1393b2fbd54SBarry Smith 1403b2fbd54SBarry Smith PetscFree(*ia); 1413b2fbd54SBarry Smith PetscFree(*ja); 1423b2fbd54SBarry Smith 1433a40ed3dSBarry Smith PetscFunctionReturn(0); 1443b2fbd54SBarry Smith } 1453b2fbd54SBarry Smith 146227d817aSBarry Smith #define CHUNKSIZE 15 14717ab2063SBarry Smith 1485615d1e5SSatish Balay #undef __FUNC__ 1495615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ" 15061d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 15117ab2063SBarry Smith { 152416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 153416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 1544b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 155d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 156416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 15717ab2063SBarry Smith 1583a40ed3dSBarry Smith PetscFunctionBegin; 15917ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 160416022c9SBarry Smith row = im[k]; 1613a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 162a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 163a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 1643b2fbd54SBarry Smith #endif 16517ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 16617ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 167416022c9SBarry Smith low = 0; 16817ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 1693a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 170a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 171a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 1723b2fbd54SBarry Smith #endif 1734b0e389bSBarry Smith col = in[l] - shift; 1744b0e389bSBarry Smith if (roworiented) { 1754b0e389bSBarry Smith value = *v++; 176bef8e0ddSBarry Smith } else { 1774b0e389bSBarry Smith value = v[k + l*m]; 1784b0e389bSBarry Smith } 179416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 180416022c9SBarry Smith while (high-low > 5) { 181416022c9SBarry Smith t = (low+high)/2; 182416022c9SBarry Smith if (rp[t] > col) high = t; 183416022c9SBarry Smith else low = t; 18417ab2063SBarry Smith } 185416022c9SBarry Smith for ( i=low; i<high; i++ ) { 18617ab2063SBarry Smith if (rp[i] > col) break; 18717ab2063SBarry Smith if (rp[i] == col) { 188416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 18917ab2063SBarry Smith else ap[i] = value; 19017ab2063SBarry Smith goto noinsert; 19117ab2063SBarry Smith } 19217ab2063SBarry Smith } 193c2653b3dSLois Curfman McInnes if (nonew == 1) goto noinsert; 194a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 19517ab2063SBarry Smith if (nrow >= rmax) { 19617ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 197416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 19817ab2063SBarry Smith Scalar *new_a; 19917ab2063SBarry Smith 200a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 20196854ed6SLois Curfman McInnes 20217ab2063SBarry Smith /* malloc new storage space */ 203416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 2040452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 20517ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 20617ab2063SBarry Smith new_i = new_j + new_nz; 20717ab2063SBarry Smith 20817ab2063SBarry Smith /* copy over old data into new slots */ 20917ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 210416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 211416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 212416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 213416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 21417ab2063SBarry Smith len*sizeof(int)); 215416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 216416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 21717ab2063SBarry Smith len*sizeof(Scalar)); 21817ab2063SBarry Smith /* free up old matrix storage */ 2190452661fSBarry Smith PetscFree(a->a); 2200452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 221416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 222416022c9SBarry Smith a->singlemalloc = 1; 22317ab2063SBarry Smith 22417ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 225416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 226416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 227416022c9SBarry Smith a->maxnz += CHUNKSIZE; 228b810aeb4SBarry Smith a->reallocs++; 22917ab2063SBarry Smith } 230416022c9SBarry Smith N = nrow++ - 1; a->nz++; 231416022c9SBarry Smith /* shift up all the later entries in this row */ 232416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 23317ab2063SBarry Smith rp[ii+1] = rp[ii]; 23417ab2063SBarry Smith ap[ii+1] = ap[ii]; 23517ab2063SBarry Smith } 23617ab2063SBarry Smith rp[i] = col; 23717ab2063SBarry Smith ap[i] = value; 23817ab2063SBarry Smith noinsert:; 239416022c9SBarry Smith low = i + 1; 24017ab2063SBarry Smith } 24117ab2063SBarry Smith ailen[row] = nrow; 24217ab2063SBarry Smith } 2433a40ed3dSBarry Smith PetscFunctionReturn(0); 24417ab2063SBarry Smith } 24517ab2063SBarry Smith 2465615d1e5SSatish Balay #undef __FUNC__ 2475615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ" 2488f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 2497eb43aa7SLois Curfman McInnes { 2507eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 251b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 2527eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 2537eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 2547eb43aa7SLois Curfman McInnes 2553a40ed3dSBarry Smith PetscFunctionBegin; 2567eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 2577eb43aa7SLois Curfman McInnes row = im[k]; 258a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 259a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 2607eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 2617eb43aa7SLois Curfman McInnes nrow = ailen[row]; 2627eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 263a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 264a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 2657eb43aa7SLois Curfman McInnes col = in[l] - shift; 2667eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 2677eb43aa7SLois Curfman McInnes while (high-low > 5) { 2687eb43aa7SLois Curfman McInnes t = (low+high)/2; 2697eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 2707eb43aa7SLois Curfman McInnes else low = t; 2717eb43aa7SLois Curfman McInnes } 2727eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 2737eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2747eb43aa7SLois Curfman McInnes if (rp[i] == col) { 275b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2767eb43aa7SLois Curfman McInnes goto finished; 2777eb43aa7SLois Curfman McInnes } 2787eb43aa7SLois Curfman McInnes } 279b49de8d1SLois Curfman McInnes *v++ = zero; 2807eb43aa7SLois Curfman McInnes finished:; 2817eb43aa7SLois Curfman McInnes } 2827eb43aa7SLois Curfman McInnes } 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 2847eb43aa7SLois Curfman McInnes } 2857eb43aa7SLois Curfman McInnes 28617ab2063SBarry Smith 2875615d1e5SSatish Balay #undef __FUNC__ 2885615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary" 289*480ef9eaSBarry Smith int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 29017ab2063SBarry Smith { 291416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 292416022c9SBarry Smith int i, fd, *col_lens, ierr; 29317ab2063SBarry Smith 2943a40ed3dSBarry Smith PetscFunctionBegin; 29590ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2960452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 297416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 298416022c9SBarry Smith col_lens[1] = a->m; 299416022c9SBarry Smith col_lens[2] = a->n; 300416022c9SBarry Smith col_lens[3] = a->nz; 301416022c9SBarry Smith 302416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 303416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 304416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 30517ab2063SBarry Smith } 3060752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3070452661fSBarry Smith PetscFree(col_lens); 308416022c9SBarry Smith 309416022c9SBarry Smith /* store column indices (zero start index) */ 310416022c9SBarry Smith if (a->indexshift) { 311416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 31217ab2063SBarry Smith } 3130752156aSBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0); CHKERRQ(ierr); 314416022c9SBarry Smith if (a->indexshift) { 315416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 31617ab2063SBarry Smith } 317416022c9SBarry Smith 318416022c9SBarry Smith /* store nonzero values */ 3190752156aSBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 3203a40ed3dSBarry Smith PetscFunctionReturn(0); 32117ab2063SBarry Smith } 322416022c9SBarry Smith 3235615d1e5SSatish Balay #undef __FUNC__ 3245615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII" 325*480ef9eaSBarry Smith int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 326416022c9SBarry Smith { 327416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 328496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 32917ab2063SBarry Smith FILE *fd; 33017ab2063SBarry Smith char *outputname; 33117ab2063SBarry Smith 3323a40ed3dSBarry Smith PetscFunctionBegin; 33390ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 334416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 33590ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 336a93ec695SBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 3373a40ed3dSBarry Smith PetscFunctionReturn(0); 3383a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 339496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 340496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 341496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 34295e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 34395e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 3443a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 345d00d2cf4SBarry Smith int nofinalvalue = 0; 346d00d2cf4SBarry Smith if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) { 347d00d2cf4SBarry Smith nofinalvalue = 1; 348d00d2cf4SBarry Smith } 349416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3504e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 351d00d2cf4SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue); 35217ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 35317ab2063SBarry Smith 35417ab2063SBarry Smith for (i=0; i<m; i++) { 355416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 3563a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 357d64ed03dSBarry Smith fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]),imag(a->a[j])); 35817ab2063SBarry Smith #else 3597a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 36017ab2063SBarry Smith #endif 36117ab2063SBarry Smith } 36217ab2063SBarry Smith } 363d00d2cf4SBarry Smith if (nofinalvalue) { 364d00d2cf4SBarry Smith fprintf(fd,"%d %d %18.16e\n", m, a->n, 0.0); 365d00d2cf4SBarry Smith } 36617ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 3673a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 36844cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 36944cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 37044cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 3713a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 372766eeae4SLois Curfman McInnes if (imag(a->a[j]) > 0.0 && real(a->a[j]) != 0.0) 37344cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 374766eeae4SLois Curfman McInnes else if (imag(a->a[j]) < 0.0 && real(a->a[j]) != 0.0) 375766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j])); 37644cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 37744cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 37844cd7ae7SLois Curfman McInnes #else 37944cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 38044cd7ae7SLois Curfman McInnes #endif 38144cd7ae7SLois Curfman McInnes } 38244cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 38344cd7ae7SLois Curfman McInnes } 38444cd7ae7SLois Curfman McInnes } 385496be53dSLois Curfman McInnes else if (format == VIEWER_FORMAT_ASCII_SYMMODU) { 386496be53dSLois Curfman McInnes int nzd=0, fshift=1, *sptr; 3872e44a96cSLois Curfman McInnes sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr); 388496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 389496be53dSLois Curfman McInnes sptr[i] = nzd+1; 390496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 391496be53dSLois Curfman McInnes if (a->j[j] >= i) { 3923a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 393496be53dSLois Curfman McInnes if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) nzd++; 394496be53dSLois Curfman McInnes #else 395496be53dSLois Curfman McInnes if (a->a[j] != 0.0) nzd++; 396496be53dSLois Curfman McInnes #endif 397496be53dSLois Curfman McInnes } 398496be53dSLois Curfman McInnes } 399496be53dSLois Curfman McInnes } 4002e44a96cSLois Curfman McInnes sptr[m] = nzd+1; 401496be53dSLois Curfman McInnes fprintf(fd," %d %d\n\n",m,nzd); 4022e44a96cSLois Curfman McInnes for ( i=0; i<m+1; i+=6 ) { 4032e44a96cSLois Curfman McInnes if (i+4<m) fprintf(fd," %d %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4],sptr[i+5]); 4042e44a96cSLois Curfman McInnes else if (i+3<m) fprintf(fd," %d %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3],sptr[i+4]); 4052e44a96cSLois Curfman McInnes else if (i+2<m) fprintf(fd," %d %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2],sptr[i+3]); 4062e44a96cSLois Curfman McInnes else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]); 4072e44a96cSLois Curfman McInnes else if (i<m) fprintf(fd," %d %d\n",sptr[i],sptr[i+1]); 4087272d637SLois Curfman McInnes else fprintf(fd," %d\n",sptr[i]); 409496be53dSLois Curfman McInnes } 410496be53dSLois Curfman McInnes fprintf(fd,"\n"); 411496be53dSLois Curfman McInnes PetscFree(sptr); 412496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 413496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 414496be53dSLois Curfman McInnes if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift); 415496be53dSLois Curfman McInnes } 416496be53dSLois Curfman McInnes fprintf(fd,"\n"); 417496be53dSLois Curfman McInnes } 418496be53dSLois Curfman McInnes fprintf(fd,"\n"); 419496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 420496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 421496be53dSLois Curfman McInnes if (a->j[j] >= i) { 4223a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 423496be53dSLois Curfman McInnes if (imag(a->a[j]) != 0.0 || real(a->a[j]) != 0.0) 424496be53dSLois Curfman McInnes fprintf(fd," %18.16e %18.16e ",real(a->a[j]),imag(a->a[j])); 425496be53dSLois Curfman McInnes #else 426496be53dSLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]); 427496be53dSLois Curfman McInnes #endif 428496be53dSLois Curfman McInnes } 429496be53dSLois Curfman McInnes } 430496be53dSLois Curfman McInnes fprintf(fd,"\n"); 431496be53dSLois Curfman McInnes } 4323a40ed3dSBarry Smith } else { 43317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 43417ab2063SBarry Smith fprintf(fd,"row %d:",i); 435416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 4363a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 437766eeae4SLois Curfman McInnes if (imag(a->a[j]) > 0.0) { 438416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 439766eeae4SLois Curfman McInnes } else if (imag(a->a[j]) < 0.0) { 440766eeae4SLois Curfman McInnes fprintf(fd," %d %g - %g i",a->j[j]+shift,real(a->a[j]),-imag(a->a[j])); 4413a40ed3dSBarry Smith } else { 442416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 44317ab2063SBarry Smith } 44417ab2063SBarry Smith #else 445416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 44617ab2063SBarry Smith #endif 44717ab2063SBarry Smith } 44817ab2063SBarry Smith fprintf(fd,"\n"); 44917ab2063SBarry Smith } 45017ab2063SBarry Smith } 45117ab2063SBarry Smith fflush(fd); 4523a40ed3dSBarry Smith PetscFunctionReturn(0); 453416022c9SBarry Smith } 454416022c9SBarry Smith 4555615d1e5SSatish Balay #undef __FUNC__ 456*480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw_Zoom" 457*480ef9eaSBarry Smith int MatView_SeqAIJ_Draw_Zoom(Draw draw,void *Aa) 458416022c9SBarry Smith { 459*480ef9eaSBarry Smith Mat A = (Mat) Aa; 460416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 461*480ef9eaSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,color; 4620513a670SBarry Smith int format; 463*480ef9eaSBarry Smith double xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0; 464cddf8d76SBarry Smith DrawButton button; 46519bcc07fSBarry Smith PetscTruth isnull; 466*480ef9eaSBarry Smith Viewer viewer; 467cddf8d76SBarry Smith 4683a40ed3dSBarry Smith PetscFunctionBegin; 469*480ef9eaSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 4700513a670SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 47119bcc07fSBarry Smith 472*480ef9eaSBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr); 473416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 4740513a670SBarry Smith 4750513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 4760513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 477cddf8d76SBarry Smith color = DRAW_BLUE; 478416022c9SBarry Smith for ( i=0; i<m; i++ ) { 479cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 480416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 481cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 4823a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 483cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 484cddf8d76SBarry Smith #else 485cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 486cddf8d76SBarry Smith #endif 487*480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 488cddf8d76SBarry Smith } 489cddf8d76SBarry Smith } 490cddf8d76SBarry Smith color = DRAW_CYAN; 491cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 492cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 493cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 494cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 495cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 496*480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 497cddf8d76SBarry Smith } 498cddf8d76SBarry Smith } 499cddf8d76SBarry Smith color = DRAW_RED; 500cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 501cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 502cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 503cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5043a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 505cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 506cddf8d76SBarry Smith #else 507cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 508cddf8d76SBarry Smith #endif 509*480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 510416022c9SBarry Smith } 511416022c9SBarry Smith } 5120513a670SBarry Smith } else { 5130513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 5140513a670SBarry Smith /* first determine max of all nonzero values */ 5150513a670SBarry Smith int nz = a->nz,count; 5160513a670SBarry Smith Draw popup; 517*480ef9eaSBarry Smith double scale; 5180513a670SBarry Smith 5190513a670SBarry Smith for ( i=0; i<nz; i++ ) { 5200513a670SBarry Smith if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 5210513a670SBarry Smith } 522*480ef9eaSBarry Smith scale = (245.0 - DRAW_BASIC_COLORS)/maxv; 523522c5e43SBarry Smith ierr = DrawGetPopup(draw,&popup); CHKERRQ(ierr); 5240513a670SBarry Smith ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr); 5250513a670SBarry Smith count = 0; 5260513a670SBarry Smith for ( i=0; i<m; i++ ) { 5270513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 5280513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 5290513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 530*480ef9eaSBarry Smith color = DRAW_BASIC_COLORS + scale*PetscAbsScalar(a->a[count]); 531*480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 5320513a670SBarry Smith count++; 5330513a670SBarry Smith } 5340513a670SBarry Smith } 5350513a670SBarry Smith } 536*480ef9eaSBarry Smith PetscFunctionReturn(0); 537*480ef9eaSBarry Smith } 538cddf8d76SBarry Smith 539*480ef9eaSBarry Smith #undef __FUNC__ 540*480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw" 541*480ef9eaSBarry Smith int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 542*480ef9eaSBarry Smith { 543*480ef9eaSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 544*480ef9eaSBarry Smith int ierr; 545*480ef9eaSBarry Smith Draw draw; 546*480ef9eaSBarry Smith double xr,yr,xl,yl,h,w; 547*480ef9eaSBarry Smith 548*480ef9eaSBarry Smith PetscTruth isnull; 549*480ef9eaSBarry Smith 550*480ef9eaSBarry Smith PetscFunctionBegin; 551*480ef9eaSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 552*480ef9eaSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); 553*480ef9eaSBarry Smith if (isnull) PetscFunctionReturn(0); 554*480ef9eaSBarry Smith 555*480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 556*480ef9eaSBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 557*480ef9eaSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 558cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 559*480ef9eaSBarry Smith ierr = DrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A); CHKERRQ(ierr); 560*480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5613a40ed3dSBarry Smith PetscFunctionReturn(0); 562416022c9SBarry Smith } 563416022c9SBarry Smith 5645615d1e5SSatish Balay #undef __FUNC__ 565d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ" 566e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer) 567416022c9SBarry Smith { 568416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 569bcd2baecSBarry Smith ViewerType vtype; 570bcd2baecSBarry Smith int ierr; 571416022c9SBarry Smith 5723a40ed3dSBarry Smith PetscFunctionBegin; 573bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 574bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 5753a40ed3dSBarry Smith ierr = ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); CHKERRQ(ierr); 5765cd90555SBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 5773a40ed3dSBarry Smith ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr); 5785cd90555SBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 5793a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr); 5805cd90555SBarry Smith } else if (vtype == DRAW_VIEWER) { 5813a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr); 5825cd90555SBarry Smith } else { 5835cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 58417ab2063SBarry Smith } 5853a40ed3dSBarry Smith PetscFunctionReturn(0); 58617ab2063SBarry Smith } 58719bcc07fSBarry Smith 588c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 5895615d1e5SSatish Balay #undef __FUNC__ 5905615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ" 5918f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 59217ab2063SBarry Smith { 593416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 59441c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 59543ee02c3SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0; 596416022c9SBarry Smith Scalar *aa = a->a, *ap; 59717ab2063SBarry Smith 5983a40ed3dSBarry Smith PetscFunctionBegin; 5993a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 60017ab2063SBarry Smith 60143ee02c3SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 60217ab2063SBarry Smith for ( i=1; i<m; i++ ) { 603416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 60417ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 60594a9d846SBarry Smith rmax = PetscMax(rmax,ailen[i]); 60617ab2063SBarry Smith if (fshift) { 607416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 60817ab2063SBarry Smith N = ailen[i]; 60917ab2063SBarry Smith for ( j=0; j<N; j++ ) { 61017ab2063SBarry Smith ip[j-fshift] = ip[j]; 61117ab2063SBarry Smith ap[j-fshift] = ap[j]; 61217ab2063SBarry Smith } 61317ab2063SBarry Smith } 61417ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 61517ab2063SBarry Smith } 61617ab2063SBarry Smith if (m) { 61717ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 61817ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 61917ab2063SBarry Smith } 62017ab2063SBarry Smith /* reset ilen and imax for each row */ 62117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 62217ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 62317ab2063SBarry Smith } 624416022c9SBarry Smith a->nz = ai[m] + shift; 62517ab2063SBarry Smith 62617ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 627416022c9SBarry Smith if (fshift && a->diag) { 6280452661fSBarry Smith PetscFree(a->diag); 629416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 630416022c9SBarry Smith a->diag = 0; 63117ab2063SBarry Smith } 6324e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 6334e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 6344e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 635b810aeb4SBarry Smith a->reallocs); 63694a9d846SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax); 637dd5f02e7SSatish Balay a->reallocs = 0; 6384e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 6394e220ebcSLois Curfman McInnes 64076dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 64141c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 6423a40ed3dSBarry Smith PetscFunctionReturn(0); 64317ab2063SBarry Smith } 64417ab2063SBarry Smith 6455615d1e5SSatish Balay #undef __FUNC__ 6465615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ" 6478f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A) 64817ab2063SBarry Smith { 649416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6503a40ed3dSBarry Smith 6513a40ed3dSBarry Smith PetscFunctionBegin; 652cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 6533a40ed3dSBarry Smith PetscFunctionReturn(0); 65417ab2063SBarry Smith } 655416022c9SBarry Smith 6565615d1e5SSatish Balay #undef __FUNC__ 6575615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ" 658e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A) 65917ab2063SBarry Smith { 660416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 66182bf6240SBarry Smith int ierr; 662d5d45c9bSBarry Smith 6633a40ed3dSBarry Smith PetscFunctionBegin; 6643a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 665e1311b90SBarry Smith PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 66617ab2063SBarry Smith #endif 6670452661fSBarry Smith PetscFree(a->a); 6680452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 6690452661fSBarry Smith if (a->diag) PetscFree(a->diag); 6700452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 6710452661fSBarry Smith if (a->imax) PetscFree(a->imax); 6720452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 67376dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 67482bf6240SBarry Smith if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 6750452661fSBarry Smith PetscFree(a); 676eed86810SBarry Smith 677f2655603SLois Curfman McInnes PLogObjectDestroy(A); 678f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 6793a40ed3dSBarry Smith PetscFunctionReturn(0); 68017ab2063SBarry Smith } 68117ab2063SBarry Smith 6825615d1e5SSatish Balay #undef __FUNC__ 6835615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ" 6848f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A) 68517ab2063SBarry Smith { 6863a40ed3dSBarry Smith PetscFunctionBegin; 6873a40ed3dSBarry Smith PetscFunctionReturn(0); 68817ab2063SBarry Smith } 68917ab2063SBarry Smith 6905615d1e5SSatish Balay #undef __FUNC__ 6915615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ" 6928f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op) 69317ab2063SBarry Smith { 694416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6953a40ed3dSBarry Smith 6963a40ed3dSBarry Smith PetscFunctionBegin; 6976d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 6986d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 6996d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 700219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 7016d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 702c2653b3dSLois Curfman McInnes else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 70396854ed6SLois Curfman McInnes else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 7046d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 7056d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 706219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 7076d4a8577SBarry Smith op == MAT_SYMMETRIC || 7086d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 70990f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 710b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES|| 711b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) 712981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n"); 7133a40ed3dSBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) { 7143a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 7153a40ed3dSBarry Smith } else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 7166d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 7176d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 7186d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 7196d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 7203a40ed3dSBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 7213a40ed3dSBarry Smith PetscFunctionReturn(0); 72217ab2063SBarry Smith } 72317ab2063SBarry Smith 7245615d1e5SSatish Balay #undef __FUNC__ 7255615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ" 7268f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 72717ab2063SBarry Smith { 728416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7293a40ed3dSBarry Smith int i,j, n,shift = a->indexshift,ierr; 73017ab2063SBarry Smith Scalar *x, zero = 0.0; 73117ab2063SBarry Smith 7323a40ed3dSBarry Smith PetscFunctionBegin; 7333a40ed3dSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 734e1311b90SBarry Smith ierr = VecGetArray(v,&x);;CHKERRQ(ierr); 735e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr); 736a8c6a408SBarry Smith if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector"); 737416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 738416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 739416022c9SBarry Smith if (a->j[j]+shift == i) { 740416022c9SBarry Smith x[i] = a->a[j]; 74117ab2063SBarry Smith break; 74217ab2063SBarry Smith } 74317ab2063SBarry Smith } 74417ab2063SBarry Smith } 745e1311b90SBarry Smith ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr); 7463a40ed3dSBarry Smith PetscFunctionReturn(0); 74717ab2063SBarry Smith } 74817ab2063SBarry Smith 74917ab2063SBarry Smith /* -------------------------------------------------------*/ 75017ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 75117ab2063SBarry Smith /* -------------------------------------------------------*/ 7525615d1e5SSatish Balay #undef __FUNC__ 7535615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ" 75444cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 75517ab2063SBarry Smith { 756416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 75717ab2063SBarry Smith Scalar *x, *y, *v, alpha; 758e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx, shift = a->indexshift; 75917ab2063SBarry Smith 7603a40ed3dSBarry Smith PetscFunctionBegin; 761e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 762e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 763cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 76417ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 76517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 766416022c9SBarry Smith idx = a->j + a->i[i] + shift; 767416022c9SBarry Smith v = a->a + a->i[i] + shift; 768416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 76917ab2063SBarry Smith alpha = x[i]; 77017ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 77117ab2063SBarry Smith } 772416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 773e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 774e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 7753a40ed3dSBarry Smith PetscFunctionReturn(0); 77617ab2063SBarry Smith } 777d5d45c9bSBarry Smith 7785615d1e5SSatish Balay #undef __FUNC__ 7795615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ" 78044cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 78117ab2063SBarry Smith { 782416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 78317ab2063SBarry Smith Scalar *x, *y, *v, alpha; 784e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx,shift = a->indexshift; 78517ab2063SBarry Smith 7863a40ed3dSBarry Smith PetscFunctionBegin; 787e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 788e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 78917ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 79017ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 79117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 792416022c9SBarry Smith idx = a->j + a->i[i] + shift; 793416022c9SBarry Smith v = a->a + a->i[i] + shift; 794416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 79517ab2063SBarry Smith alpha = x[i]; 79617ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 79717ab2063SBarry Smith } 79890f02eecSBarry Smith PLogFlops(2*a->nz); 799e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 800e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8013a40ed3dSBarry Smith PetscFunctionReturn(0); 80217ab2063SBarry Smith } 80317ab2063SBarry Smith 8045615d1e5SSatish Balay #undef __FUNC__ 8055615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ" 80644cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 80717ab2063SBarry Smith { 808416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 80917ab2063SBarry Smith Scalar *x, *y, *v, sum; 810e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 811e36a17ebSSatish Balay #if !defined(USE_FORTRAN_KERNELS) 812e36a17ebSSatish Balay int n, i, jrow,j; 813e36a17ebSSatish Balay #endif 81417ab2063SBarry Smith 815fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT) 816fee21e36SBarry Smith #pragma disjoint(*x,*y,*v) 817fee21e36SBarry Smith #endif 818fee21e36SBarry Smith 8193a40ed3dSBarry Smith PetscFunctionBegin; 820e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 821e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 82217ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 823416022c9SBarry Smith idx = a->j; 824d64ed03dSBarry Smith v = a->a; 825416022c9SBarry Smith ii = a->i; 8268d195f9aSBarry Smith #if defined(USE_FORTRAN_KERNELS) 82729b5ca8aSSatish Balay fortranmultaij_(&m,x,ii,idx+shift,v+shift,y); 8288d195f9aSBarry Smith #else 829d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 830d64ed03dSBarry Smith idx += shift; 83117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8329ea0dfa2SSatish Balay jrow = ii[i]; 8339ea0dfa2SSatish Balay n = ii[i+1] - jrow; 83417ab2063SBarry Smith sum = 0.0; 8359ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8369ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8379ea0dfa2SSatish Balay } 83817ab2063SBarry Smith y[i] = sum; 83917ab2063SBarry Smith } 8408d195f9aSBarry Smith #endif 841416022c9SBarry Smith PLogFlops(2*a->nz - m); 842e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 843e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 8443a40ed3dSBarry Smith PetscFunctionReturn(0); 84517ab2063SBarry Smith } 84617ab2063SBarry Smith 8475615d1e5SSatish Balay #undef __FUNC__ 8485615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ" 84944cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 85017ab2063SBarry Smith { 851416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 85217ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 853e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 854e36a17ebSSatish Balay #if !defined(USE_FORTRAN_KERNELS) 855e36a17ebSSatish Balay int n,i,jrow,j; 856e36a17ebSSatish Balay #endif 8579ea0dfa2SSatish Balay 8583a40ed3dSBarry Smith PetscFunctionBegin; 859e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 860e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 861e1311b90SBarry Smith ierr = VecGetArray(zz,&z); CHKERRQ(ierr); 86217ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 863cddf8d76SBarry Smith idx = a->j; 864d64ed03dSBarry Smith v = a->a; 865cddf8d76SBarry Smith ii = a->i; 86602ab625aSSatish Balay #if defined(USE_FORTRAN_KERNELS) 86729b5ca8aSSatish Balay fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z); 86802ab625aSSatish Balay #else 869d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 870d64ed03dSBarry Smith idx += shift; 87117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8729ea0dfa2SSatish Balay jrow = ii[i]; 8739ea0dfa2SSatish Balay n = ii[i+1] - jrow; 87417ab2063SBarry Smith sum = y[i]; 8759ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8769ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8779ea0dfa2SSatish Balay } 87817ab2063SBarry Smith z[i] = sum; 87917ab2063SBarry Smith } 88002ab625aSSatish Balay #endif 881416022c9SBarry Smith PLogFlops(2*a->nz); 882e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 883e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 884e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr); 8853a40ed3dSBarry Smith PetscFunctionReturn(0); 88617ab2063SBarry Smith } 88717ab2063SBarry Smith 88817ab2063SBarry Smith /* 88917ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 89017ab2063SBarry Smith */ 89117ab2063SBarry Smith 8925615d1e5SSatish Balay #undef __FUNC__ 8935615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ" 894416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 89517ab2063SBarry Smith { 896416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 897416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 89817ab2063SBarry Smith 8993a40ed3dSBarry Smith PetscFunctionBegin; 9000452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 901416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 902416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 903416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 904416022c9SBarry Smith if (a->j[j]+shift == i) { 90517ab2063SBarry Smith diag[i] = j - shift; 90617ab2063SBarry Smith break; 90717ab2063SBarry Smith } 90817ab2063SBarry Smith } 90917ab2063SBarry Smith } 910416022c9SBarry Smith a->diag = diag; 9113a40ed3dSBarry Smith PetscFunctionReturn(0); 91217ab2063SBarry Smith } 91317ab2063SBarry Smith 9145615d1e5SSatish Balay #undef __FUNC__ 9155615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ" 91644cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 91717ab2063SBarry Smith double fshift,int its,Vec xx) 91817ab2063SBarry Smith { 919416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 920416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 921d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 92217ab2063SBarry Smith 9233a40ed3dSBarry Smith PetscFunctionBegin; 924e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 925e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 926d00d2cf4SBarry Smith if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);} 927416022c9SBarry Smith diag = a->diag; 928416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 92917ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 93017ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 93117ab2063SBarry Smith bs = b + shift; 93217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 933416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 934416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 935416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 936416022c9SBarry Smith v = a->a + diag[i] + (!shift); 93717ab2063SBarry Smith sum = b[i]*d/omega; 93817ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 93917ab2063SBarry Smith x[i] = sum; 94017ab2063SBarry Smith } 9413a40ed3dSBarry Smith PetscFunctionReturn(0); 94217ab2063SBarry Smith } 94317ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 944a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done"); 9453a40ed3dSBarry Smith } else if (flag & SOR_EISENSTAT) { 94617ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 94717ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 94817ab2063SBarry Smith 94917ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 95017ab2063SBarry Smith 95117ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 95217ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 95317ab2063SBarry Smith is the relaxation factor. 95417ab2063SBarry Smith */ 9550452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 95617ab2063SBarry Smith scale = (2.0/omega) - 1.0; 95717ab2063SBarry Smith 95817ab2063SBarry Smith /* x = (E + U)^{-1} b */ 95917ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 960416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 961416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 962416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 963416022c9SBarry Smith v = a->a + diag[i] + (!shift); 96417ab2063SBarry Smith sum = b[i]; 96517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 96617ab2063SBarry Smith x[i] = omega*(sum/d); 96717ab2063SBarry Smith } 96817ab2063SBarry Smith 96917ab2063SBarry Smith /* t = b - (2*E - D)x */ 970416022c9SBarry Smith v = a->a; 97117ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 97217ab2063SBarry Smith 97317ab2063SBarry Smith /* t = (E + L)^{-1}t */ 974416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 975416022c9SBarry Smith diag = a->diag; 97617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 977416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 978416022c9SBarry Smith n = diag[i] - a->i[i]; 979416022c9SBarry Smith idx = a->j + a->i[i] + shift; 980416022c9SBarry Smith v = a->a + a->i[i] + shift; 98117ab2063SBarry Smith sum = t[i]; 98217ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 98317ab2063SBarry Smith t[i] = omega*(sum/d); 98417ab2063SBarry Smith } 98517ab2063SBarry Smith 98617ab2063SBarry Smith /* x = x + t */ 98717ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 9880452661fSBarry Smith PetscFree(t); 9893a40ed3dSBarry Smith PetscFunctionReturn(0); 99017ab2063SBarry Smith } 99117ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 99217ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 99317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 994416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 995416022c9SBarry Smith n = diag[i] - a->i[i]; 996416022c9SBarry Smith idx = a->j + a->i[i] + shift; 997416022c9SBarry Smith v = a->a + a->i[i] + shift; 99817ab2063SBarry Smith sum = b[i]; 99917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 100017ab2063SBarry Smith x[i] = omega*(sum/d); 100117ab2063SBarry Smith } 100217ab2063SBarry Smith xb = x; 10033a40ed3dSBarry Smith } else xb = b; 100417ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 100517ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 100617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1007416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 100817ab2063SBarry Smith } 100917ab2063SBarry Smith } 101017ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 101117ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1012416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1013416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1014416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1015416022c9SBarry Smith v = a->a + diag[i] + (!shift); 101617ab2063SBarry Smith sum = xb[i]; 101717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 101817ab2063SBarry Smith x[i] = omega*(sum/d); 101917ab2063SBarry Smith } 102017ab2063SBarry Smith } 102117ab2063SBarry Smith its--; 102217ab2063SBarry Smith } 102317ab2063SBarry Smith while (its--) { 102417ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 102517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1026416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1027416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1028416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1029416022c9SBarry Smith v = a->a + a->i[i] + shift; 103017ab2063SBarry Smith sum = b[i]; 103117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 10327e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 103317ab2063SBarry Smith } 103417ab2063SBarry Smith } 103517ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 103617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1037416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1038416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1039416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1040416022c9SBarry Smith v = a->a + a->i[i] + shift; 104117ab2063SBarry Smith sum = b[i]; 104217ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 10437e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 104417ab2063SBarry Smith } 104517ab2063SBarry Smith } 104617ab2063SBarry Smith } 10473a40ed3dSBarry Smith PetscFunctionReturn(0); 104817ab2063SBarry Smith } 104917ab2063SBarry Smith 10505615d1e5SSatish Balay #undef __FUNC__ 10515615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ" 10528f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 105317ab2063SBarry Smith { 1054416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 10554e220ebcSLois Curfman McInnes 10563a40ed3dSBarry Smith PetscFunctionBegin; 10574e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 10584e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 10594e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 10604e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 10614e220ebcSLois Curfman McInnes info->block_size = 1.0; 10624e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 10634e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 10644e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 10654e220ebcSLois Curfman McInnes /* if (info->nz_unneeded != A->info.nz_unneeded) 10664e220ebcSLois Curfman McInnes printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 10674e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 10684e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 10694e220ebcSLois Curfman McInnes info->memory = A->mem; 10704e220ebcSLois Curfman McInnes if (A->factor) { 10714e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 10724e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 10734e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 10744e220ebcSLois Curfman McInnes } else { 10754e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 10764e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10774e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10784e220ebcSLois Curfman McInnes } 10793a40ed3dSBarry Smith PetscFunctionReturn(0); 108017ab2063SBarry Smith } 108117ab2063SBarry Smith 108217ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 108317ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 108417ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 108517ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 108617ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 108717ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 108817ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 108917ab2063SBarry Smith 10905615d1e5SSatish Balay #undef __FUNC__ 10915615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ" 10928f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 109317ab2063SBarry Smith { 1094416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1095416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 109617ab2063SBarry Smith 10973a40ed3dSBarry Smith PetscFunctionBegin; 109877c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 109917ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 110017ab2063SBarry Smith if (diag) { 110117ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1102a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1103416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 1104416022c9SBarry Smith a->ilen[rows[i]] = 1; 1105416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 1106416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 11073a40ed3dSBarry Smith } else { 1108d64ed03dSBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr); 110917ab2063SBarry Smith } 111017ab2063SBarry Smith } 11113a40ed3dSBarry Smith } else { 111217ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1113a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1114416022c9SBarry Smith a->ilen[rows[i]] = 0; 111517ab2063SBarry Smith } 111617ab2063SBarry Smith } 111717ab2063SBarry Smith ISRestoreIndices(is,&rows); 111843a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11193a40ed3dSBarry Smith PetscFunctionReturn(0); 112017ab2063SBarry Smith } 112117ab2063SBarry Smith 11225615d1e5SSatish Balay #undef __FUNC__ 11235615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ" 11248f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 112517ab2063SBarry Smith { 1126416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11273a40ed3dSBarry Smith 11283a40ed3dSBarry Smith PetscFunctionBegin; 11290752156aSBarry Smith if (m) *m = a->m; 11300752156aSBarry Smith if (n) *n = a->n; 11313a40ed3dSBarry Smith PetscFunctionReturn(0); 113217ab2063SBarry Smith } 113317ab2063SBarry Smith 11345615d1e5SSatish Balay #undef __FUNC__ 11355615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ" 11368f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 113717ab2063SBarry Smith { 1138416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11393a40ed3dSBarry Smith 11403a40ed3dSBarry Smith PetscFunctionBegin; 1141416022c9SBarry Smith *m = 0; *n = a->m; 11423a40ed3dSBarry Smith PetscFunctionReturn(0); 114317ab2063SBarry Smith } 1144026e39d0SSatish Balay 11455615d1e5SSatish Balay #undef __FUNC__ 11465615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ" 11474e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 114817ab2063SBarry Smith { 1149416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1150c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 115117ab2063SBarry Smith 11523a40ed3dSBarry Smith PetscFunctionBegin; 1153a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 115417ab2063SBarry Smith 1155416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 1156416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 115717ab2063SBarry Smith if (idx) { 1158416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 11594e093b46SBarry Smith if (*nz && shift) { 11600452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 116117ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 11624e093b46SBarry Smith } else if (*nz) { 11634e093b46SBarry Smith *idx = itmp; 116417ab2063SBarry Smith } 116517ab2063SBarry Smith else *idx = 0; 116617ab2063SBarry Smith } 11673a40ed3dSBarry Smith PetscFunctionReturn(0); 116817ab2063SBarry Smith } 116917ab2063SBarry Smith 11705615d1e5SSatish Balay #undef __FUNC__ 11715615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ" 11724e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 117317ab2063SBarry Smith { 11744e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11753a40ed3dSBarry Smith 11763a40ed3dSBarry Smith PetscFunctionBegin; 11774e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 11783a40ed3dSBarry Smith PetscFunctionReturn(0); 117917ab2063SBarry Smith } 118017ab2063SBarry Smith 11815615d1e5SSatish Balay #undef __FUNC__ 11825615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ" 11838f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 118417ab2063SBarry Smith { 1185416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1186416022c9SBarry Smith Scalar *v = a->a; 118717ab2063SBarry Smith double sum = 0.0; 1188416022c9SBarry Smith int i, j,shift = a->indexshift; 118917ab2063SBarry Smith 11903a40ed3dSBarry Smith PetscFunctionBegin; 119117ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1192416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 11933a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 119417ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 119517ab2063SBarry Smith #else 119617ab2063SBarry Smith sum += (*v)*(*v); v++; 119717ab2063SBarry Smith #endif 119817ab2063SBarry Smith } 119917ab2063SBarry Smith *norm = sqrt(sum); 12003a40ed3dSBarry Smith } else if (type == NORM_1) { 120117ab2063SBarry Smith double *tmp; 1202416022c9SBarry Smith int *jj = a->j; 120366963ce1SSatish Balay tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp); 1204cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 120517ab2063SBarry Smith *norm = 0.0; 1206416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1207a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 120817ab2063SBarry Smith } 1209416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 121017ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 121117ab2063SBarry Smith } 12120452661fSBarry Smith PetscFree(tmp); 12133a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 121417ab2063SBarry Smith *norm = 0.0; 1215416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1216416022c9SBarry Smith v = a->a + a->i[j] + shift; 121717ab2063SBarry Smith sum = 0.0; 1218416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1219cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 122017ab2063SBarry Smith } 122117ab2063SBarry Smith if (sum > *norm) *norm = sum; 122217ab2063SBarry Smith } 12233a40ed3dSBarry Smith } else { 1224a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 122517ab2063SBarry Smith } 12263a40ed3dSBarry Smith PetscFunctionReturn(0); 122717ab2063SBarry Smith } 122817ab2063SBarry Smith 12295615d1e5SSatish Balay #undef __FUNC__ 12305615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ" 12318f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B) 123217ab2063SBarry Smith { 1233416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1234416022c9SBarry Smith Mat C; 1235416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1236416022c9SBarry Smith int shift = a->indexshift; 1237d5d45c9bSBarry Smith Scalar *array = a->a; 123817ab2063SBarry Smith 12393a40ed3dSBarry Smith PetscFunctionBegin; 1240a8c6a408SBarry Smith if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place"); 12410452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1242cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 124317ab2063SBarry Smith if (shift) { 124417ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 124517ab2063SBarry Smith } 124617ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1247416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 12480452661fSBarry Smith PetscFree(col); 124917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 125017ab2063SBarry Smith len = ai[i+1]-ai[i]; 1251416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 125217ab2063SBarry Smith array += len; aj += len; 125317ab2063SBarry Smith } 125417ab2063SBarry Smith if (shift) { 125517ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 125617ab2063SBarry Smith } 125717ab2063SBarry Smith 12586d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 12596d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 126017ab2063SBarry Smith 12613638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1262416022c9SBarry Smith *B = C; 126317ab2063SBarry Smith } else { 1264f830108cSBarry Smith PetscOps *Abops; 1265f830108cSBarry Smith struct _MatOps *Aops; 1266f830108cSBarry Smith 1267416022c9SBarry Smith /* This isn't really an in-place transpose */ 12680452661fSBarry Smith PetscFree(a->a); 12690452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 12700452661fSBarry Smith if (a->diag) PetscFree(a->diag); 12710452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 12720452661fSBarry Smith if (a->imax) PetscFree(a->imax); 12730452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 12741073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 12750452661fSBarry Smith PetscFree(a); 1276f830108cSBarry Smith 1277f830108cSBarry Smith /* 1278f830108cSBarry Smith This is horrible, horrible code. We need to keep the 12798f0f457cSSatish Balay the bops and ops Structures, copy everything from C 12808f0f457cSSatish Balay including the function pointers.. 1281f830108cSBarry Smith */ 12828f0f457cSSatish Balay PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps)); 12838f0f457cSSatish Balay PetscMemcpy(A->bops,C->bops,sizeof(PetscOps)); 1284f830108cSBarry Smith Abops = A->bops; 1285f830108cSBarry Smith Aops = A->ops; 1286f09e8eb9SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 1287f830108cSBarry Smith A->bops = Abops; 1288f830108cSBarry Smith A->ops = Aops; 128927a8da17SBarry Smith A->qlist = 0; 1290f830108cSBarry Smith 12910452661fSBarry Smith PetscHeaderDestroy(C); 129217ab2063SBarry Smith } 12933a40ed3dSBarry Smith PetscFunctionReturn(0); 129417ab2063SBarry Smith } 129517ab2063SBarry Smith 12965615d1e5SSatish Balay #undef __FUNC__ 12975615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ" 12988f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 129917ab2063SBarry Smith { 1300416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 130117ab2063SBarry Smith Scalar *l,*r,x,*v; 1302e1311b90SBarry Smith int ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 130317ab2063SBarry Smith 13043a40ed3dSBarry Smith PetscFunctionBegin; 130517ab2063SBarry Smith if (ll) { 13063ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 13073ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 1308e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr); 1309a8c6a408SBarry Smith if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length"); 1310e1311b90SBarry Smith ierr = VecGetArray(ll,&l); CHKERRQ(ierr); 1311416022c9SBarry Smith v = a->a; 131217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 131317ab2063SBarry Smith x = l[i]; 1314416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 131517ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 131617ab2063SBarry Smith } 1317e1311b90SBarry Smith ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr); 131844cd7ae7SLois Curfman McInnes PLogFlops(nz); 131917ab2063SBarry Smith } 132017ab2063SBarry Smith if (rr) { 1321e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr); 1322a8c6a408SBarry Smith if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length"); 1323e1311b90SBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1324416022c9SBarry Smith v = a->a; jj = a->j; 132517ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 132617ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 132717ab2063SBarry Smith } 1328e1311b90SBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 132944cd7ae7SLois Curfman McInnes PLogFlops(nz); 133017ab2063SBarry Smith } 13313a40ed3dSBarry Smith PetscFunctionReturn(0); 133217ab2063SBarry Smith } 133317ab2063SBarry Smith 13345615d1e5SSatish Balay #undef __FUNC__ 13355615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ" 13366a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B) 133717ab2063SBarry Smith { 1338db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 1339d5db1b72SBarry Smith int *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 134099141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1341a2744918SBarry Smith register int sum,lensi; 134299141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 134399141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 134499141d43SSatish Balay Scalar *a_new,*mat_a; 1345416022c9SBarry Smith Mat C; 1346fee21e36SBarry Smith PetscTruth stride; 134717ab2063SBarry Smith 13483a40ed3dSBarry Smith PetscFunctionBegin; 1349d64ed03dSBarry Smith ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr); 1350a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted"); 1351d64ed03dSBarry Smith ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 1352a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted"); 135399141d43SSatish Balay 135417ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 135517ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 135617ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 135717ab2063SBarry Smith 1358fee21e36SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr); 1359fee21e36SBarry Smith ierr = ISStride(iscol,&stride); CHKERRQ(ierr); 1360fee21e36SBarry Smith if (stride && step == 1) { 136102834360SBarry Smith /* special case of contiguous rows */ 136257faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 136302834360SBarry Smith starts = lens + ncols; 136402834360SBarry Smith /* loop over new rows determining lens and starting points */ 136502834360SBarry Smith for (i=0; i<nrows; i++) { 1366a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1367a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 136802834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1369d8ced48eSBarry Smith if (aj[k]+shift >= first) { 137002834360SBarry Smith starts[i] = k; 137102834360SBarry Smith break; 137202834360SBarry Smith } 137302834360SBarry Smith } 1374a2744918SBarry Smith sum = 0; 137502834360SBarry Smith while (k < kend) { 1376d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1377a2744918SBarry Smith sum++; 137802834360SBarry Smith } 1379a2744918SBarry Smith lens[i] = sum; 138002834360SBarry Smith } 138102834360SBarry Smith /* create submatrix */ 1382cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 138308480c60SBarry Smith int n_cols,n_rows; 138408480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1385a8c6a408SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); 1386d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 138708480c60SBarry Smith C = *B; 13883a40ed3dSBarry Smith } else { 138902834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 139008480c60SBarry Smith } 1391db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1392db02288aSLois Curfman McInnes 139302834360SBarry Smith /* loop over rows inserting into submatrix */ 1394db02288aSLois Curfman McInnes a_new = c->a; 1395db02288aSLois Curfman McInnes j_new = c->j; 1396db02288aSLois Curfman McInnes i_new = c->i; 1397db02288aSLois Curfman McInnes i_new[0] = -shift; 139802834360SBarry Smith for (i=0; i<nrows; i++) { 1399a2744918SBarry Smith ii = starts[i]; 1400a2744918SBarry Smith lensi = lens[i]; 1401a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1402a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 140302834360SBarry Smith } 1404a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1405a2744918SBarry Smith a_new += lensi; 1406a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1407a2744918SBarry Smith c->ilen[i] = lensi; 140802834360SBarry Smith } 14090452661fSBarry Smith PetscFree(lens); 14103a40ed3dSBarry Smith } else { 141102834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 14120452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 141302834360SBarry Smith ssmap = smap + shift; 141499141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1415cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 141617ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 141702834360SBarry Smith /* determine lens of each row */ 141802834360SBarry Smith for (i=0; i<nrows; i++) { 1419d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 142002834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 142102834360SBarry Smith lens[i] = 0; 142202834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1423d8ced48eSBarry Smith if (ssmap[aj[k]]) { 142402834360SBarry Smith lens[i]++; 142502834360SBarry Smith } 142602834360SBarry Smith } 142702834360SBarry Smith } 142817ab2063SBarry Smith /* Create and fill new matrix */ 1429a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 143099141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 1431a8c6a408SBarry Smith if (c->m != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size"); 143299141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1433a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros"); 143499141d43SSatish Balay } 143599141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 143608480c60SBarry Smith C = *B; 14373a40ed3dSBarry Smith } else { 143802834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 143908480c60SBarry Smith } 144099141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 144117ab2063SBarry Smith for (i=0; i<nrows; i++) { 144299141d43SSatish Balay row = irow[i]; 144399141d43SSatish Balay kstart = ai[row]+shift; 144499141d43SSatish Balay kend = kstart + a->ilen[row]; 144599141d43SSatish Balay mat_i = c->i[i]+shift; 144699141d43SSatish Balay mat_j = c->j + mat_i; 144799141d43SSatish Balay mat_a = c->a + mat_i; 144899141d43SSatish Balay mat_ilen = c->ilen + i; 144917ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 145099141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 145199141d43SSatish Balay *mat_j++ = tcol - (!shift); 145299141d43SSatish Balay *mat_a++ = a->a[k]; 145399141d43SSatish Balay (*mat_ilen)++; 145499141d43SSatish Balay 145517ab2063SBarry Smith } 145617ab2063SBarry Smith } 145717ab2063SBarry Smith } 145802834360SBarry Smith /* Free work space */ 145902834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 146099141d43SSatish Balay PetscFree(smap); PetscFree(lens); 146102834360SBarry Smith } 14626d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 14636d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 146417ab2063SBarry Smith 146517ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1466416022c9SBarry Smith *B = C; 14673a40ed3dSBarry Smith PetscFunctionReturn(0); 146817ab2063SBarry Smith } 146917ab2063SBarry Smith 1470a871dcd8SBarry Smith /* 147163b91edcSBarry Smith note: This can only work for identity for row and col. It would 147263b91edcSBarry Smith be good to check this and otherwise generate an error. 1473a871dcd8SBarry Smith */ 14745615d1e5SSatish Balay #undef __FUNC__ 14755615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ" 14768f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1477a871dcd8SBarry Smith { 147863b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 147908480c60SBarry Smith int ierr; 148063b91edcSBarry Smith Mat outA; 148163b91edcSBarry Smith 14823a40ed3dSBarry Smith PetscFunctionBegin; 1483a8c6a408SBarry Smith if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 1484a871dcd8SBarry Smith 148563b91edcSBarry Smith outA = inA; 148663b91edcSBarry Smith inA->factor = FACTOR_LU; 148763b91edcSBarry Smith a->row = row; 148863b91edcSBarry Smith a->col = col; 148963b91edcSBarry Smith 1490f0ec6fceSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1491f0ec6fceSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 14921e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1493f0ec6fceSSatish Balay 149494a9d846SBarry Smith if (!a->solve_work) { /* this matrix may have been factored before */ 14950452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work); 149694a9d846SBarry Smith } 149763b91edcSBarry Smith 149808480c60SBarry Smith if (!a->diag) { 149908480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 150063b91edcSBarry Smith } 150163b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 15023a40ed3dSBarry Smith PetscFunctionReturn(0); 1503a871dcd8SBarry Smith } 1504a871dcd8SBarry Smith 150574cdf0dfSBarry Smith #include "pinclude/blaslapack.h" 15065615d1e5SSatish Balay #undef __FUNC__ 15075615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ" 15088f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1509f0b747eeSBarry Smith { 1510f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1511f0b747eeSBarry Smith int one = 1; 15123a40ed3dSBarry Smith 15133a40ed3dSBarry Smith PetscFunctionBegin; 1514f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1515f0b747eeSBarry Smith PLogFlops(a->nz); 15163a40ed3dSBarry Smith PetscFunctionReturn(0); 1517f0b747eeSBarry Smith } 1518f0b747eeSBarry Smith 15195615d1e5SSatish Balay #undef __FUNC__ 15205615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ" 15216a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B) 1522cddf8d76SBarry Smith { 1523cddf8d76SBarry Smith int ierr,i; 1524cddf8d76SBarry Smith 15253a40ed3dSBarry Smith PetscFunctionBegin; 1526cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 15270452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1528cddf8d76SBarry Smith } 1529cddf8d76SBarry Smith 1530cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 15316a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1532cddf8d76SBarry Smith } 15333a40ed3dSBarry Smith PetscFunctionReturn(0); 1534cddf8d76SBarry Smith } 1535cddf8d76SBarry Smith 15365615d1e5SSatish Balay #undef __FUNC__ 15375615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ" 15388f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 15395a838052SSatish Balay { 1540f830108cSBarry Smith PetscFunctionBegin; 15415a838052SSatish Balay *bs = 1; 15423a40ed3dSBarry Smith PetscFunctionReturn(0); 15435a838052SSatish Balay } 15445a838052SSatish Balay 15455615d1e5SSatish Balay #undef __FUNC__ 15465615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ" 15478f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 15484dcbc457SBarry Smith { 1549e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 155006763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 15518a047759SSatish Balay int start, end, *ai, *aj; 1552bbd702dbSSatish Balay BT table; 1553bbd702dbSSatish Balay 15543a40ed3dSBarry Smith PetscFunctionBegin; 15558a047759SSatish Balay shift = a->indexshift; 1556e4d965acSSatish Balay m = a->m; 1557e4d965acSSatish Balay ai = a->i; 15588a047759SSatish Balay aj = a->j+shift; 15598a047759SSatish Balay 1560a8c6a408SBarry Smith if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used"); 156106763907SSatish Balay 156206763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 1563bbd702dbSSatish Balay ierr = BTCreate(m,table); CHKERRQ(ierr); 156406763907SSatish Balay 1565e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1566b97fc60eSLois Curfman McInnes /* Initialize the two local arrays */ 1567e4d965acSSatish Balay isz = 0; 1568bbd702dbSSatish Balay BTMemzero(m,table); 1569e4d965acSSatish Balay 1570e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 15714dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 157277c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1573e4d965acSSatish Balay 1574dd097bc3SLois Curfman McInnes /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 1575e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 1576bbd702dbSSatish Balay if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];} 15774dcbc457SBarry Smith } 157806763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 157906763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1580e4d965acSSatish Balay 158104a348a9SBarry Smith k = 0; 158204a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap */ 158304a348a9SBarry Smith n = isz; 158406763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1585e4d965acSSatish Balay row = nidx[k]; 1586e4d965acSSatish Balay start = ai[row]; 1587e4d965acSSatish Balay end = ai[row+1]; 158804a348a9SBarry Smith for ( l = start; l<end ; l++){ 15898a047759SSatish Balay val = aj[l] + shift; 15902a8f2162SSatish Balay if (!BTLookupSet(table,val)) {nidx[isz++] = val;} 1591e4d965acSSatish Balay } 1592e4d965acSSatish Balay } 1593e4d965acSSatish Balay } 1594029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1595e4d965acSSatish Balay } 1596bbd702dbSSatish Balay BTDestroy(table); 159706763907SSatish Balay PetscFree(nidx); 15983a40ed3dSBarry Smith PetscFunctionReturn(0); 15994dcbc457SBarry Smith } 160017ab2063SBarry Smith 16010513a670SBarry Smith /* -------------------------------------------------------------- */ 16020513a670SBarry Smith #undef __FUNC__ 16030513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ" 16040513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 16050513a670SBarry Smith { 16060513a670SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 16070513a670SBarry Smith Scalar *vwork; 16080513a670SBarry Smith int i, ierr, nz, m = a->m, n = a->n, *cwork; 16090513a670SBarry Smith int *row,*col,*cnew,j,*lens; 161056cd22aeSBarry Smith IS icolp,irowp; 16110513a670SBarry Smith 16123a40ed3dSBarry Smith PetscFunctionBegin; 161356cd22aeSBarry Smith ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr); 161456cd22aeSBarry Smith ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr); 161556cd22aeSBarry Smith ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr); 161656cd22aeSBarry Smith ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr); 16170513a670SBarry Smith 16180513a670SBarry Smith /* determine lengths of permuted rows */ 16190513a670SBarry Smith lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens); 16200513a670SBarry Smith for (i=0; i<m; i++ ) { 16210513a670SBarry Smith lens[row[i]] = a->i[i+1] - a->i[i]; 16220513a670SBarry Smith } 16230513a670SBarry Smith ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr); 16240513a670SBarry Smith PetscFree(lens); 16250513a670SBarry Smith 16260513a670SBarry Smith cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew); 16270513a670SBarry Smith for (i=0; i<m; i++) { 16280513a670SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 16290513a670SBarry Smith for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];} 16300513a670SBarry Smith ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr); 16310513a670SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 16320513a670SBarry Smith } 16330513a670SBarry Smith PetscFree(cnew); 16340513a670SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 16350513a670SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 163656cd22aeSBarry Smith ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr); 163756cd22aeSBarry Smith ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr); 163856cd22aeSBarry Smith ierr = ISDestroy(irowp); CHKERRQ(ierr); 163956cd22aeSBarry Smith ierr = ISDestroy(icolp); CHKERRQ(ierr); 16403a40ed3dSBarry Smith PetscFunctionReturn(0); 16410513a670SBarry Smith } 16420513a670SBarry Smith 16435615d1e5SSatish Balay #undef __FUNC__ 16445615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ" 1645682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1646682d7d0cSBarry Smith { 1647682d7d0cSBarry Smith static int called = 0; 1648682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1649682d7d0cSBarry Smith 16503a40ed3dSBarry Smith PetscFunctionBegin; 16513a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 165276be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 165376be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n"); 165476be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n"); 165576be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_no_inode: Do not use inodes\n"); 165676be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n"); 1657682d7d0cSBarry Smith #if defined(HAVE_ESSL) 165876be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.\n"); 1659682d7d0cSBarry Smith #endif 16603a40ed3dSBarry Smith PetscFunctionReturn(0); 1661682d7d0cSBarry Smith } 16628f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1663a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1664a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1665a93ec695SBarry Smith 1666682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 166717ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 166817ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1669416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1670416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 167117ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 167217ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 167317ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 167417ab2063SBarry Smith MatRelax_SeqAIJ, 167517ab2063SBarry Smith MatTranspose_SeqAIJ, 16767264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1677f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 167817ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 167917ab2063SBarry Smith MatCompress_SeqAIJ, 168017ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 168117ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 168217ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 168317ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 168494a9d846SBarry Smith 0,0, 16853d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1686cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 16877eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1688682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1689f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 16905a838052SSatish Balay MatScale_SeqAIJ,0,0, 16916945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 16926945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 16933b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 16943b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 16953b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1696a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1697a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 16980513a670SBarry Smith MatColoringPatch_SeqAIJ, 16990513a670SBarry Smith 0, 17000513a670SBarry Smith MatPermute_SeqAIJ}; 170117ab2063SBarry Smith 170217ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 170317ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 170417ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 170517ab2063SBarry Smith 17065615d1e5SSatish Balay #undef __FUNC__ 1707bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ" 1708bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices) 1709bef8e0ddSBarry Smith { 1710bef8e0ddSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1711bef8e0ddSBarry Smith int i,nz,n; 1712bef8e0ddSBarry Smith 1713bef8e0ddSBarry Smith PetscFunctionBegin; 1714bef8e0ddSBarry Smith if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing"); 1715bef8e0ddSBarry Smith 1716bef8e0ddSBarry Smith nz = aij->maxnz; 1717bef8e0ddSBarry Smith n = aij->n; 1718bef8e0ddSBarry Smith for (i=0; i<nz; i++) { 1719bef8e0ddSBarry Smith aij->j[i] = indices[i]; 1720bef8e0ddSBarry Smith } 1721bef8e0ddSBarry Smith aij->nz = nz; 1722bef8e0ddSBarry Smith for ( i=0; i<n; i++ ) { 1723bef8e0ddSBarry Smith aij->ilen[i] = aij->imax[i]; 1724bef8e0ddSBarry Smith } 1725bef8e0ddSBarry Smith 1726bef8e0ddSBarry Smith PetscFunctionReturn(0); 1727bef8e0ddSBarry Smith } 1728bef8e0ddSBarry Smith 1729bef8e0ddSBarry Smith #undef __FUNC__ 1730bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices" 1731bef8e0ddSBarry Smith /*@ 1732bef8e0ddSBarry Smith MatSeqAIJSetColumnIndices - Set the column indices for all the rows 1733bef8e0ddSBarry Smith in the matrix. 1734bef8e0ddSBarry Smith 1735bef8e0ddSBarry Smith Input Parameters: 1736bef8e0ddSBarry Smith + mat - the SeqAIJ matrix 1737bef8e0ddSBarry Smith - indices - the column indices 1738bef8e0ddSBarry Smith 1739bef8e0ddSBarry Smith Notes: 1740bef8e0ddSBarry Smith This can be called if you have precomputed the nonzero structure of the 1741bef8e0ddSBarry Smith matrix and want to provide it to the matrix object to improve the performance 1742bef8e0ddSBarry Smith of the MatSetValues() operation. 1743bef8e0ddSBarry Smith 1744bef8e0ddSBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 1745bef8e0ddSBarry Smith MatCreateSeqAIJ(). 1746bef8e0ddSBarry Smith 1747bef8e0ddSBarry Smith MUST be called before any calls to MatSetValues(); 1748bef8e0ddSBarry Smith 1749bef8e0ddSBarry Smith @*/ 1750bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices) 1751bef8e0ddSBarry Smith { 1752bef8e0ddSBarry Smith int ierr,(*f)(Mat,int *); 1753bef8e0ddSBarry Smith 1754bef8e0ddSBarry Smith PetscFunctionBegin; 1755bef8e0ddSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1756bef8e0ddSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f); 1757bef8e0ddSBarry Smith CHKERRQ(ierr); 1758bef8e0ddSBarry Smith if (f) { 1759bef8e0ddSBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 1760bef8e0ddSBarry Smith } else { 1761bef8e0ddSBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 1762bef8e0ddSBarry Smith } 1763bef8e0ddSBarry Smith PetscFunctionReturn(0); 1764bef8e0ddSBarry Smith } 1765bef8e0ddSBarry Smith 1766bef8e0ddSBarry Smith #undef __FUNC__ 17675615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ" 176817ab2063SBarry Smith /*@C 1769682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 17700d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 17716e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 17722bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 17732bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 177417ab2063SBarry Smith 1775db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1776db81eaa0SLois Curfman McInnes 177717ab2063SBarry Smith Input Parameters: 1778db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 177917ab2063SBarry Smith . m - number of rows 178017ab2063SBarry Smith . n - number of columns 178117ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 1782db81eaa0SLois Curfman McInnes - nzz - array containing the number of nonzeros in the various rows 17832bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 178417ab2063SBarry Smith 178517ab2063SBarry Smith Output Parameter: 1786416022c9SBarry Smith . A - the matrix 178717ab2063SBarry Smith 1788b259b22eSLois Curfman McInnes Notes: 178917ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 179017ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 17910002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 179244cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 179317ab2063SBarry Smith 179417ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1795a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 17963d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 17976da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 179817ab2063SBarry Smith 1799682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 18004fca80b9SLois Curfman McInnes improve numerical efficiency of matrix-vector products and solves. We 1801682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 18026c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 18036c7ebb05SLois Curfman McInnes 18046c7ebb05SLois Curfman McInnes Options Database Keys: 1805db81eaa0SLois Curfman McInnes + -mat_aij_no_inode - Do not use inodes 1806db81eaa0SLois Curfman McInnes . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5) 1807db81eaa0SLois Curfman McInnes - -mat_aij_oneindex - Internally use indexing starting at 1 1808db81eaa0SLois Curfman McInnes rather than 0. Note that when calling MatSetValues(), 1809db81eaa0SLois Curfman McInnes the user still MUST index entries starting at 0! 181017ab2063SBarry Smith 1811bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices() 181217ab2063SBarry Smith @*/ 1813416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 181417ab2063SBarry Smith { 1815416022c9SBarry Smith Mat B; 1816416022c9SBarry Smith Mat_SeqAIJ *b; 18176945ee14SBarry Smith int i, len, ierr, flg,size; 18186945ee14SBarry Smith 18193a40ed3dSBarry Smith PetscFunctionBegin; 18206945ee14SBarry Smith MPI_Comm_size(comm,&size); 1821a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1"); 1822d5d45c9bSBarry Smith 1823416022c9SBarry Smith *A = 0; 1824f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView); 1825416022c9SBarry Smith PLogObjectCreate(B); 18260452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 182744cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1828f830108cSBarry Smith PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps)); 1829e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqAIJ; 1830e1311b90SBarry Smith B->ops->view = MatView_SeqAIJ; 1831416022c9SBarry Smith B->factor = 0; 1832416022c9SBarry Smith B->lupivotthreshold = 1.0; 183390f02eecSBarry Smith B->mapping = 0; 1834e1311b90SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg); CHKERRQ(ierr); 18357a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 1836e1311b90SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*) &b->ilu_preserve_row_sums);CHKERRQ(ierr); 1837416022c9SBarry Smith b->row = 0; 1838416022c9SBarry Smith b->col = 0; 183982bf6240SBarry Smith b->icol = 0; 1840416022c9SBarry Smith b->indexshift = 0; 1841b810aeb4SBarry Smith b->reallocs = 0; 184269957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 184369957df2SSatish Balay if (flg) b->indexshift = -1; 184417ab2063SBarry Smith 184544cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 184644cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 18470452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1848b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 18497b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 18507b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1851416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 185217ab2063SBarry Smith nz = nz*m; 18533a40ed3dSBarry Smith } else { 185417ab2063SBarry Smith nz = 0; 1855416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 185617ab2063SBarry Smith } 185717ab2063SBarry Smith 185817ab2063SBarry Smith /* allocate the matrix space */ 1859416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 18600452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1861416022c9SBarry Smith b->j = (int *) (b->a + nz); 1862cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1863416022c9SBarry Smith b->i = b->j + nz; 1864416022c9SBarry Smith b->singlemalloc = 1; 186517ab2063SBarry Smith 1866416022c9SBarry Smith b->i[0] = -b->indexshift; 186717ab2063SBarry Smith for (i=1; i<m+1; i++) { 1868416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 186917ab2063SBarry Smith } 187017ab2063SBarry Smith 1871416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 18720452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1873f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 1874416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 187517ab2063SBarry Smith 1876416022c9SBarry Smith b->nz = 0; 1877416022c9SBarry Smith b->maxnz = nz; 1878416022c9SBarry Smith b->sorted = 0; 1879416022c9SBarry Smith b->roworiented = 1; 1880416022c9SBarry Smith b->nonew = 0; 1881416022c9SBarry Smith b->diag = 0; 1882416022c9SBarry Smith b->solve_work = 0; 188371bd300dSLois Curfman McInnes b->spptr = 0; 1884754ec7b1SSatish Balay b->inode.node_count = 0; 1885754ec7b1SSatish Balay b->inode.size = 0; 18866c7ebb05SLois Curfman McInnes b->inode.limit = 5; 18876c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 18884e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 188917ab2063SBarry Smith 1890416022c9SBarry Smith *A = B; 18914e220ebcSLois Curfman McInnes 18924b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 18934b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 189469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 189569957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 18964b14c69eSBarry Smith #endif 189769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 189869957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 189969957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 190069957df2SSatish Balay if (flg) { 1901a8c6a408SBarry Smith if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml"); 1902416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 190317ab2063SBarry Smith } 190469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 190569957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 1906bef8e0ddSBarry Smith 1907bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C", 1908bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 1909bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 19103a40ed3dSBarry Smith PetscFunctionReturn(0); 191117ab2063SBarry Smith } 191217ab2063SBarry Smith 19135615d1e5SSatish Balay #undef __FUNC__ 19145615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ" 19153d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 191617ab2063SBarry Smith { 1917416022c9SBarry Smith Mat C; 1918416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 1919bef8e0ddSBarry Smith int i,len, m = a->m,shift = a->indexshift,ierr; 192017ab2063SBarry Smith 19213a40ed3dSBarry Smith PetscFunctionBegin; 19224043dd9cSLois Curfman McInnes *B = 0; 1923f830108cSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView); 1924416022c9SBarry Smith PLogObjectCreate(C); 19250452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 1926f830108cSBarry Smith PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 1927e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqAIJ; 1928e1311b90SBarry Smith C->ops->view = MatView_SeqAIJ; 1929416022c9SBarry Smith C->factor = A->factor; 1930416022c9SBarry Smith c->row = 0; 1931416022c9SBarry Smith c->col = 0; 193282bf6240SBarry Smith c->icol = 0; 1933416022c9SBarry Smith c->indexshift = shift; 1934c456f294SBarry Smith C->assembled = PETSC_TRUE; 193517ab2063SBarry Smith 193644cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 193744cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 193844cd7ae7SLois Curfman McInnes C->M = a->m; 193944cd7ae7SLois Curfman McInnes C->N = a->n; 194017ab2063SBarry Smith 19410452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 19420452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 194317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1944416022c9SBarry Smith c->imax[i] = a->imax[i]; 1945416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 194617ab2063SBarry Smith } 194717ab2063SBarry Smith 194817ab2063SBarry Smith /* allocate the matrix space */ 1949416022c9SBarry Smith c->singlemalloc = 1; 1950416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 19510452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1952416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1953416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1954416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 195517ab2063SBarry Smith if (m > 0) { 1956416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 195708480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1958416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 195917ab2063SBarry Smith } 196008480c60SBarry Smith } 196117ab2063SBarry Smith 1962f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 1963416022c9SBarry Smith c->sorted = a->sorted; 1964416022c9SBarry Smith c->roworiented = a->roworiented; 1965416022c9SBarry Smith c->nonew = a->nonew; 19667a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 196717ab2063SBarry Smith 1968416022c9SBarry Smith if (a->diag) { 19690452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1970416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 197117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1972416022c9SBarry Smith c->diag[i] = a->diag[i]; 197317ab2063SBarry Smith } 19743a40ed3dSBarry Smith } else c->diag = 0; 19756c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 19766c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1977754ec7b1SSatish Balay if (a->inode.size){ 1978daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 1979754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1980daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 1981754ec7b1SSatish Balay } else { 1982754ec7b1SSatish Balay c->inode.size = 0; 1983754ec7b1SSatish Balay c->inode.node_count = 0; 1984754ec7b1SSatish Balay } 1985416022c9SBarry Smith c->nz = a->nz; 1986416022c9SBarry Smith c->maxnz = a->maxnz; 1987416022c9SBarry Smith c->solve_work = 0; 198876dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1989754ec7b1SSatish Balay 1990416022c9SBarry Smith *B = C; 1991bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C", 1992bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 1993bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 19943a40ed3dSBarry Smith PetscFunctionReturn(0); 199517ab2063SBarry Smith } 199617ab2063SBarry Smith 19975615d1e5SSatish Balay #undef __FUNC__ 19985615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ" 199919bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 200017ab2063SBarry Smith { 2001416022c9SBarry Smith Mat_SeqAIJ *a; 2002416022c9SBarry Smith Mat B; 200317699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 2004bcd2baecSBarry Smith MPI_Comm comm; 200517ab2063SBarry Smith 20063a40ed3dSBarry Smith PetscFunctionBegin; 200719bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 200817699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 2009a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor"); 201090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 20110752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 2012a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file"); 201317ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 201417ab2063SBarry Smith 2015d64ed03dSBarry Smith if (nz < 0) { 2016a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ"); 2017d64ed03dSBarry Smith } 2018d64ed03dSBarry Smith 201917ab2063SBarry Smith /* read in row lengths */ 20200452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 20210752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 202217ab2063SBarry Smith 202317ab2063SBarry Smith /* create our matrix */ 2024416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 2025416022c9SBarry Smith B = *A; 2026416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 2027416022c9SBarry Smith shift = a->indexshift; 202817ab2063SBarry Smith 202917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 20300752156aSBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr); 203117ab2063SBarry Smith if (shift) { 203217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 2033416022c9SBarry Smith a->j[i] += 1; 203417ab2063SBarry Smith } 203517ab2063SBarry Smith } 203617ab2063SBarry Smith 203717ab2063SBarry Smith /* read in nonzero values */ 20380752156aSBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr); 203917ab2063SBarry Smith 204017ab2063SBarry Smith /* set matrix "i" values */ 2041416022c9SBarry Smith a->i[0] = -shift; 204217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 2043416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 2044416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 204517ab2063SBarry Smith } 20460452661fSBarry Smith PetscFree(rowlengths); 204717ab2063SBarry Smith 20486d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 20496d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 20503a40ed3dSBarry Smith PetscFunctionReturn(0); 205117ab2063SBarry Smith } 205217ab2063SBarry Smith 20535615d1e5SSatish Balay #undef __FUNC__ 20545615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ" 20558f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 20567264ac53SSatish Balay { 20577264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 20587264ac53SSatish Balay 20593a40ed3dSBarry Smith PetscFunctionBegin; 2060a8c6a408SBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 20617264ac53SSatish Balay 20627264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 20637264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 2064bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 20653a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2066bcd2baecSBarry Smith } 20677264ac53SSatish Balay 20687264ac53SSatish Balay /* if the a->i are the same */ 20698108c231SLois Curfman McInnes if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) { 20703a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 20717264ac53SSatish Balay } 20727264ac53SSatish Balay 20737264ac53SSatish Balay /* if a->j are the same */ 2074bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 20753a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2076bcd2baecSBarry Smith } 2077bcd2baecSBarry Smith 2078bcd2baecSBarry Smith /* if a->a are the same */ 207919bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 20803a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 20817264ac53SSatish Balay } 208277c4ece6SBarry Smith *flg = PETSC_TRUE; 20833a40ed3dSBarry Smith PetscFunctionReturn(0); 20847264ac53SSatish Balay 20857264ac53SSatish Balay } 2086