1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*02594712SBarry Smith static char vcid[] = "$Id: aij.c,v 1.285 1998/10/19 22:17:52 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" 16eeb667a2SSatish Balay #include "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" 289480ef9eaSBarry 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" 325480ef9eaSBarry 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) 357e20fef11SSatish Balay fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,PetscReal(a->a[j]),PetscImaginary(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) 372e20fef11SSatish Balay if (PetscImaginary(a->a[j]) > 0.0 && PetscReal(a->a[j]) != 0.0) 373e20fef11SSatish Balay fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j])); 374e20fef11SSatish Balay else if (PetscImaginary(a->a[j]) < 0.0 && PetscReal(a->a[j]) != 0.0) 375e20fef11SSatish Balay fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j])); 376e20fef11SSatish Balay else if (PetscReal(a->a[j]) != 0.0) 377e20fef11SSatish Balay fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(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 } 384*02594712SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_SYMMODU) { 385496be53dSLois Curfman McInnes int nzd=0, fshift=1, *sptr; 3862e44a96cSLois Curfman McInnes sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr); 387496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 388496be53dSLois Curfman McInnes sptr[i] = nzd+1; 389496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 390496be53dSLois Curfman McInnes if (a->j[j] >= i) { 3913a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 392e20fef11SSatish Balay if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) nzd++; 393496be53dSLois Curfman McInnes #else 394496be53dSLois Curfman McInnes if (a->a[j] != 0.0) nzd++; 395496be53dSLois Curfman McInnes #endif 396496be53dSLois Curfman McInnes } 397496be53dSLois Curfman McInnes } 398496be53dSLois Curfman McInnes } 3992e44a96cSLois Curfman McInnes sptr[m] = nzd+1; 400496be53dSLois Curfman McInnes fprintf(fd," %d %d\n\n",m,nzd); 4012e44a96cSLois Curfman McInnes for ( i=0; i<m+1; i+=6 ) { 4022e44a96cSLois 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]); 4032e44a96cSLois 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]); 4042e44a96cSLois 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]); 4052e44a96cSLois Curfman McInnes else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]); 4062e44a96cSLois Curfman McInnes else if (i<m) fprintf(fd," %d %d\n",sptr[i],sptr[i+1]); 4077272d637SLois Curfman McInnes else fprintf(fd," %d\n",sptr[i]); 408496be53dSLois Curfman McInnes } 409496be53dSLois Curfman McInnes fprintf(fd,"\n"); 410496be53dSLois Curfman McInnes PetscFree(sptr); 411496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 412496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 413496be53dSLois Curfman McInnes if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift); 414496be53dSLois Curfman McInnes } 415496be53dSLois Curfman McInnes fprintf(fd,"\n"); 416496be53dSLois Curfman McInnes } 417496be53dSLois Curfman McInnes fprintf(fd,"\n"); 418496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 419496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 420496be53dSLois Curfman McInnes if (a->j[j] >= i) { 4213a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 422e20fef11SSatish Balay if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) 423e20fef11SSatish Balay fprintf(fd," %18.16e %18.16e ",PetscReal(a->a[j]),PetscImaginary(a->a[j])); 424496be53dSLois Curfman McInnes #else 425496be53dSLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]); 426496be53dSLois Curfman McInnes #endif 427496be53dSLois Curfman McInnes } 428496be53dSLois Curfman McInnes } 429496be53dSLois Curfman McInnes fprintf(fd,"\n"); 430496be53dSLois Curfman McInnes } 431*02594712SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_DENSE) { 432*02594712SBarry Smith int cnt = 0,jcnt; 433*02594712SBarry Smith Scalar value; 434*02594712SBarry Smith 435*02594712SBarry Smith for ( i=0; i<m; i++ ) { 436*02594712SBarry Smith jcnt = 0; 437*02594712SBarry Smith for ( j=0; j<a->n; j++ ) { 438*02594712SBarry Smith if ( jcnt++ < a->i[i+1]-a->i[i+1] && j == a->j[cnt]) { 439*02594712SBarry Smith value = a->a[cnt++]; 440*02594712SBarry Smith } else { 441*02594712SBarry Smith value = 0.0; 442*02594712SBarry Smith } 443*02594712SBarry Smith #if defined(USE_PETSC_COMPLEX) 444*02594712SBarry Smith fprintf(fd," %7.5e+%7.5e i ",PetscReal(value),PetscImaginary(value)); 445*02594712SBarry Smith #else 446*02594712SBarry Smith fprintf(fd," %7.5e ",value); 447*02594712SBarry Smith #endif 448*02594712SBarry Smith } 449*02594712SBarry Smith fprintf(fd,"\n"); 450*02594712SBarry Smith } 4513a40ed3dSBarry Smith } else { 45217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 45317ab2063SBarry Smith fprintf(fd,"row %d:",i); 454416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 4553a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 456e20fef11SSatish Balay if (PetscImaginary(a->a[j]) > 0.0) { 457e20fef11SSatish Balay fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j])); 458e20fef11SSatish Balay } else if (PetscImaginary(a->a[j]) < 0.0) { 459e20fef11SSatish Balay fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j])); 4603a40ed3dSBarry Smith } else { 461e20fef11SSatish Balay fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j])); 46217ab2063SBarry Smith } 46317ab2063SBarry Smith #else 464416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 46517ab2063SBarry Smith #endif 46617ab2063SBarry Smith } 46717ab2063SBarry Smith fprintf(fd,"\n"); 46817ab2063SBarry Smith } 46917ab2063SBarry Smith } 47017ab2063SBarry Smith fflush(fd); 4713a40ed3dSBarry Smith PetscFunctionReturn(0); 472416022c9SBarry Smith } 473416022c9SBarry Smith 4745615d1e5SSatish Balay #undef __FUNC__ 475480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw_Zoom" 476480ef9eaSBarry Smith int MatView_SeqAIJ_Draw_Zoom(Draw draw,void *Aa) 477416022c9SBarry Smith { 478480ef9eaSBarry Smith Mat A = (Mat) Aa; 479416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 480480ef9eaSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,color; 4810513a670SBarry Smith int format; 482480ef9eaSBarry Smith double xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0; 483480ef9eaSBarry Smith Viewer viewer; 484cddf8d76SBarry Smith 4853a40ed3dSBarry Smith PetscFunctionBegin; 486480ef9eaSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 4870513a670SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 48819bcc07fSBarry Smith 489480ef9eaSBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr); 490416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 4910513a670SBarry Smith 4920513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 4930513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 494cddf8d76SBarry Smith color = DRAW_BLUE; 495416022c9SBarry Smith for ( i=0; i<m; i++ ) { 496cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 497416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 498cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 4993a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 500e20fef11SSatish Balay if (PetscReal(a->a[j]) >= 0.) continue; 501cddf8d76SBarry Smith #else 502cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 503cddf8d76SBarry Smith #endif 504480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 505cddf8d76SBarry Smith } 506cddf8d76SBarry Smith } 507cddf8d76SBarry Smith color = DRAW_CYAN; 508cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 509cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 510cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 511cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 512cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 513480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 514cddf8d76SBarry Smith } 515cddf8d76SBarry Smith } 516cddf8d76SBarry Smith color = DRAW_RED; 517cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 518cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 519cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 520cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5213a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 522e20fef11SSatish Balay if (PetscReal(a->a[j]) <= 0.) continue; 523cddf8d76SBarry Smith #else 524cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 525cddf8d76SBarry Smith #endif 526480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 527416022c9SBarry Smith } 528416022c9SBarry Smith } 5290513a670SBarry Smith } else { 5300513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 5310513a670SBarry Smith /* first determine max of all nonzero values */ 5320513a670SBarry Smith int nz = a->nz,count; 5330513a670SBarry Smith Draw popup; 534480ef9eaSBarry Smith double scale; 5350513a670SBarry Smith 5360513a670SBarry Smith for ( i=0; i<nz; i++ ) { 5370513a670SBarry Smith if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 5380513a670SBarry Smith } 539480ef9eaSBarry Smith scale = (245.0 - DRAW_BASIC_COLORS)/maxv; 540522c5e43SBarry Smith ierr = DrawGetPopup(draw,&popup); CHKERRQ(ierr); 5410513a670SBarry Smith ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr); 5420513a670SBarry Smith count = 0; 5430513a670SBarry Smith for ( i=0; i<m; i++ ) { 5440513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 5450513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 5460513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5476d6b6c30SSatish Balay color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count])); 548480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 5490513a670SBarry Smith count++; 5500513a670SBarry Smith } 5510513a670SBarry Smith } 5520513a670SBarry Smith } 553480ef9eaSBarry Smith PetscFunctionReturn(0); 554480ef9eaSBarry Smith } 555cddf8d76SBarry Smith 556480ef9eaSBarry Smith #undef __FUNC__ 557480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw" 558480ef9eaSBarry Smith int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 559480ef9eaSBarry Smith { 560480ef9eaSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 561480ef9eaSBarry Smith int ierr; 562480ef9eaSBarry Smith Draw draw; 563480ef9eaSBarry Smith double xr,yr,xl,yl,h,w; 564480ef9eaSBarry Smith 565480ef9eaSBarry Smith PetscTruth isnull; 566480ef9eaSBarry Smith 567480ef9eaSBarry Smith PetscFunctionBegin; 568480ef9eaSBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 569480ef9eaSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); 570480ef9eaSBarry Smith if (isnull) PetscFunctionReturn(0); 571480ef9eaSBarry Smith 572480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 573480ef9eaSBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 574480ef9eaSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 575cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 576480ef9eaSBarry Smith ierr = DrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A); CHKERRQ(ierr); 577480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5783a40ed3dSBarry Smith PetscFunctionReturn(0); 579416022c9SBarry Smith } 580416022c9SBarry Smith 5815615d1e5SSatish Balay #undef __FUNC__ 582d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ" 583e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer) 584416022c9SBarry Smith { 585416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 586bcd2baecSBarry Smith ViewerType vtype; 587bcd2baecSBarry Smith int ierr; 588416022c9SBarry Smith 5893a40ed3dSBarry Smith PetscFunctionBegin; 590bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 591bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 5923a40ed3dSBarry Smith ierr = ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); CHKERRQ(ierr); 5935cd90555SBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 5943a40ed3dSBarry Smith ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr); 5955cd90555SBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 5963a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr); 5975cd90555SBarry Smith } else if (vtype == DRAW_VIEWER) { 5983a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr); 5995cd90555SBarry Smith } else { 6005cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 60117ab2063SBarry Smith } 6023a40ed3dSBarry Smith PetscFunctionReturn(0); 60317ab2063SBarry Smith } 60419bcc07fSBarry Smith 605c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 6065615d1e5SSatish Balay #undef __FUNC__ 6075615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ" 6088f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 60917ab2063SBarry Smith { 610416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 61141c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 61243ee02c3SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0; 613416022c9SBarry Smith Scalar *aa = a->a, *ap; 61417ab2063SBarry Smith 6153a40ed3dSBarry Smith PetscFunctionBegin; 6163a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 61717ab2063SBarry Smith 61843ee02c3SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 61917ab2063SBarry Smith for ( i=1; i<m; i++ ) { 620416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 62117ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 62294a9d846SBarry Smith rmax = PetscMax(rmax,ailen[i]); 62317ab2063SBarry Smith if (fshift) { 624416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 62517ab2063SBarry Smith N = ailen[i]; 62617ab2063SBarry Smith for ( j=0; j<N; j++ ) { 62717ab2063SBarry Smith ip[j-fshift] = ip[j]; 62817ab2063SBarry Smith ap[j-fshift] = ap[j]; 62917ab2063SBarry Smith } 63017ab2063SBarry Smith } 63117ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 63217ab2063SBarry Smith } 63317ab2063SBarry Smith if (m) { 63417ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 63517ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 63617ab2063SBarry Smith } 63717ab2063SBarry Smith /* reset ilen and imax for each row */ 63817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 63917ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 64017ab2063SBarry Smith } 641416022c9SBarry Smith a->nz = ai[m] + shift; 64217ab2063SBarry Smith 64317ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 644416022c9SBarry Smith if (fshift && a->diag) { 6450452661fSBarry Smith PetscFree(a->diag); 646416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 647416022c9SBarry Smith a->diag = 0; 64817ab2063SBarry Smith } 6494e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 6504e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 6514e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 652b810aeb4SBarry Smith a->reallocs); 65394a9d846SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax); 654dd5f02e7SSatish Balay a->reallocs = 0; 6554e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 6564e220ebcSLois Curfman McInnes 65776dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 65841c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 6593a40ed3dSBarry Smith PetscFunctionReturn(0); 66017ab2063SBarry Smith } 66117ab2063SBarry Smith 6625615d1e5SSatish Balay #undef __FUNC__ 6635615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ" 6648f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A) 66517ab2063SBarry Smith { 666416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6673a40ed3dSBarry Smith 6683a40ed3dSBarry Smith PetscFunctionBegin; 669cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 6703a40ed3dSBarry Smith PetscFunctionReturn(0); 67117ab2063SBarry Smith } 672416022c9SBarry Smith 6735615d1e5SSatish Balay #undef __FUNC__ 6745615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ" 675e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A) 67617ab2063SBarry Smith { 677416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 67882bf6240SBarry Smith int ierr; 679d5d45c9bSBarry Smith 6803a40ed3dSBarry Smith PetscFunctionBegin; 68194d884c6SBarry Smith if (--A->refct > 0) PetscFunctionReturn(0); 68294d884c6SBarry Smith 68394d884c6SBarry Smith if (A->mapping) { 68494d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 68594d884c6SBarry Smith } 68694d884c6SBarry Smith if (A->bmapping) { 68794d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr); 68894d884c6SBarry Smith } 68961b13de0SBarry Smith if (A->rmap) { 69061b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 69161b13de0SBarry Smith } 69261b13de0SBarry Smith if (A->cmap) { 69361b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 69461b13de0SBarry Smith } 6953a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 696e1311b90SBarry Smith PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 69717ab2063SBarry Smith #endif 6980452661fSBarry Smith PetscFree(a->a); 6990452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 7000452661fSBarry Smith if (a->diag) PetscFree(a->diag); 7010452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 7020452661fSBarry Smith if (a->imax) PetscFree(a->imax); 7030452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 70476dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 70582bf6240SBarry Smith if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 706be6bf707SBarry Smith if (a->saved_values) PetscFree(a->saved_values); 7070452661fSBarry Smith PetscFree(a); 708eed86810SBarry Smith 709f2655603SLois Curfman McInnes PLogObjectDestroy(A); 710f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 7113a40ed3dSBarry Smith PetscFunctionReturn(0); 71217ab2063SBarry Smith } 71317ab2063SBarry Smith 7145615d1e5SSatish Balay #undef __FUNC__ 7155615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ" 7168f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A) 71717ab2063SBarry Smith { 7183a40ed3dSBarry Smith PetscFunctionBegin; 7193a40ed3dSBarry Smith PetscFunctionReturn(0); 72017ab2063SBarry Smith } 72117ab2063SBarry Smith 7225615d1e5SSatish Balay #undef __FUNC__ 7235615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ" 7248f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op) 72517ab2063SBarry Smith { 726416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7273a40ed3dSBarry Smith 7283a40ed3dSBarry Smith PetscFunctionBegin; 7296d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 7306d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 7316d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 732219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 7336d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 734c2653b3dSLois Curfman McInnes else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 73596854ed6SLois Curfman McInnes else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 7366d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 7376d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 738219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 7396d4a8577SBarry Smith op == MAT_SYMMETRIC || 7406d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 74190f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 742b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES|| 743b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) 744981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n"); 7453a40ed3dSBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) { 7463a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 7473a40ed3dSBarry Smith } else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 7486d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 7496d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 7506d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 7516d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 7523a40ed3dSBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 7533a40ed3dSBarry Smith PetscFunctionReturn(0); 75417ab2063SBarry Smith } 75517ab2063SBarry Smith 7565615d1e5SSatish Balay #undef __FUNC__ 7575615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ" 7588f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 75917ab2063SBarry Smith { 760416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7613a40ed3dSBarry Smith int i,j, n,shift = a->indexshift,ierr; 76217ab2063SBarry Smith Scalar *x, zero = 0.0; 76317ab2063SBarry Smith 7643a40ed3dSBarry Smith PetscFunctionBegin; 7653a40ed3dSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 766e1311b90SBarry Smith ierr = VecGetArray(v,&x);;CHKERRQ(ierr); 767e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr); 768a8c6a408SBarry Smith if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector"); 769416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 770416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 771416022c9SBarry Smith if (a->j[j]+shift == i) { 772416022c9SBarry Smith x[i] = a->a[j]; 77317ab2063SBarry Smith break; 77417ab2063SBarry Smith } 77517ab2063SBarry Smith } 77617ab2063SBarry Smith } 777e1311b90SBarry Smith ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr); 7783a40ed3dSBarry Smith PetscFunctionReturn(0); 77917ab2063SBarry Smith } 78017ab2063SBarry Smith 78117ab2063SBarry Smith /* -------------------------------------------------------*/ 78217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 78317ab2063SBarry Smith /* -------------------------------------------------------*/ 7845615d1e5SSatish Balay #undef __FUNC__ 7855615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ" 78644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 78717ab2063SBarry Smith { 788416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 78917ab2063SBarry Smith Scalar *x, *y, *v, alpha; 790e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx, shift = a->indexshift; 79117ab2063SBarry Smith 7923a40ed3dSBarry Smith PetscFunctionBegin; 793e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 794e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 795cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 79617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 79717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 798416022c9SBarry Smith idx = a->j + a->i[i] + shift; 799416022c9SBarry Smith v = a->a + a->i[i] + shift; 800416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 80117ab2063SBarry Smith alpha = x[i]; 80217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 80317ab2063SBarry Smith } 804416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 805e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 806e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8073a40ed3dSBarry Smith PetscFunctionReturn(0); 80817ab2063SBarry Smith } 809d5d45c9bSBarry Smith 8105615d1e5SSatish Balay #undef __FUNC__ 8115615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ" 81244cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 81317ab2063SBarry Smith { 814416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 81517ab2063SBarry Smith Scalar *x, *y, *v, alpha; 816e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx,shift = a->indexshift; 81717ab2063SBarry Smith 8183a40ed3dSBarry Smith PetscFunctionBegin; 8192e8a6d31SBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 820e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 821e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 82217ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 82317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 824416022c9SBarry Smith idx = a->j + a->i[i] + shift; 825416022c9SBarry Smith v = a->a + a->i[i] + shift; 826416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 82717ab2063SBarry Smith alpha = x[i]; 82817ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 82917ab2063SBarry Smith } 83090f02eecSBarry Smith PLogFlops(2*a->nz); 831e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 832e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8333a40ed3dSBarry Smith PetscFunctionReturn(0); 83417ab2063SBarry Smith } 83517ab2063SBarry Smith 8365615d1e5SSatish Balay #undef __FUNC__ 8375615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ" 83844cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 83917ab2063SBarry Smith { 840416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 84117ab2063SBarry Smith Scalar *x, *y, *v, sum; 842e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 8432e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTAIJ) 844e36a17ebSSatish Balay int n, i, jrow,j; 845e36a17ebSSatish Balay #endif 84617ab2063SBarry Smith 847fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT) 848fee21e36SBarry Smith #pragma disjoint(*x,*y,*v) 849fee21e36SBarry Smith #endif 850fee21e36SBarry Smith 8513a40ed3dSBarry Smith PetscFunctionBegin; 852e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 853e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 85417ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 855416022c9SBarry Smith idx = a->j; 856d64ed03dSBarry Smith v = a->a; 857416022c9SBarry Smith ii = a->i; 858acc4d009SSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTAIJ) 85929b5ca8aSSatish Balay fortranmultaij_(&m,x,ii,idx+shift,v+shift,y); 8608d195f9aSBarry Smith #else 861d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 862d64ed03dSBarry Smith idx += shift; 86317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8649ea0dfa2SSatish Balay jrow = ii[i]; 8659ea0dfa2SSatish Balay n = ii[i+1] - jrow; 86617ab2063SBarry Smith sum = 0.0; 8679ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8689ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8699ea0dfa2SSatish Balay } 87017ab2063SBarry Smith y[i] = sum; 87117ab2063SBarry Smith } 8728d195f9aSBarry Smith #endif 873416022c9SBarry Smith PLogFlops(2*a->nz - m); 874e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 875e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 8763a40ed3dSBarry Smith PetscFunctionReturn(0); 87717ab2063SBarry Smith } 87817ab2063SBarry Smith 8795615d1e5SSatish Balay #undef __FUNC__ 8805615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ" 88144cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 88217ab2063SBarry Smith { 883416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 88417ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 885e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 8862e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTADDAIJ) 887e36a17ebSSatish Balay int n,i,jrow,j; 888e36a17ebSSatish Balay #endif 8899ea0dfa2SSatish Balay 8903a40ed3dSBarry Smith PetscFunctionBegin; 891e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 892e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 8932e8a6d31SBarry Smith if (zz != yy) { 894e1311b90SBarry Smith ierr = VecGetArray(zz,&z); CHKERRQ(ierr); 8952e8a6d31SBarry Smith } else { 8962e8a6d31SBarry Smith z = y; 8972e8a6d31SBarry Smith } 89817ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 899cddf8d76SBarry Smith idx = a->j; 900d64ed03dSBarry Smith v = a->a; 901cddf8d76SBarry Smith ii = a->i; 9022e40c62eSSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTADDAIJ) 90329b5ca8aSSatish Balay fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z); 90402ab625aSSatish Balay #else 905d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 906d64ed03dSBarry Smith idx += shift; 90717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 9089ea0dfa2SSatish Balay jrow = ii[i]; 9099ea0dfa2SSatish Balay n = ii[i+1] - jrow; 91017ab2063SBarry Smith sum = y[i]; 9119ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 9129ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 9139ea0dfa2SSatish Balay } 91417ab2063SBarry Smith z[i] = sum; 91517ab2063SBarry Smith } 91602ab625aSSatish Balay #endif 917416022c9SBarry Smith PLogFlops(2*a->nz); 918e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 919e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 9202e8a6d31SBarry Smith if (zz != yy) { 921e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr); 9222e8a6d31SBarry Smith } 9233a40ed3dSBarry Smith PetscFunctionReturn(0); 92417ab2063SBarry Smith } 92517ab2063SBarry Smith 92617ab2063SBarry Smith /* 92717ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 92817ab2063SBarry Smith */ 92917ab2063SBarry Smith 9305615d1e5SSatish Balay #undef __FUNC__ 9315615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ" 932416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 93317ab2063SBarry Smith { 934416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 935416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 93617ab2063SBarry Smith 9373a40ed3dSBarry Smith PetscFunctionBegin; 9380452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 939416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 940416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 94135b0346bSBarry Smith diag[i] = a->i[i+1]; 942416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 943416022c9SBarry Smith if (a->j[j]+shift == i) { 94417ab2063SBarry Smith diag[i] = j - shift; 94517ab2063SBarry Smith break; 94617ab2063SBarry Smith } 94717ab2063SBarry Smith } 94817ab2063SBarry Smith } 949416022c9SBarry Smith a->diag = diag; 9503a40ed3dSBarry Smith PetscFunctionReturn(0); 95117ab2063SBarry Smith } 95217ab2063SBarry Smith 9535615d1e5SSatish Balay #undef __FUNC__ 9545615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ" 95544cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 95617ab2063SBarry Smith double fshift,int its,Vec xx) 95717ab2063SBarry Smith { 958416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 959416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 960d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 96117ab2063SBarry Smith 9623a40ed3dSBarry Smith PetscFunctionBegin; 963e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 964fb2e594dSBarry Smith if (xx != bb) { 965e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 966fb2e594dSBarry Smith } else { 967fb2e594dSBarry Smith b = x; 968fb2e594dSBarry Smith } 969fb2e594dSBarry Smith 970d00d2cf4SBarry Smith if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);} 971416022c9SBarry Smith diag = a->diag; 972416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 97317ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 97417ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 97517ab2063SBarry Smith bs = b + shift; 97617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 977416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 978416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 979488ecbafSBarry Smith PLogFlops(2*n-1); 980416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 981416022c9SBarry Smith v = a->a + diag[i] + (!shift); 98217ab2063SBarry Smith sum = b[i]*d/omega; 98317ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 98417ab2063SBarry Smith x[i] = sum; 98517ab2063SBarry Smith } 986cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 987fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 9883a40ed3dSBarry Smith PetscFunctionReturn(0); 98917ab2063SBarry Smith } 99017ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 991a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done"); 9923a40ed3dSBarry Smith } else if (flag & SOR_EISENSTAT) { 99317ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 99417ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 99517ab2063SBarry Smith 99617ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 99717ab2063SBarry Smith 99817ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 99917ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 100017ab2063SBarry Smith is the relaxation factor. 100117ab2063SBarry Smith */ 10020452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 100317ab2063SBarry Smith scale = (2.0/omega) - 1.0; 100417ab2063SBarry Smith 100517ab2063SBarry Smith /* x = (E + U)^{-1} b */ 100617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1007416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1008416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1009488ecbafSBarry Smith PLogFlops(2*n-1); 1010416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1011416022c9SBarry Smith v = a->a + diag[i] + (!shift); 101217ab2063SBarry Smith sum = b[i]; 101317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 101417ab2063SBarry Smith x[i] = omega*(sum/d); 101517ab2063SBarry Smith } 101617ab2063SBarry Smith 101717ab2063SBarry Smith /* t = b - (2*E - D)x */ 1018416022c9SBarry Smith v = a->a; 101917ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 1020488ecbafSBarry Smith PLogFlops(2*m); 1021488ecbafSBarry Smith 102217ab2063SBarry Smith 102317ab2063SBarry Smith /* t = (E + L)^{-1}t */ 1024416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 1025416022c9SBarry Smith diag = a->diag; 102617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1027416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1028416022c9SBarry Smith n = diag[i] - a->i[i]; 1029488ecbafSBarry Smith PLogFlops(2*n-1); 1030416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1031416022c9SBarry Smith v = a->a + a->i[i] + shift; 103217ab2063SBarry Smith sum = t[i]; 103317ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 103417ab2063SBarry Smith t[i] = omega*(sum/d); 103517ab2063SBarry Smith } 103617ab2063SBarry Smith 103717ab2063SBarry Smith /* x = x + t */ 103817ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 10390452661fSBarry Smith PetscFree(t); 1040cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1041fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 10423a40ed3dSBarry Smith PetscFunctionReturn(0); 104317ab2063SBarry Smith } 104417ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 104517ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 104617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1047416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1048416022c9SBarry Smith n = diag[i] - a->i[i]; 1049488ecbafSBarry Smith PLogFlops(2*n-1); 1050416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1051416022c9SBarry Smith v = a->a + a->i[i] + shift; 105217ab2063SBarry Smith sum = b[i]; 105317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 105417ab2063SBarry Smith x[i] = omega*(sum/d); 105517ab2063SBarry Smith } 105617ab2063SBarry Smith xb = x; 10573a40ed3dSBarry Smith } else xb = b; 105817ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 105917ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 106017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1061416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 106217ab2063SBarry Smith } 1063488ecbafSBarry Smith PLogFlops(m); 106417ab2063SBarry Smith } 106517ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 106617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1067416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1068416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1069488ecbafSBarry Smith PLogFlops(2*n-1); 1070416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1071416022c9SBarry Smith v = a->a + diag[i] + (!shift); 107217ab2063SBarry Smith sum = xb[i]; 107317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 107417ab2063SBarry Smith x[i] = omega*(sum/d); 107517ab2063SBarry Smith } 107617ab2063SBarry Smith } 107717ab2063SBarry Smith its--; 107817ab2063SBarry Smith } 107917ab2063SBarry Smith while (its--) { 108017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 108117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1082416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1083416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1084488ecbafSBarry Smith PLogFlops(2*n-1); 1085416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1086416022c9SBarry Smith v = a->a + a->i[i] + shift; 108717ab2063SBarry Smith sum = b[i]; 108817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 10897e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 109017ab2063SBarry Smith } 109117ab2063SBarry Smith } 109217ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 109317ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1094416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1095416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1096488ecbafSBarry Smith PLogFlops(2*n-1); 1097416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1098416022c9SBarry Smith v = a->a + a->i[i] + shift; 109917ab2063SBarry Smith sum = b[i]; 110017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 11017e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 110217ab2063SBarry Smith } 110317ab2063SBarry Smith } 110417ab2063SBarry Smith } 1105cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1106fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 11073a40ed3dSBarry Smith PetscFunctionReturn(0); 110817ab2063SBarry Smith } 110917ab2063SBarry Smith 11105615d1e5SSatish Balay #undef __FUNC__ 11115615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ" 11128f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 111317ab2063SBarry Smith { 1114416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11154e220ebcSLois Curfman McInnes 11163a40ed3dSBarry Smith PetscFunctionBegin; 11174e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 11184e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 11194e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 11204e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 11214e220ebcSLois Curfman McInnes info->block_size = 1.0; 11224e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 11234e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 11244e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 11254e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 11264e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 11274e220ebcSLois Curfman McInnes info->memory = A->mem; 11284e220ebcSLois Curfman McInnes if (A->factor) { 11294e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 11304e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 11314e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 11324e220ebcSLois Curfman McInnes } else { 11334e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 11344e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 11354e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 11364e220ebcSLois Curfman McInnes } 11373a40ed3dSBarry Smith PetscFunctionReturn(0); 113817ab2063SBarry Smith } 113917ab2063SBarry Smith 114017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 114117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 114217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 114317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 114417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 114517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 114617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 114717ab2063SBarry Smith 11485615d1e5SSatish Balay #undef __FUNC__ 11495615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ" 11508f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 115117ab2063SBarry Smith { 1152416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1153416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 115417ab2063SBarry Smith 11553a40ed3dSBarry Smith PetscFunctionBegin; 115677c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 115717ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 115817ab2063SBarry Smith if (diag) { 115917ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1160a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1161416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 1162416022c9SBarry Smith a->ilen[rows[i]] = 1; 1163416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 1164416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 11653a40ed3dSBarry Smith } else { 1166d64ed03dSBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr); 116717ab2063SBarry Smith } 116817ab2063SBarry Smith } 11693a40ed3dSBarry Smith } else { 117017ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1171a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1172416022c9SBarry Smith a->ilen[rows[i]] = 0; 117317ab2063SBarry Smith } 117417ab2063SBarry Smith } 117517ab2063SBarry Smith ISRestoreIndices(is,&rows); 117643a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11773a40ed3dSBarry Smith PetscFunctionReturn(0); 117817ab2063SBarry Smith } 117917ab2063SBarry Smith 11805615d1e5SSatish Balay #undef __FUNC__ 11815615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ" 11828f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 118317ab2063SBarry Smith { 1184416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11853a40ed3dSBarry Smith 11863a40ed3dSBarry Smith PetscFunctionBegin; 11870752156aSBarry Smith if (m) *m = a->m; 11880752156aSBarry Smith if (n) *n = a->n; 11893a40ed3dSBarry Smith PetscFunctionReturn(0); 119017ab2063SBarry Smith } 119117ab2063SBarry Smith 11925615d1e5SSatish Balay #undef __FUNC__ 11935615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ" 11948f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 119517ab2063SBarry Smith { 1196416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11973a40ed3dSBarry Smith 11983a40ed3dSBarry Smith PetscFunctionBegin; 1199416022c9SBarry Smith *m = 0; *n = a->m; 12003a40ed3dSBarry Smith PetscFunctionReturn(0); 120117ab2063SBarry Smith } 1202026e39d0SSatish Balay 12035615d1e5SSatish Balay #undef __FUNC__ 12045615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ" 12054e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 120617ab2063SBarry Smith { 1207416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1208c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 120917ab2063SBarry Smith 12103a40ed3dSBarry Smith PetscFunctionBegin; 1211a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 121217ab2063SBarry Smith 1213416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 1214416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 121517ab2063SBarry Smith if (idx) { 1216416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 12174e093b46SBarry Smith if (*nz && shift) { 12180452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 121917ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 12204e093b46SBarry Smith } else if (*nz) { 12214e093b46SBarry Smith *idx = itmp; 122217ab2063SBarry Smith } 122317ab2063SBarry Smith else *idx = 0; 122417ab2063SBarry Smith } 12253a40ed3dSBarry Smith PetscFunctionReturn(0); 122617ab2063SBarry Smith } 122717ab2063SBarry Smith 12285615d1e5SSatish Balay #undef __FUNC__ 12295615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ" 12304e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 123117ab2063SBarry Smith { 12324e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 12333a40ed3dSBarry Smith 12343a40ed3dSBarry Smith PetscFunctionBegin; 12354e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 12363a40ed3dSBarry Smith PetscFunctionReturn(0); 123717ab2063SBarry Smith } 123817ab2063SBarry Smith 12395615d1e5SSatish Balay #undef __FUNC__ 12405615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ" 12418f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 124217ab2063SBarry Smith { 1243416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1244416022c9SBarry Smith Scalar *v = a->a; 124517ab2063SBarry Smith double sum = 0.0; 1246416022c9SBarry Smith int i, j,shift = a->indexshift; 124717ab2063SBarry Smith 12483a40ed3dSBarry Smith PetscFunctionBegin; 124917ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1250416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 12513a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 1252e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 125317ab2063SBarry Smith #else 125417ab2063SBarry Smith sum += (*v)*(*v); v++; 125517ab2063SBarry Smith #endif 125617ab2063SBarry Smith } 125717ab2063SBarry Smith *norm = sqrt(sum); 12583a40ed3dSBarry Smith } else if (type == NORM_1) { 125917ab2063SBarry Smith double *tmp; 1260416022c9SBarry Smith int *jj = a->j; 126166963ce1SSatish Balay tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp); 1262cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 126317ab2063SBarry Smith *norm = 0.0; 1264416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1265a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 126617ab2063SBarry Smith } 1267416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 126817ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 126917ab2063SBarry Smith } 12700452661fSBarry Smith PetscFree(tmp); 12713a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 127217ab2063SBarry Smith *norm = 0.0; 1273416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1274416022c9SBarry Smith v = a->a + a->i[j] + shift; 127517ab2063SBarry Smith sum = 0.0; 1276416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1277cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 127817ab2063SBarry Smith } 127917ab2063SBarry Smith if (sum > *norm) *norm = sum; 128017ab2063SBarry Smith } 12813a40ed3dSBarry Smith } else { 1282a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 128317ab2063SBarry Smith } 12843a40ed3dSBarry Smith PetscFunctionReturn(0); 128517ab2063SBarry Smith } 128617ab2063SBarry Smith 12875615d1e5SSatish Balay #undef __FUNC__ 12885615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ" 12898f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B) 129017ab2063SBarry Smith { 1291416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1292416022c9SBarry Smith Mat C; 1293416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1294416022c9SBarry Smith int shift = a->indexshift; 1295d5d45c9bSBarry Smith Scalar *array = a->a; 129617ab2063SBarry Smith 12973a40ed3dSBarry Smith PetscFunctionBegin; 1298a8c6a408SBarry Smith if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place"); 12990452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1300cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 130117ab2063SBarry Smith if (shift) { 130217ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 130317ab2063SBarry Smith } 130417ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1305416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 13060452661fSBarry Smith PetscFree(col); 130717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 130817ab2063SBarry Smith len = ai[i+1]-ai[i]; 1309416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 131017ab2063SBarry Smith array += len; aj += len; 131117ab2063SBarry Smith } 131217ab2063SBarry Smith if (shift) { 131317ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 131417ab2063SBarry Smith } 131517ab2063SBarry Smith 13166d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13176d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 131817ab2063SBarry Smith 13193638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1320416022c9SBarry Smith *B = C; 132117ab2063SBarry Smith } else { 1322f830108cSBarry Smith PetscOps *Abops; 13230a6ffc59SBarry Smith MatOps Aops; 1324f830108cSBarry Smith 1325416022c9SBarry Smith /* This isn't really an in-place transpose */ 13260452661fSBarry Smith PetscFree(a->a); 13270452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 13280452661fSBarry Smith if (a->diag) PetscFree(a->diag); 13290452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 13300452661fSBarry Smith if (a->imax) PetscFree(a->imax); 13310452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 13321073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 13330452661fSBarry Smith PetscFree(a); 1334f830108cSBarry Smith 1335488ecbafSBarry Smith 1336488ecbafSBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 1337488ecbafSBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 1338488ecbafSBarry Smith 1339f830108cSBarry Smith /* 1340f830108cSBarry Smith This is horrible, horrible code. We need to keep the 13418f0f457cSSatish Balay the bops and ops Structures, copy everything from C 13428f0f457cSSatish Balay including the function pointers.. 1343f830108cSBarry Smith */ 13448f0f457cSSatish Balay PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps)); 13458f0f457cSSatish Balay PetscMemcpy(A->bops,C->bops,sizeof(PetscOps)); 1346f830108cSBarry Smith Abops = A->bops; 1347f830108cSBarry Smith Aops = A->ops; 1348f09e8eb9SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 1349f830108cSBarry Smith A->bops = Abops; 1350f830108cSBarry Smith A->ops = Aops; 135127a8da17SBarry Smith A->qlist = 0; 1352f830108cSBarry Smith 13530452661fSBarry Smith PetscHeaderDestroy(C); 135417ab2063SBarry Smith } 13553a40ed3dSBarry Smith PetscFunctionReturn(0); 135617ab2063SBarry Smith } 135717ab2063SBarry Smith 13585615d1e5SSatish Balay #undef __FUNC__ 13595615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ" 13608f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 136117ab2063SBarry Smith { 1362416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 136317ab2063SBarry Smith Scalar *l,*r,x,*v; 1364e1311b90SBarry Smith int ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 136517ab2063SBarry Smith 13663a40ed3dSBarry Smith PetscFunctionBegin; 136717ab2063SBarry Smith if (ll) { 13683ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 13693ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 1370e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr); 1371a8c6a408SBarry Smith if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length"); 1372e1311b90SBarry Smith ierr = VecGetArray(ll,&l); CHKERRQ(ierr); 1373416022c9SBarry Smith v = a->a; 137417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 137517ab2063SBarry Smith x = l[i]; 1376416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 137717ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 137817ab2063SBarry Smith } 1379e1311b90SBarry Smith ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr); 138044cd7ae7SLois Curfman McInnes PLogFlops(nz); 138117ab2063SBarry Smith } 138217ab2063SBarry Smith if (rr) { 1383e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr); 1384a8c6a408SBarry Smith if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length"); 1385e1311b90SBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1386416022c9SBarry Smith v = a->a; jj = a->j; 138717ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 138817ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 138917ab2063SBarry Smith } 1390e1311b90SBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 139144cd7ae7SLois Curfman McInnes PLogFlops(nz); 139217ab2063SBarry Smith } 13933a40ed3dSBarry Smith PetscFunctionReturn(0); 139417ab2063SBarry Smith } 139517ab2063SBarry Smith 13965615d1e5SSatish Balay #undef __FUNC__ 13975615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ" 13986a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B) 139917ab2063SBarry Smith { 1400db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 1401d5db1b72SBarry Smith int *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 140299141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1403a2744918SBarry Smith register int sum,lensi; 140499141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 140599141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 140699141d43SSatish Balay Scalar *a_new,*mat_a; 1407416022c9SBarry Smith Mat C; 1408fee21e36SBarry Smith PetscTruth stride; 140917ab2063SBarry Smith 14103a40ed3dSBarry Smith PetscFunctionBegin; 1411d64ed03dSBarry Smith ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr); 1412a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted"); 1413d64ed03dSBarry Smith ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 1414a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted"); 141599141d43SSatish Balay 141617ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 141717ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 141817ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 141917ab2063SBarry Smith 1420fee21e36SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr); 1421fee21e36SBarry Smith ierr = ISStride(iscol,&stride); CHKERRQ(ierr); 1422fee21e36SBarry Smith if (stride && step == 1) { 142302834360SBarry Smith /* special case of contiguous rows */ 142457faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 142502834360SBarry Smith starts = lens + ncols; 142602834360SBarry Smith /* loop over new rows determining lens and starting points */ 142702834360SBarry Smith for (i=0; i<nrows; i++) { 1428a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1429a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 143002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1431d8ced48eSBarry Smith if (aj[k]+shift >= first) { 143202834360SBarry Smith starts[i] = k; 143302834360SBarry Smith break; 143402834360SBarry Smith } 143502834360SBarry Smith } 1436a2744918SBarry Smith sum = 0; 143702834360SBarry Smith while (k < kend) { 1438d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1439a2744918SBarry Smith sum++; 144002834360SBarry Smith } 1441a2744918SBarry Smith lens[i] = sum; 144202834360SBarry Smith } 144302834360SBarry Smith /* create submatrix */ 1444cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 144508480c60SBarry Smith int n_cols,n_rows; 144608480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1447a8c6a408SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); 1448d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 144908480c60SBarry Smith C = *B; 14503a40ed3dSBarry Smith } else { 145102834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 145208480c60SBarry Smith } 1453db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1454db02288aSLois Curfman McInnes 145502834360SBarry Smith /* loop over rows inserting into submatrix */ 1456db02288aSLois Curfman McInnes a_new = c->a; 1457db02288aSLois Curfman McInnes j_new = c->j; 1458db02288aSLois Curfman McInnes i_new = c->i; 1459db02288aSLois Curfman McInnes i_new[0] = -shift; 146002834360SBarry Smith for (i=0; i<nrows; i++) { 1461a2744918SBarry Smith ii = starts[i]; 1462a2744918SBarry Smith lensi = lens[i]; 1463a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1464a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 146502834360SBarry Smith } 1466a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1467a2744918SBarry Smith a_new += lensi; 1468a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1469a2744918SBarry Smith c->ilen[i] = lensi; 147002834360SBarry Smith } 14710452661fSBarry Smith PetscFree(lens); 14723a40ed3dSBarry Smith } else { 147302834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 14740452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 147502834360SBarry Smith ssmap = smap + shift; 147699141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1477cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 147817ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 147902834360SBarry Smith /* determine lens of each row */ 148002834360SBarry Smith for (i=0; i<nrows; i++) { 1481d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 148202834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 148302834360SBarry Smith lens[i] = 0; 148402834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1485d8ced48eSBarry Smith if (ssmap[aj[k]]) { 148602834360SBarry Smith lens[i]++; 148702834360SBarry Smith } 148802834360SBarry Smith } 148902834360SBarry Smith } 149017ab2063SBarry Smith /* Create and fill new matrix */ 1491a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 149299141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 1493a8c6a408SBarry Smith if (c->m != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size"); 149499141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1495a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros"); 149699141d43SSatish Balay } 149799141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 149808480c60SBarry Smith C = *B; 14993a40ed3dSBarry Smith } else { 150002834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 150108480c60SBarry Smith } 150299141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 150317ab2063SBarry Smith for (i=0; i<nrows; i++) { 150499141d43SSatish Balay row = irow[i]; 150599141d43SSatish Balay kstart = ai[row]+shift; 150699141d43SSatish Balay kend = kstart + a->ilen[row]; 150799141d43SSatish Balay mat_i = c->i[i]+shift; 150899141d43SSatish Balay mat_j = c->j + mat_i; 150999141d43SSatish Balay mat_a = c->a + mat_i; 151099141d43SSatish Balay mat_ilen = c->ilen + i; 151117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 151299141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 151399141d43SSatish Balay *mat_j++ = tcol - (!shift); 151499141d43SSatish Balay *mat_a++ = a->a[k]; 151599141d43SSatish Balay (*mat_ilen)++; 151699141d43SSatish Balay 151717ab2063SBarry Smith } 151817ab2063SBarry Smith } 151917ab2063SBarry Smith } 152002834360SBarry Smith /* Free work space */ 152102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 152299141d43SSatish Balay PetscFree(smap); PetscFree(lens); 152302834360SBarry Smith } 15246d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 15256d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 152617ab2063SBarry Smith 152717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1528416022c9SBarry Smith *B = C; 15293a40ed3dSBarry Smith PetscFunctionReturn(0); 153017ab2063SBarry Smith } 153117ab2063SBarry Smith 1532a871dcd8SBarry Smith /* 153363b91edcSBarry Smith note: This can only work for identity for row and col. It would 153463b91edcSBarry Smith be good to check this and otherwise generate an error. 1535a871dcd8SBarry Smith */ 15365615d1e5SSatish Balay #undef __FUNC__ 15375615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ" 15388f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1539a871dcd8SBarry Smith { 154063b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 154108480c60SBarry Smith int ierr; 154263b91edcSBarry Smith Mat outA; 154363b91edcSBarry Smith 15443a40ed3dSBarry Smith PetscFunctionBegin; 1545a8c6a408SBarry Smith if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 1546a871dcd8SBarry Smith 154763b91edcSBarry Smith outA = inA; 154863b91edcSBarry Smith inA->factor = FACTOR_LU; 154963b91edcSBarry Smith a->row = row; 155063b91edcSBarry Smith a->col = col; 155163b91edcSBarry Smith 1552f0ec6fceSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1553f0ec6fceSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 15541e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1555f0ec6fceSSatish Balay 155694a9d846SBarry Smith if (!a->solve_work) { /* this matrix may have been factored before */ 15570452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work); 155894a9d846SBarry Smith } 155963b91edcSBarry Smith 156008480c60SBarry Smith if (!a->diag) { 156108480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 156263b91edcSBarry Smith } 156363b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 15643a40ed3dSBarry Smith PetscFunctionReturn(0); 1565a871dcd8SBarry Smith } 1566a871dcd8SBarry Smith 156774cdf0dfSBarry Smith #include "pinclude/blaslapack.h" 15685615d1e5SSatish Balay #undef __FUNC__ 15695615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ" 15708f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1571f0b747eeSBarry Smith { 1572f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1573f0b747eeSBarry Smith int one = 1; 15743a40ed3dSBarry Smith 15753a40ed3dSBarry Smith PetscFunctionBegin; 1576f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1577f0b747eeSBarry Smith PLogFlops(a->nz); 15783a40ed3dSBarry Smith PetscFunctionReturn(0); 1579f0b747eeSBarry Smith } 1580f0b747eeSBarry Smith 15815615d1e5SSatish Balay #undef __FUNC__ 15825615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ" 15836a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B) 1584cddf8d76SBarry Smith { 1585cddf8d76SBarry Smith int ierr,i; 1586cddf8d76SBarry Smith 15873a40ed3dSBarry Smith PetscFunctionBegin; 1588cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 15890452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1590cddf8d76SBarry Smith } 1591cddf8d76SBarry Smith 1592cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 15936a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1594cddf8d76SBarry Smith } 15953a40ed3dSBarry Smith PetscFunctionReturn(0); 1596cddf8d76SBarry Smith } 1597cddf8d76SBarry Smith 15985615d1e5SSatish Balay #undef __FUNC__ 15995615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ" 16008f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 16015a838052SSatish Balay { 1602f830108cSBarry Smith PetscFunctionBegin; 16035a838052SSatish Balay *bs = 1; 16043a40ed3dSBarry Smith PetscFunctionReturn(0); 16055a838052SSatish Balay } 16065a838052SSatish Balay 16075615d1e5SSatish Balay #undef __FUNC__ 16085615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ" 16098f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 16104dcbc457SBarry Smith { 1611e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 161206763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 16138a047759SSatish Balay int start, end, *ai, *aj; 1614bbd702dbSSatish Balay BT table; 1615bbd702dbSSatish Balay 16163a40ed3dSBarry Smith PetscFunctionBegin; 16178a047759SSatish Balay shift = a->indexshift; 1618e4d965acSSatish Balay m = a->m; 1619e4d965acSSatish Balay ai = a->i; 16208a047759SSatish Balay aj = a->j+shift; 16218a047759SSatish Balay 1622a8c6a408SBarry Smith if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used"); 162306763907SSatish Balay 162406763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 1625bbd702dbSSatish Balay ierr = BTCreate(m,table); CHKERRQ(ierr); 162606763907SSatish Balay 1627e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1628b97fc60eSLois Curfman McInnes /* Initialize the two local arrays */ 1629e4d965acSSatish Balay isz = 0; 1630bbd702dbSSatish Balay BTMemzero(m,table); 1631e4d965acSSatish Balay 1632e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 16334dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 163477c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1635e4d965acSSatish Balay 1636dd097bc3SLois Curfman McInnes /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 1637e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 1638bbd702dbSSatish Balay if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];} 16394dcbc457SBarry Smith } 164006763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 164106763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1642e4d965acSSatish Balay 164304a348a9SBarry Smith k = 0; 164404a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap */ 164504a348a9SBarry Smith n = isz; 164606763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1647e4d965acSSatish Balay row = nidx[k]; 1648e4d965acSSatish Balay start = ai[row]; 1649e4d965acSSatish Balay end = ai[row+1]; 165004a348a9SBarry Smith for ( l = start; l<end ; l++){ 16518a047759SSatish Balay val = aj[l] + shift; 16522a8f2162SSatish Balay if (!BTLookupSet(table,val)) {nidx[isz++] = val;} 1653e4d965acSSatish Balay } 1654e4d965acSSatish Balay } 1655e4d965acSSatish Balay } 1656029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1657e4d965acSSatish Balay } 1658bbd702dbSSatish Balay BTDestroy(table); 165906763907SSatish Balay PetscFree(nidx); 16603a40ed3dSBarry Smith PetscFunctionReturn(0); 16614dcbc457SBarry Smith } 166217ab2063SBarry Smith 16630513a670SBarry Smith /* -------------------------------------------------------------- */ 16640513a670SBarry Smith #undef __FUNC__ 16650513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ" 16660513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 16670513a670SBarry Smith { 16680513a670SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 16690513a670SBarry Smith Scalar *vwork; 16700513a670SBarry Smith int i, ierr, nz, m = a->m, n = a->n, *cwork; 16710513a670SBarry Smith int *row,*col,*cnew,j,*lens; 167256cd22aeSBarry Smith IS icolp,irowp; 16730513a670SBarry Smith 16743a40ed3dSBarry Smith PetscFunctionBegin; 167556cd22aeSBarry Smith ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr); 167656cd22aeSBarry Smith ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr); 167756cd22aeSBarry Smith ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr); 167856cd22aeSBarry Smith ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr); 16790513a670SBarry Smith 16800513a670SBarry Smith /* determine lengths of permuted rows */ 16810513a670SBarry Smith lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens); 16820513a670SBarry Smith for (i=0; i<m; i++ ) { 16830513a670SBarry Smith lens[row[i]] = a->i[i+1] - a->i[i]; 16840513a670SBarry Smith } 16850513a670SBarry Smith ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr); 16860513a670SBarry Smith PetscFree(lens); 16870513a670SBarry Smith 16880513a670SBarry Smith cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew); 16890513a670SBarry Smith for (i=0; i<m; i++) { 16900513a670SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 16910513a670SBarry Smith for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];} 16920513a670SBarry Smith ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr); 16930513a670SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 16940513a670SBarry Smith } 16950513a670SBarry Smith PetscFree(cnew); 16960513a670SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 16970513a670SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 169856cd22aeSBarry Smith ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr); 169956cd22aeSBarry Smith ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr); 170056cd22aeSBarry Smith ierr = ISDestroy(irowp); CHKERRQ(ierr); 170156cd22aeSBarry Smith ierr = ISDestroy(icolp); CHKERRQ(ierr); 17023a40ed3dSBarry Smith PetscFunctionReturn(0); 17030513a670SBarry Smith } 17040513a670SBarry Smith 17055615d1e5SSatish Balay #undef __FUNC__ 17065615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ" 1707682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1708682d7d0cSBarry Smith { 1709682d7d0cSBarry Smith static int called = 0; 1710682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1711682d7d0cSBarry Smith 17123a40ed3dSBarry Smith PetscFunctionBegin; 17133a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 171476be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 171576be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n"); 171676be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n"); 171776be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_no_inode: Do not use inodes\n"); 171876be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n"); 1719682d7d0cSBarry Smith #if defined(HAVE_ESSL) 172076be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.\n"); 1721682d7d0cSBarry Smith #endif 17223a40ed3dSBarry Smith PetscFunctionReturn(0); 1723682d7d0cSBarry Smith } 17248f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1725a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1726a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1727a93ec695SBarry Smith 1728cb5b572fSBarry Smith #undef __FUNC__ 1729cb5b572fSBarry Smith #define __FUNC__ "MatCopy_SeqAIJ" 1730cb5b572fSBarry Smith int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str) 1731cb5b572fSBarry Smith { 1732be6bf707SBarry Smith int ierr; 1733cb5b572fSBarry Smith 1734cb5b572fSBarry Smith PetscFunctionBegin; 1735be6bf707SBarry Smith if (str == SAME_NONZERO_PATTERN && B->type == MATSEQAIJ) { 1736be6bf707SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1737be6bf707SBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data; 1738be6bf707SBarry Smith 1739be6bf707SBarry Smith if (a->nonew != 1) SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1740be6bf707SBarry Smith if (b->nonew != 1) SETERRQ(1,1,"Must call MatSetOption(B,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1741be6bf707SBarry Smith if (a->i[a->m]+a->indexshift != b->i[b->m]+a->indexshift) { 1742be6bf707SBarry Smith SETERRQ(1,1,"Number of nonzeros in two matrices are different"); 1743cb5b572fSBarry Smith } 1744be6bf707SBarry Smith PetscMemcpy(b->a,a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 1745cb5b572fSBarry Smith } else { 1746cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1747cb5b572fSBarry Smith } 1748cb5b572fSBarry Smith PetscFunctionReturn(0); 1749cb5b572fSBarry Smith } 1750cb5b572fSBarry Smith 1751682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 17520a6ffc59SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ, 1753cb5b572fSBarry Smith MatGetRow_SeqAIJ, 1754cb5b572fSBarry Smith MatRestoreRow_SeqAIJ, 1755cb5b572fSBarry Smith MatMult_SeqAIJ, 1756cb5b572fSBarry Smith MatMultAdd_SeqAIJ, 1757cb5b572fSBarry Smith MatMultTrans_SeqAIJ, 1758cb5b572fSBarry Smith MatMultTransAdd_SeqAIJ, 1759cb5b572fSBarry Smith MatSolve_SeqAIJ, 1760cb5b572fSBarry Smith MatSolveAdd_SeqAIJ, 1761cb5b572fSBarry Smith MatSolveTrans_SeqAIJ, 1762cb5b572fSBarry Smith MatSolveTransAdd_SeqAIJ, 1763cb5b572fSBarry Smith MatLUFactor_SeqAIJ, 1764cb5b572fSBarry Smith 0, 176517ab2063SBarry Smith MatRelax_SeqAIJ, 176617ab2063SBarry Smith MatTranspose_SeqAIJ, 1767cb5b572fSBarry Smith MatGetInfo_SeqAIJ, 1768cb5b572fSBarry Smith MatEqual_SeqAIJ, 1769cb5b572fSBarry Smith MatGetDiagonal_SeqAIJ, 1770cb5b572fSBarry Smith MatDiagonalScale_SeqAIJ, 1771cb5b572fSBarry Smith MatNorm_SeqAIJ, 1772cb5b572fSBarry Smith 0, 1773cb5b572fSBarry Smith MatAssemblyEnd_SeqAIJ, 177417ab2063SBarry Smith MatCompress_SeqAIJ, 1775cb5b572fSBarry Smith MatSetOption_SeqAIJ, 1776cb5b572fSBarry Smith MatZeroEntries_SeqAIJ, 1777cb5b572fSBarry Smith MatZeroRows_SeqAIJ, 1778cb5b572fSBarry Smith MatLUFactorSymbolic_SeqAIJ, 1779cb5b572fSBarry Smith MatLUFactorNumeric_SeqAIJ, 1780cb5b572fSBarry Smith 0, 1781cb5b572fSBarry Smith 0, 1782cb5b572fSBarry Smith MatGetSize_SeqAIJ, 1783cb5b572fSBarry Smith MatGetSize_SeqAIJ, 1784cb5b572fSBarry Smith MatGetOwnershipRange_SeqAIJ, 1785cb5b572fSBarry Smith MatILUFactorSymbolic_SeqAIJ, 1786cb5b572fSBarry Smith 0, 1787cb5b572fSBarry Smith 0, 1788cb5b572fSBarry Smith 0, 1789be6bf707SBarry Smith MatDuplicate_SeqAIJ, 1790cb5b572fSBarry Smith 0, 1791cb5b572fSBarry Smith 0, 1792cb5b572fSBarry Smith MatILUFactor_SeqAIJ, 1793cb5b572fSBarry Smith 0, 1794cb5b572fSBarry Smith 0, 1795cb5b572fSBarry Smith MatGetSubMatrices_SeqAIJ, 1796cb5b572fSBarry Smith MatIncreaseOverlap_SeqAIJ, 1797cb5b572fSBarry Smith MatGetValues_SeqAIJ, 1798cb5b572fSBarry Smith MatCopy_SeqAIJ, 1799f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1800cb5b572fSBarry Smith MatScale_SeqAIJ, 1801cb5b572fSBarry Smith 0, 1802cb5b572fSBarry Smith 0, 18036945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 18046945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 18053b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 18063b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 18073b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1808a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1809a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 18100513a670SBarry Smith MatColoringPatch_SeqAIJ, 18110513a670SBarry Smith 0, 1812cda55fadSBarry Smith MatPermute_SeqAIJ, 1813cda55fadSBarry Smith 0, 1814cda55fadSBarry Smith 0, 1815cda55fadSBarry Smith 0, 1816cda55fadSBarry Smith 0, 1817cda55fadSBarry Smith MatGetMaps_Petsc}; 181817ab2063SBarry Smith 181917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 182017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 182117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 182217ab2063SBarry Smith 1823fb2e594dSBarry Smith EXTERN_C_BEGIN 18245615d1e5SSatish Balay #undef __FUNC__ 1825bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ" 1826bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices) 1827bef8e0ddSBarry Smith { 1828bef8e0ddSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1829bef8e0ddSBarry Smith int i,nz,n; 1830bef8e0ddSBarry Smith 1831bef8e0ddSBarry Smith PetscFunctionBegin; 1832bef8e0ddSBarry Smith if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing"); 1833bef8e0ddSBarry Smith 1834bef8e0ddSBarry Smith nz = aij->maxnz; 1835bef8e0ddSBarry Smith n = aij->n; 1836bef8e0ddSBarry Smith for (i=0; i<nz; i++) { 1837bef8e0ddSBarry Smith aij->j[i] = indices[i]; 1838bef8e0ddSBarry Smith } 1839bef8e0ddSBarry Smith aij->nz = nz; 1840bef8e0ddSBarry Smith for ( i=0; i<n; i++ ) { 1841bef8e0ddSBarry Smith aij->ilen[i] = aij->imax[i]; 1842bef8e0ddSBarry Smith } 1843bef8e0ddSBarry Smith 1844bef8e0ddSBarry Smith PetscFunctionReturn(0); 1845bef8e0ddSBarry Smith } 1846fb2e594dSBarry Smith EXTERN_C_END 1847bef8e0ddSBarry Smith 1848bef8e0ddSBarry Smith #undef __FUNC__ 1849bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices" 1850bef8e0ddSBarry Smith /*@ 1851bef8e0ddSBarry Smith MatSeqAIJSetColumnIndices - Set the column indices for all the rows 1852bef8e0ddSBarry Smith in the matrix. 1853bef8e0ddSBarry Smith 1854bef8e0ddSBarry Smith Input Parameters: 1855bef8e0ddSBarry Smith + mat - the SeqAIJ matrix 1856bef8e0ddSBarry Smith - indices - the column indices 1857bef8e0ddSBarry Smith 1858bef8e0ddSBarry Smith Notes: 1859bef8e0ddSBarry Smith This can be called if you have precomputed the nonzero structure of the 1860bef8e0ddSBarry Smith matrix and want to provide it to the matrix object to improve the performance 1861bef8e0ddSBarry Smith of the MatSetValues() operation. 1862bef8e0ddSBarry Smith 1863bef8e0ddSBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 1864bef8e0ddSBarry Smith MatCreateSeqAIJ(). 1865bef8e0ddSBarry Smith 1866bef8e0ddSBarry Smith MUST be called before any calls to MatSetValues(); 1867bef8e0ddSBarry Smith 1868bef8e0ddSBarry Smith @*/ 1869bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices) 1870bef8e0ddSBarry Smith { 1871bef8e0ddSBarry Smith int ierr,(*f)(Mat,int *); 1872bef8e0ddSBarry Smith 1873bef8e0ddSBarry Smith PetscFunctionBegin; 1874bef8e0ddSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1875bef8e0ddSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f); 1876bef8e0ddSBarry Smith CHKERRQ(ierr); 1877bef8e0ddSBarry Smith if (f) { 1878bef8e0ddSBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 1879bef8e0ddSBarry Smith } else { 1880bef8e0ddSBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 1881bef8e0ddSBarry Smith } 1882bef8e0ddSBarry Smith PetscFunctionReturn(0); 1883bef8e0ddSBarry Smith } 1884bef8e0ddSBarry Smith 1885be6bf707SBarry Smith /* ----------------------------------------------------------------------------------------*/ 1886be6bf707SBarry Smith 1887fb2e594dSBarry Smith EXTERN_C_BEGIN 1888be6bf707SBarry Smith #undef __FUNC__ 1889be6bf707SBarry Smith #define __FUNC__ "MatStoreValues_SeqAIJ" 1890be6bf707SBarry Smith int MatStoreValues_SeqAIJ(Mat mat) 1891be6bf707SBarry Smith { 1892be6bf707SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1893be6bf707SBarry Smith int nz = aij->i[aij->m]+aij->indexshift; 1894be6bf707SBarry Smith 1895be6bf707SBarry Smith PetscFunctionBegin; 1896be6bf707SBarry Smith if (aij->nonew != 1) { 1897be6bf707SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1898be6bf707SBarry Smith } 1899be6bf707SBarry Smith 1900be6bf707SBarry Smith /* allocate space for values if not already there */ 1901be6bf707SBarry Smith if (!aij->saved_values) { 1902be6bf707SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 1903be6bf707SBarry Smith } 1904be6bf707SBarry Smith 1905be6bf707SBarry Smith /* copy values over */ 1906be6bf707SBarry Smith PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar)); 1907be6bf707SBarry Smith PetscFunctionReturn(0); 1908be6bf707SBarry Smith } 1909fb2e594dSBarry Smith EXTERN_C_END 1910be6bf707SBarry Smith 1911be6bf707SBarry Smith #undef __FUNC__ 1912be6bf707SBarry Smith #define __FUNC__ "MatStoreValues" 1913be6bf707SBarry Smith /*@ 1914be6bf707SBarry Smith MatStoreValues - Stashes a copy of the matrix values; this allows, for 1915be6bf707SBarry Smith example, reuse of the linear part of a Jacobian, while recomputing the 1916be6bf707SBarry Smith nonlinear portion. 1917be6bf707SBarry Smith 1918be6bf707SBarry Smith Collect on Mat 1919be6bf707SBarry Smith 1920be6bf707SBarry Smith Input Parameters: 1921be6bf707SBarry Smith . mat - the matrix (currently on AIJ matrices support this option) 1922be6bf707SBarry Smith 1923be6bf707SBarry Smith Common Usage, with SNESSolve(): 1924be6bf707SBarry Smith $ Create Jacobian matrix 1925be6bf707SBarry Smith $ Set linear terms into matrix 1926be6bf707SBarry Smith $ Apply boundary conditions to matrix, at this time matrix must have 1927be6bf707SBarry Smith $ final nonzero structure (i.e. setting the nonlinear terms and applying 1928be6bf707SBarry Smith $ boundary conditions again will not change the nonzero structure 1929be6bf707SBarry Smith $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); 1930be6bf707SBarry Smith $ ierr = MatStoreValues(mat); 1931be6bf707SBarry Smith $ Call SNESSetJacobian() with matrix 1932be6bf707SBarry Smith $ In your Jacobian routine 1933be6bf707SBarry Smith $ ierr = MatRetrieveValues(mat); 1934be6bf707SBarry Smith $ Set nonlinear terms in matrix 1935be6bf707SBarry Smith 1936be6bf707SBarry Smith Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself: 1937be6bf707SBarry Smith $ // build linear portion of Jacobian 1938be6bf707SBarry Smith $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); 1939be6bf707SBarry Smith $ ierr = MatStoreValues(mat); 1940be6bf707SBarry Smith $ loop over nonlinear iterations 1941be6bf707SBarry Smith $ ierr = MatRetrieveValues(mat); 1942be6bf707SBarry Smith $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian 1943be6bf707SBarry Smith $ // call MatAssemblyBegin/End() on matrix 1944be6bf707SBarry Smith $ Solve linear system with Jacobian 1945be6bf707SBarry Smith $ endloop 1946be6bf707SBarry Smith 1947be6bf707SBarry Smith Notes: 1948be6bf707SBarry Smith Matrix must already be assemblied before calling this routine 1949be6bf707SBarry Smith Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 1950be6bf707SBarry Smith calling this routine. 1951be6bf707SBarry Smith 1952be6bf707SBarry Smith .seealso: MatRetrieveValues() 1953be6bf707SBarry Smith 1954be6bf707SBarry Smith @*/ 1955be6bf707SBarry Smith int MatStoreValues(Mat mat) 1956be6bf707SBarry Smith { 1957be6bf707SBarry Smith int ierr,(*f)(Mat); 1958be6bf707SBarry Smith 1959be6bf707SBarry Smith PetscFunctionBegin; 1960be6bf707SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1961be6bf707SBarry Smith if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1962be6bf707SBarry Smith if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1963be6bf707SBarry Smith 1964be6bf707SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void **)&f); 1965be6bf707SBarry Smith CHKERRQ(ierr); 1966be6bf707SBarry Smith if (f) { 1967be6bf707SBarry Smith ierr = (*f)(mat);CHKERRQ(ierr); 1968be6bf707SBarry Smith } else { 1969be6bf707SBarry Smith SETERRQ(1,1,"Wrong type of matrix to store values"); 1970be6bf707SBarry Smith } 1971be6bf707SBarry Smith PetscFunctionReturn(0); 1972be6bf707SBarry Smith } 1973be6bf707SBarry Smith 1974fb2e594dSBarry Smith EXTERN_C_BEGIN 1975be6bf707SBarry Smith #undef __FUNC__ 1976be6bf707SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqAIJ" 1977be6bf707SBarry Smith int MatRetrieveValues_SeqAIJ(Mat mat) 1978be6bf707SBarry Smith { 1979be6bf707SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1980be6bf707SBarry Smith int nz = aij->i[aij->m]+aij->indexshift; 1981be6bf707SBarry Smith 1982be6bf707SBarry Smith PetscFunctionBegin; 1983be6bf707SBarry Smith if (aij->nonew != 1) { 1984be6bf707SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1985be6bf707SBarry Smith } 1986be6bf707SBarry Smith if (!aij->saved_values) { 1987be6bf707SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 1988be6bf707SBarry Smith } 1989be6bf707SBarry Smith 1990be6bf707SBarry Smith /* copy values over */ 1991be6bf707SBarry Smith PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar)); 1992be6bf707SBarry Smith PetscFunctionReturn(0); 1993be6bf707SBarry Smith } 1994fb2e594dSBarry Smith EXTERN_C_END 1995be6bf707SBarry Smith 1996be6bf707SBarry Smith #undef __FUNC__ 1997be6bf707SBarry Smith #define __FUNC__ "MatRetrieveValues" 1998be6bf707SBarry Smith /*@ 1999be6bf707SBarry Smith MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 2000be6bf707SBarry Smith example, reuse of the linear part of a Jacobian, while recomputing the 2001be6bf707SBarry Smith nonlinear portion. 2002be6bf707SBarry Smith 2003be6bf707SBarry Smith Collect on Mat 2004be6bf707SBarry Smith 2005be6bf707SBarry Smith Input Parameters: 2006be6bf707SBarry Smith . mat - the matrix (currently on AIJ matrices support this option) 2007be6bf707SBarry Smith 2008be6bf707SBarry Smith .seealso: MatStoreValues() 2009be6bf707SBarry Smith 2010be6bf707SBarry Smith @*/ 2011be6bf707SBarry Smith int MatRetrieveValues(Mat mat) 2012be6bf707SBarry Smith { 2013be6bf707SBarry Smith int ierr,(*f)(Mat); 2014be6bf707SBarry Smith 2015be6bf707SBarry Smith PetscFunctionBegin; 2016be6bf707SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 2017be6bf707SBarry Smith if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2018be6bf707SBarry Smith if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2019be6bf707SBarry Smith 2020be6bf707SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void **)&f); 2021be6bf707SBarry Smith CHKERRQ(ierr); 2022be6bf707SBarry Smith if (f) { 2023be6bf707SBarry Smith ierr = (*f)(mat);CHKERRQ(ierr); 2024be6bf707SBarry Smith } else { 2025be6bf707SBarry Smith SETERRQ(1,1,"Wrong type of matrix to retrieve values"); 2026be6bf707SBarry Smith } 2027be6bf707SBarry Smith PetscFunctionReturn(0); 2028be6bf707SBarry Smith } 2029be6bf707SBarry Smith 2030be6bf707SBarry Smith /* --------------------------------------------------------------------------------*/ 2031cb5b572fSBarry Smith 2032bef8e0ddSBarry Smith #undef __FUNC__ 20335615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ" 203417ab2063SBarry Smith /*@C 2035682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 20360d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 20376e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 20382bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 20392bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 204017ab2063SBarry Smith 2041db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2042db81eaa0SLois Curfman McInnes 204317ab2063SBarry Smith Input Parameters: 2044db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 204517ab2063SBarry Smith . m - number of rows 204617ab2063SBarry Smith . n - number of columns 204717ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 2048db81eaa0SLois Curfman McInnes - nzz - array containing the number of nonzeros in the various rows 20492bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 205017ab2063SBarry Smith 205117ab2063SBarry Smith Output Parameter: 2052416022c9SBarry Smith . A - the matrix 205317ab2063SBarry Smith 2054b259b22eSLois Curfman McInnes Notes: 205517ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 205617ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 20570002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 205844cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 205917ab2063SBarry Smith 206017ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2061a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 20623d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 20636da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 206417ab2063SBarry Smith 2065682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 20664fca80b9SLois Curfman McInnes improve numerical efficiency of matrix-vector products and solves. We 2067682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 20686c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 20696c7ebb05SLois Curfman McInnes 20706c7ebb05SLois Curfman McInnes Options Database Keys: 2071db81eaa0SLois Curfman McInnes + -mat_aij_no_inode - Do not use inodes 2072db81eaa0SLois Curfman McInnes . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5) 2073db81eaa0SLois Curfman McInnes - -mat_aij_oneindex - Internally use indexing starting at 1 2074db81eaa0SLois Curfman McInnes rather than 0. Note that when calling MatSetValues(), 2075db81eaa0SLois Curfman McInnes the user still MUST index entries starting at 0! 207617ab2063SBarry Smith 2077bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices() 207817ab2063SBarry Smith @*/ 2079416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 208017ab2063SBarry Smith { 2081416022c9SBarry Smith Mat B; 2082416022c9SBarry Smith Mat_SeqAIJ *b; 20836945ee14SBarry Smith int i, len, ierr, flg,size; 20846945ee14SBarry Smith 20853a40ed3dSBarry Smith PetscFunctionBegin; 20866945ee14SBarry Smith MPI_Comm_size(comm,&size); 2087a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1"); 2088d5d45c9bSBarry Smith 2089416022c9SBarry Smith *A = 0; 2090f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,comm,MatDestroy,MatView); 2091416022c9SBarry Smith PLogObjectCreate(B); 20920452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 209344cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 20940a6ffc59SBarry Smith PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)); 2095e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqAIJ; 2096e1311b90SBarry Smith B->ops->view = MatView_SeqAIJ; 2097416022c9SBarry Smith B->factor = 0; 2098416022c9SBarry Smith B->lupivotthreshold = 1.0; 209990f02eecSBarry Smith B->mapping = 0; 2100e1311b90SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg);CHKERRQ(ierr); 21017a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 2102e1311b90SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*)&b->ilu_preserve_row_sums);CHKERRQ(ierr); 2103416022c9SBarry Smith b->row = 0; 2104416022c9SBarry Smith b->col = 0; 210582bf6240SBarry Smith b->icol = 0; 2106416022c9SBarry Smith b->indexshift = 0; 2107b810aeb4SBarry Smith b->reallocs = 0; 210869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 210969957df2SSatish Balay if (flg) b->indexshift = -1; 211017ab2063SBarry Smith 211144cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 211244cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 2113a5ae1ecdSBarry Smith 21144d197716SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 21154d197716SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 2116a5ae1ecdSBarry Smith 21170452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 2118b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 21197b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 21207b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 2121416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 212217ab2063SBarry Smith nz = nz*m; 21233a40ed3dSBarry Smith } else { 212417ab2063SBarry Smith nz = 0; 2125416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 212617ab2063SBarry Smith } 212717ab2063SBarry Smith 212817ab2063SBarry Smith /* allocate the matrix space */ 2129416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 21300452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 2131416022c9SBarry Smith b->j = (int *) (b->a + nz); 2132cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 2133416022c9SBarry Smith b->i = b->j + nz; 2134416022c9SBarry Smith b->singlemalloc = 1; 213517ab2063SBarry Smith 2136416022c9SBarry Smith b->i[0] = -b->indexshift; 213717ab2063SBarry Smith for (i=1; i<m+1; i++) { 2138416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 213917ab2063SBarry Smith } 214017ab2063SBarry Smith 2141416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 21420452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 2143f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 2144416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 214517ab2063SBarry Smith 2146416022c9SBarry Smith b->nz = 0; 2147416022c9SBarry Smith b->maxnz = nz; 2148416022c9SBarry Smith b->sorted = 0; 2149416022c9SBarry Smith b->roworiented = 1; 2150416022c9SBarry Smith b->nonew = 0; 2151416022c9SBarry Smith b->diag = 0; 2152416022c9SBarry Smith b->solve_work = 0; 215371bd300dSLois Curfman McInnes b->spptr = 0; 2154754ec7b1SSatish Balay b->inode.node_count = 0; 2155754ec7b1SSatish Balay b->inode.size = 0; 21566c7ebb05SLois Curfman McInnes b->inode.limit = 5; 21576c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 2158be6bf707SBarry Smith b->saved_values = 0; 21594e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 216017ab2063SBarry Smith 2161416022c9SBarry Smith *A = B; 21624e220ebcSLois Curfman McInnes 21634b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 21644b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 216569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 216669957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 21674b14c69eSBarry Smith #endif 216869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 216969957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 217069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 217169957df2SSatish Balay if (flg) { 2172a8c6a408SBarry Smith if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml"); 2173416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 217417ab2063SBarry Smith } 217569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 217669957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 2177bef8e0ddSBarry Smith 2178bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C", 2179bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 2180bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 2181be6bf707SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 2182be6bf707SBarry Smith "MatStoreValues_SeqAIJ", 2183be6bf707SBarry Smith (void*)MatStoreValues_SeqAIJ);CHKERRQ(ierr); 2184be6bf707SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 2185be6bf707SBarry Smith "MatRetrieveValues_SeqAIJ", 2186be6bf707SBarry Smith (void*)MatRetrieveValues_SeqAIJ);CHKERRQ(ierr); 21873a40ed3dSBarry Smith PetscFunctionReturn(0); 218817ab2063SBarry Smith } 218917ab2063SBarry Smith 21905615d1e5SSatish Balay #undef __FUNC__ 2191be6bf707SBarry Smith #define __FUNC__ "MatDuplicate_SeqAIJ" 2192be6bf707SBarry Smith int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 219317ab2063SBarry Smith { 2194416022c9SBarry Smith Mat C; 2195416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 2196bef8e0ddSBarry Smith int i,len, m = a->m,shift = a->indexshift,ierr; 219717ab2063SBarry Smith 21983a40ed3dSBarry Smith PetscFunctionBegin; 21994043dd9cSLois Curfman McInnes *B = 0; 2200f830108cSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,A->comm,MatDestroy,MatView); 2201416022c9SBarry Smith PLogObjectCreate(C); 22020452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 2203f830108cSBarry Smith PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 2204e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqAIJ; 2205e1311b90SBarry Smith C->ops->view = MatView_SeqAIJ; 2206416022c9SBarry Smith C->factor = A->factor; 2207416022c9SBarry Smith c->row = 0; 2208416022c9SBarry Smith c->col = 0; 220982bf6240SBarry Smith c->icol = 0; 2210416022c9SBarry Smith c->indexshift = shift; 2211c456f294SBarry Smith C->assembled = PETSC_TRUE; 221217ab2063SBarry Smith 221344cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 221444cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 221544cd7ae7SLois Curfman McInnes C->M = a->m; 221644cd7ae7SLois Curfman McInnes C->N = a->n; 221717ab2063SBarry Smith 22180452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 22190452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 222017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 2221416022c9SBarry Smith c->imax[i] = a->imax[i]; 2222416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 222317ab2063SBarry Smith } 222417ab2063SBarry Smith 222517ab2063SBarry Smith /* allocate the matrix space */ 2226416022c9SBarry Smith c->singlemalloc = 1; 2227416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 22280452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 2229416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 2230416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 2231416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 223217ab2063SBarry Smith if (m > 0) { 2233416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 2234be6bf707SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 2235416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 2236be6bf707SBarry Smith } else { 2237be6bf707SBarry Smith PetscMemzero(c->a,(a->i[m]+shift)*sizeof(Scalar)); 223817ab2063SBarry Smith } 223908480c60SBarry Smith } 224017ab2063SBarry Smith 2241f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 2242416022c9SBarry Smith c->sorted = a->sorted; 2243416022c9SBarry Smith c->roworiented = a->roworiented; 2244416022c9SBarry Smith c->nonew = a->nonew; 22457a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 2246be6bf707SBarry Smith c->saved_values = 0; 224717ab2063SBarry Smith 2248416022c9SBarry Smith if (a->diag) { 22490452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 2250416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 225117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 2252416022c9SBarry Smith c->diag[i] = a->diag[i]; 225317ab2063SBarry Smith } 22543a40ed3dSBarry Smith } else c->diag = 0; 22556c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 22566c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 2257754ec7b1SSatish Balay if (a->inode.size){ 2258daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 2259754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 2260daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 2261754ec7b1SSatish Balay } else { 2262754ec7b1SSatish Balay c->inode.size = 0; 2263754ec7b1SSatish Balay c->inode.node_count = 0; 2264754ec7b1SSatish Balay } 2265416022c9SBarry Smith c->nz = a->nz; 2266416022c9SBarry Smith c->maxnz = a->maxnz; 2267416022c9SBarry Smith c->solve_work = 0; 226876dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 2269754ec7b1SSatish Balay 2270416022c9SBarry Smith *B = C; 2271bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C", 2272bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 2273bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 22743a40ed3dSBarry Smith PetscFunctionReturn(0); 227517ab2063SBarry Smith } 227617ab2063SBarry Smith 22775615d1e5SSatish Balay #undef __FUNC__ 22785615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ" 227919bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 228017ab2063SBarry Smith { 2281416022c9SBarry Smith Mat_SeqAIJ *a; 2282416022c9SBarry Smith Mat B; 228317699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 2284bcd2baecSBarry Smith MPI_Comm comm; 228517ab2063SBarry Smith 22863a40ed3dSBarry Smith PetscFunctionBegin; 228719bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 228817699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 2289a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor"); 229090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 22910752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 2292a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file"); 229317ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 229417ab2063SBarry Smith 2295d64ed03dSBarry Smith if (nz < 0) { 2296a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ"); 2297d64ed03dSBarry Smith } 2298d64ed03dSBarry Smith 229917ab2063SBarry Smith /* read in row lengths */ 23000452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 23010752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 230217ab2063SBarry Smith 230317ab2063SBarry Smith /* create our matrix */ 2304416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 2305416022c9SBarry Smith B = *A; 2306416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 2307416022c9SBarry Smith shift = a->indexshift; 230817ab2063SBarry Smith 230917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 23100752156aSBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr); 231117ab2063SBarry Smith if (shift) { 231217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 2313416022c9SBarry Smith a->j[i] += 1; 231417ab2063SBarry Smith } 231517ab2063SBarry Smith } 231617ab2063SBarry Smith 231717ab2063SBarry Smith /* read in nonzero values */ 23180752156aSBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr); 231917ab2063SBarry Smith 232017ab2063SBarry Smith /* set matrix "i" values */ 2321416022c9SBarry Smith a->i[0] = -shift; 232217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 2323416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 2324416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 232517ab2063SBarry Smith } 23260452661fSBarry Smith PetscFree(rowlengths); 232717ab2063SBarry Smith 23286d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 23296d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 23303a40ed3dSBarry Smith PetscFunctionReturn(0); 233117ab2063SBarry Smith } 233217ab2063SBarry Smith 23335615d1e5SSatish Balay #undef __FUNC__ 23345615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ" 23358f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 23367264ac53SSatish Balay { 23377264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 23387264ac53SSatish Balay 23393a40ed3dSBarry Smith PetscFunctionBegin; 2340a8c6a408SBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 23417264ac53SSatish Balay 23427264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 23437264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 2344bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 23453a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2346bcd2baecSBarry Smith } 23477264ac53SSatish Balay 23487264ac53SSatish Balay /* if the a->i are the same */ 23498108c231SLois Curfman McInnes if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) { 23503a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 23517264ac53SSatish Balay } 23527264ac53SSatish Balay 23537264ac53SSatish Balay /* if a->j are the same */ 2354bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 23553a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2356bcd2baecSBarry Smith } 2357bcd2baecSBarry Smith 2358bcd2baecSBarry Smith /* if a->a are the same */ 235919bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 23603a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 23617264ac53SSatish Balay } 236277c4ece6SBarry Smith *flg = PETSC_TRUE; 23633a40ed3dSBarry Smith PetscFunctionReturn(0); 23647264ac53SSatish Balay 23657264ac53SSatish Balay } 2366