1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3f1db9ecSBarry Smith static char vcid[] = "$Id: aij.c,v 1.289 1998/12/03 03:59:57 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 "sys.h" 1170f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 12f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 13f5eb4b81SSatish Balay #include "src/inline/spops.h" 148d195f9aSBarry Smith #include "src/inline/dot.h" 15eeb667a2SSatish Balay #include "bitarray.h" 1617ab2063SBarry Smith 17a2ce50c7SBarry Smith /* 18a2ce50c7SBarry Smith Basic AIJ format ILU based on drop tolerance 19a2ce50c7SBarry Smith */ 205615d1e5SSatish Balay #undef __FUNC__ 215615d1e5SSatish Balay #define __FUNC__ "MatILUDTFactor_SeqAIJ" 22a2ce50c7SBarry Smith int MatILUDTFactor_SeqAIJ(Mat A,double dt,int maxnz,IS row,IS col,Mat *fact) 23a2ce50c7SBarry Smith { 24a2ce50c7SBarry Smith /* Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; */ 25a2ce50c7SBarry Smith 263a40ed3dSBarry Smith PetscFunctionBegin; 27*3f1db9ecSBarry Smith SETERRQ(1,0,"Not implemented"); 28d64ed03dSBarry Smith #if !defined(USE_PETSC_DEBUG) 29d64ed03dSBarry Smith PetscFunctionReturn(0); 30d64ed03dSBarry Smith #endif 31a2ce50c7SBarry Smith } 32a2ce50c7SBarry Smith 33bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 3417ab2063SBarry Smith 355615d1e5SSatish Balay #undef __FUNC__ 365615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ" 378f6be9afSLois Curfman McInnes int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja, 386945ee14SBarry Smith PetscTruth *done) 3917ab2063SBarry Smith { 40416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 416945ee14SBarry Smith int ierr,i,ishift; 4217ab2063SBarry Smith 433a40ed3dSBarry Smith PetscFunctionBegin; 4431625ec5SSatish Balay *m = A->m; 453a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 466945ee14SBarry Smith ishift = a->indexshift; 476945ee14SBarry Smith if (symmetric) { 4831625ec5SSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 496945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 5031625ec5SSatish Balay int nz = a->i[a->m]; 513b2fbd54SBarry Smith /* malloc space and subtract 1 from i and j indices */ 5231625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 533b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 543b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1; 5531625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1; 566945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 5731625ec5SSatish Balay int nz = a->i[a->m] + 1; 583b2fbd54SBarry Smith /* malloc space and add 1 to i and j indices */ 5931625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 603b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 613b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1; 6231625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1; 636945ee14SBarry Smith } else { 646945ee14SBarry Smith *ia = a->i; *ja = a->j; 65a2ce50c7SBarry Smith } 66a2ce50c7SBarry Smith 673a40ed3dSBarry Smith PetscFunctionReturn(0); 68a2744918SBarry Smith } 69a2744918SBarry Smith 705615d1e5SSatish Balay #undef __FUNC__ 715615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ" 728f6be9afSLois Curfman McInnes int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 736945ee14SBarry Smith PetscTruth *done) 746945ee14SBarry Smith { 756945ee14SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 763b2fbd54SBarry Smith int ishift = a->indexshift; 776945ee14SBarry Smith 783a40ed3dSBarry Smith PetscFunctionBegin; 793a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 803b2fbd54SBarry Smith if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) { 816945ee14SBarry Smith PetscFree(*ia); 826945ee14SBarry Smith PetscFree(*ja); 83bcd2baecSBarry Smith } 843a40ed3dSBarry Smith PetscFunctionReturn(0); 8517ab2063SBarry Smith } 8617ab2063SBarry Smith 875615d1e5SSatish Balay #undef __FUNC__ 885615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ" 8943a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 903b2fbd54SBarry Smith PetscTruth *done) 913b2fbd54SBarry Smith { 923b2fbd54SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 93a93ec695SBarry Smith int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m; 94a93ec695SBarry Smith int nz = a->i[m]+ishift,row,*jj,mr,col; 953b2fbd54SBarry Smith 963a40ed3dSBarry Smith PetscFunctionBegin; 973b2fbd54SBarry Smith *nn = A->n; 983a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 993b2fbd54SBarry Smith if (symmetric) { 100179192dfSSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 1013b2fbd54SBarry Smith } else { 10261d2ded1SBarry Smith collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths); 1033b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1043b2fbd54SBarry Smith cia = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia); 105a93ec695SBarry Smith cja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja); 1063b2fbd54SBarry Smith jj = a->j; 1073b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) { 1083b2fbd54SBarry Smith collengths[jj[i] + ishift]++; 1093b2fbd54SBarry Smith } 1103b2fbd54SBarry Smith cia[0] = oshift; 1113b2fbd54SBarry Smith for ( i=0; i<n; i++) { 1123b2fbd54SBarry Smith cia[i+1] = cia[i] + collengths[i]; 1133b2fbd54SBarry Smith } 1143b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1153b2fbd54SBarry Smith jj = a->j; 116a93ec695SBarry Smith for ( row=0; row<m; row++ ) { 117a93ec695SBarry Smith mr = a->i[row+1] - a->i[row]; 118a93ec695SBarry Smith for ( i=0; i<mr; i++ ) { 1193b2fbd54SBarry Smith col = *jj++ + ishift; 1203b2fbd54SBarry Smith cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1213b2fbd54SBarry Smith } 1223b2fbd54SBarry Smith } 1233b2fbd54SBarry Smith PetscFree(collengths); 1243b2fbd54SBarry Smith *ia = cia; *ja = cja; 1253b2fbd54SBarry Smith } 1263b2fbd54SBarry Smith 1273a40ed3dSBarry Smith PetscFunctionReturn(0); 1283b2fbd54SBarry Smith } 1293b2fbd54SBarry Smith 1305615d1e5SSatish Balay #undef __FUNC__ 1315615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ" 13243a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia, 1333b2fbd54SBarry Smith int **ja,PetscTruth *done) 1343b2fbd54SBarry Smith { 1353a40ed3dSBarry Smith PetscFunctionBegin; 1363a40ed3dSBarry Smith if (!ia) PetscFunctionReturn(0); 1373b2fbd54SBarry Smith 1383b2fbd54SBarry Smith PetscFree(*ia); 1393b2fbd54SBarry Smith PetscFree(*ja); 1403b2fbd54SBarry Smith 1413a40ed3dSBarry Smith PetscFunctionReturn(0); 1423b2fbd54SBarry Smith } 1433b2fbd54SBarry Smith 144227d817aSBarry Smith #define CHUNKSIZE 15 14517ab2063SBarry Smith 1465615d1e5SSatish Balay #undef __FUNC__ 1475615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ" 14861d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 14917ab2063SBarry Smith { 150416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 151416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 1524b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 153d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 154416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 15517ab2063SBarry Smith 1563a40ed3dSBarry Smith PetscFunctionBegin; 15717ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 158416022c9SBarry Smith row = im[k]; 1593a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 160a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 161a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 1623b2fbd54SBarry Smith #endif 16317ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 16417ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 165416022c9SBarry Smith low = 0; 16617ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 1673a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 168a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 169a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 1703b2fbd54SBarry Smith #endif 1714b0e389bSBarry Smith col = in[l] - shift; 1724b0e389bSBarry Smith if (roworiented) { 1734b0e389bSBarry Smith value = *v++; 174bef8e0ddSBarry Smith } else { 1754b0e389bSBarry Smith value = v[k + l*m]; 1764b0e389bSBarry Smith } 177416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 178416022c9SBarry Smith while (high-low > 5) { 179416022c9SBarry Smith t = (low+high)/2; 180416022c9SBarry Smith if (rp[t] > col) high = t; 181416022c9SBarry Smith else low = t; 18217ab2063SBarry Smith } 183416022c9SBarry Smith for ( i=low; i<high; i++ ) { 18417ab2063SBarry Smith if (rp[i] > col) break; 18517ab2063SBarry Smith if (rp[i] == col) { 186416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 18717ab2063SBarry Smith else ap[i] = value; 18817ab2063SBarry Smith goto noinsert; 18917ab2063SBarry Smith } 19017ab2063SBarry Smith } 191c2653b3dSLois Curfman McInnes if (nonew == 1) goto noinsert; 192a8c6a408SBarry Smith else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 19317ab2063SBarry Smith if (nrow >= rmax) { 19417ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 195416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 19617ab2063SBarry Smith Scalar *new_a; 19717ab2063SBarry Smith 198a8c6a408SBarry Smith if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix"); 19996854ed6SLois Curfman McInnes 20017ab2063SBarry Smith /* malloc new storage space */ 201416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 2020452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 20317ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 20417ab2063SBarry Smith new_i = new_j + new_nz; 20517ab2063SBarry Smith 20617ab2063SBarry Smith /* copy over old data into new slots */ 20717ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 208416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 209416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 210416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 211416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 21217ab2063SBarry Smith len*sizeof(int)); 213416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 214416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 21517ab2063SBarry Smith len*sizeof(Scalar)); 21617ab2063SBarry Smith /* free up old matrix storage */ 2170452661fSBarry Smith PetscFree(a->a); 2180452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 219416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 220416022c9SBarry Smith a->singlemalloc = 1; 22117ab2063SBarry Smith 22217ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 223416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 224416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 225416022c9SBarry Smith a->maxnz += CHUNKSIZE; 226b810aeb4SBarry Smith a->reallocs++; 22717ab2063SBarry Smith } 228416022c9SBarry Smith N = nrow++ - 1; a->nz++; 229416022c9SBarry Smith /* shift up all the later entries in this row */ 230416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 23117ab2063SBarry Smith rp[ii+1] = rp[ii]; 23217ab2063SBarry Smith ap[ii+1] = ap[ii]; 23317ab2063SBarry Smith } 23417ab2063SBarry Smith rp[i] = col; 23517ab2063SBarry Smith ap[i] = value; 23617ab2063SBarry Smith noinsert:; 237416022c9SBarry Smith low = i + 1; 23817ab2063SBarry Smith } 23917ab2063SBarry Smith ailen[row] = nrow; 24017ab2063SBarry Smith } 2413a40ed3dSBarry Smith PetscFunctionReturn(0); 24217ab2063SBarry Smith } 24317ab2063SBarry Smith 2445615d1e5SSatish Balay #undef __FUNC__ 2455615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ" 2468f6be9afSLois Curfman McInnes int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 2477eb43aa7SLois Curfman McInnes { 2487eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 249b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 2507eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 2517eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 2527eb43aa7SLois Curfman McInnes 2533a40ed3dSBarry Smith PetscFunctionBegin; 2547eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 2557eb43aa7SLois Curfman McInnes row = im[k]; 256a8c6a408SBarry Smith if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 257a8c6a408SBarry Smith if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 2587eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 2597eb43aa7SLois Curfman McInnes nrow = ailen[row]; 2607eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 261a8c6a408SBarry Smith if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 262a8c6a408SBarry Smith if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 2637eb43aa7SLois Curfman McInnes col = in[l] - shift; 2647eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 2657eb43aa7SLois Curfman McInnes while (high-low > 5) { 2667eb43aa7SLois Curfman McInnes t = (low+high)/2; 2677eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 2687eb43aa7SLois Curfman McInnes else low = t; 2697eb43aa7SLois Curfman McInnes } 2707eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 2717eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2727eb43aa7SLois Curfman McInnes if (rp[i] == col) { 273b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2747eb43aa7SLois Curfman McInnes goto finished; 2757eb43aa7SLois Curfman McInnes } 2767eb43aa7SLois Curfman McInnes } 277b49de8d1SLois Curfman McInnes *v++ = zero; 2787eb43aa7SLois Curfman McInnes finished:; 2797eb43aa7SLois Curfman McInnes } 2807eb43aa7SLois Curfman McInnes } 2813a40ed3dSBarry Smith PetscFunctionReturn(0); 2827eb43aa7SLois Curfman McInnes } 2837eb43aa7SLois Curfman McInnes 28417ab2063SBarry Smith 2855615d1e5SSatish Balay #undef __FUNC__ 2865615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary" 287480ef9eaSBarry Smith int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 28817ab2063SBarry Smith { 289416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 290416022c9SBarry Smith int i, fd, *col_lens, ierr; 29117ab2063SBarry Smith 2923a40ed3dSBarry Smith PetscFunctionBegin; 29390ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2940452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 295416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 296416022c9SBarry Smith col_lens[1] = a->m; 297416022c9SBarry Smith col_lens[2] = a->n; 298416022c9SBarry Smith col_lens[3] = a->nz; 299416022c9SBarry Smith 300416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 301416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 302416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 30317ab2063SBarry Smith } 3040752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1); CHKERRQ(ierr); 3050452661fSBarry Smith PetscFree(col_lens); 306416022c9SBarry Smith 307416022c9SBarry Smith /* store column indices (zero start index) */ 308416022c9SBarry Smith if (a->indexshift) { 309416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 31017ab2063SBarry Smith } 3110752156aSBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,0); CHKERRQ(ierr); 312416022c9SBarry Smith if (a->indexshift) { 313416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 31417ab2063SBarry Smith } 315416022c9SBarry Smith 316416022c9SBarry Smith /* store nonzero values */ 3170752156aSBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,0); CHKERRQ(ierr); 3183a40ed3dSBarry Smith PetscFunctionReturn(0); 31917ab2063SBarry Smith } 320416022c9SBarry Smith 3215615d1e5SSatish Balay #undef __FUNC__ 3225615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII" 323480ef9eaSBarry Smith int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 324416022c9SBarry Smith { 325416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 326496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 32717ab2063SBarry Smith FILE *fd; 32817ab2063SBarry Smith char *outputname; 32917ab2063SBarry Smith 3303a40ed3dSBarry Smith PetscFunctionBegin; 33190ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 332fb4b0f7fSBarry Smith ierr = ViewerGetOutputname(viewer,&outputname); CHKERRQ(ierr); 33390ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 334a93ec695SBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 3353a40ed3dSBarry Smith PetscFunctionReturn(0); 3363a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 337496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 338496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 339496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 34095e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 34195e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 3423a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 343d00d2cf4SBarry Smith int nofinalvalue = 0; 344d00d2cf4SBarry Smith if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) { 345d00d2cf4SBarry Smith nofinalvalue = 1; 346d00d2cf4SBarry Smith } 347416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3484e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 349d00d2cf4SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue); 35017ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 35117ab2063SBarry Smith 35217ab2063SBarry Smith for (i=0; i<m; i++) { 353416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 3543a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 355e20fef11SSatish 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])); 35617ab2063SBarry Smith #else 3577a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 35817ab2063SBarry Smith #endif 35917ab2063SBarry Smith } 36017ab2063SBarry Smith } 361d00d2cf4SBarry Smith if (nofinalvalue) { 362d00d2cf4SBarry Smith fprintf(fd,"%d %d %18.16e\n", m, a->n, 0.0); 363d00d2cf4SBarry Smith } 364e24b481bSBarry Smith if (outputname) fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 365e24b481bSBarry Smith else fprintf(fd,"];\n M = spconvert(zzz);\n"); 3663a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 36744cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 36844cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 36944cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 3703a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 371e20fef11SSatish Balay if (PetscImaginary(a->a[j]) > 0.0 && PetscReal(a->a[j]) != 0.0) 372e20fef11SSatish Balay fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j])); 373e20fef11SSatish Balay else if (PetscImaginary(a->a[j]) < 0.0 && PetscReal(a->a[j]) != 0.0) 374e20fef11SSatish Balay fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j])); 375e20fef11SSatish Balay else if (PetscReal(a->a[j]) != 0.0) 376e20fef11SSatish Balay fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j])); 37744cd7ae7SLois Curfman McInnes #else 37844cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 37944cd7ae7SLois Curfman McInnes #endif 38044cd7ae7SLois Curfman McInnes } 38144cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 38244cd7ae7SLois Curfman McInnes } 38302594712SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_SYMMODU) { 384496be53dSLois Curfman McInnes int nzd=0, fshift=1, *sptr; 3852e44a96cSLois Curfman McInnes sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr); 386496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 387496be53dSLois Curfman McInnes sptr[i] = nzd+1; 388496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 389496be53dSLois Curfman McInnes if (a->j[j] >= i) { 3903a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 391e20fef11SSatish Balay if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) nzd++; 392496be53dSLois Curfman McInnes #else 393496be53dSLois Curfman McInnes if (a->a[j] != 0.0) nzd++; 394496be53dSLois Curfman McInnes #endif 395496be53dSLois Curfman McInnes } 396496be53dSLois Curfman McInnes } 397496be53dSLois Curfman McInnes } 3982e44a96cSLois Curfman McInnes sptr[m] = nzd+1; 399496be53dSLois Curfman McInnes fprintf(fd," %d %d\n\n",m,nzd); 4002e44a96cSLois Curfman McInnes for ( i=0; i<m+1; i+=6 ) { 4012e44a96cSLois 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]); 4022e44a96cSLois 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]); 4032e44a96cSLois 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]); 4042e44a96cSLois Curfman McInnes else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]); 4052e44a96cSLois Curfman McInnes else if (i<m) fprintf(fd," %d %d\n",sptr[i],sptr[i+1]); 4067272d637SLois Curfman McInnes else fprintf(fd," %d\n",sptr[i]); 407496be53dSLois Curfman McInnes } 408496be53dSLois Curfman McInnes fprintf(fd,"\n"); 409496be53dSLois Curfman McInnes PetscFree(sptr); 410496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 411496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 412496be53dSLois Curfman McInnes if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift); 413496be53dSLois Curfman McInnes } 414496be53dSLois Curfman McInnes fprintf(fd,"\n"); 415496be53dSLois Curfman McInnes } 416496be53dSLois Curfman McInnes fprintf(fd,"\n"); 417496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 418496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 419496be53dSLois Curfman McInnes if (a->j[j] >= i) { 4203a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 421e20fef11SSatish Balay if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) 422e20fef11SSatish Balay fprintf(fd," %18.16e %18.16e ",PetscReal(a->a[j]),PetscImaginary(a->a[j])); 423496be53dSLois Curfman McInnes #else 424496be53dSLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]); 425496be53dSLois Curfman McInnes #endif 426496be53dSLois Curfman McInnes } 427496be53dSLois Curfman McInnes } 428496be53dSLois Curfman McInnes fprintf(fd,"\n"); 429496be53dSLois Curfman McInnes } 43002594712SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_DENSE) { 43102594712SBarry Smith int cnt = 0,jcnt; 43202594712SBarry Smith Scalar value; 43302594712SBarry Smith 43402594712SBarry Smith for ( i=0; i<m; i++ ) { 43502594712SBarry Smith jcnt = 0; 43602594712SBarry Smith for ( j=0; j<a->n; j++ ) { 437e24b481bSBarry Smith if ( jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) { 43802594712SBarry Smith value = a->a[cnt++]; 439e24b481bSBarry Smith jcnt++; 44002594712SBarry Smith } else { 44102594712SBarry Smith value = 0.0; 44202594712SBarry Smith } 44302594712SBarry Smith #if defined(USE_PETSC_COMPLEX) 44402594712SBarry Smith fprintf(fd," %7.5e+%7.5e i ",PetscReal(value),PetscImaginary(value)); 44502594712SBarry Smith #else 44602594712SBarry Smith fprintf(fd," %7.5e ",value); 44702594712SBarry Smith #endif 44802594712SBarry Smith } 44902594712SBarry Smith fprintf(fd,"\n"); 45002594712SBarry 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; 48077ed5343SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,color,rank; 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; 48477ed5343SBarry Smith MPI_Comm comm; 485cddf8d76SBarry Smith 4863a40ed3dSBarry Smith PetscFunctionBegin; 48777ed5343SBarry Smith /* 48877ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 48977ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 49077ed5343SBarry Smith rest should return immediately. 49177ed5343SBarry Smith */ 49277ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 49377ed5343SBarry Smith MPI_Comm_rank(comm,&rank); 49477ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 49577ed5343SBarry Smith 496480ef9eaSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 4970513a670SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 49819bcc07fSBarry Smith 499480ef9eaSBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr); 500416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 5010513a670SBarry Smith 5020513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 5030513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 504cddf8d76SBarry Smith color = DRAW_BLUE; 505416022c9SBarry Smith for ( i=0; i<m; i++ ) { 506cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 507416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 508cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5093a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 510e20fef11SSatish Balay if (PetscReal(a->a[j]) >= 0.) continue; 511cddf8d76SBarry Smith #else 512cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 513cddf8d76SBarry Smith #endif 514480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 515cddf8d76SBarry Smith } 516cddf8d76SBarry Smith } 517cddf8d76SBarry Smith color = DRAW_CYAN; 518cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 519cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 520cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 521cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 522cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 523480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 524cddf8d76SBarry Smith } 525cddf8d76SBarry Smith } 526cddf8d76SBarry Smith color = DRAW_RED; 527cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 528cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 529cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 530cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5313a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 532e20fef11SSatish Balay if (PetscReal(a->a[j]) <= 0.) continue; 533cddf8d76SBarry Smith #else 534cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 535cddf8d76SBarry Smith #endif 536480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 537416022c9SBarry Smith } 538416022c9SBarry Smith } 5390513a670SBarry Smith } else { 5400513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 5410513a670SBarry Smith /* first determine max of all nonzero values */ 5420513a670SBarry Smith int nz = a->nz,count; 5430513a670SBarry Smith Draw popup; 544480ef9eaSBarry Smith double scale; 5450513a670SBarry Smith 5460513a670SBarry Smith for ( i=0; i<nz; i++ ) { 5470513a670SBarry Smith if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 5480513a670SBarry Smith } 549480ef9eaSBarry Smith scale = (245.0 - DRAW_BASIC_COLORS)/maxv; 550522c5e43SBarry Smith ierr = DrawGetPopup(draw,&popup); CHKERRQ(ierr); 5510513a670SBarry Smith ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr); 5520513a670SBarry Smith count = 0; 5530513a670SBarry Smith for ( i=0; i<m; i++ ) { 5540513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 5550513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 5560513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5576d6b6c30SSatish Balay color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count])); 558480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 5590513a670SBarry Smith count++; 5600513a670SBarry Smith } 5610513a670SBarry Smith } 5620513a670SBarry Smith } 563480ef9eaSBarry Smith PetscFunctionReturn(0); 564480ef9eaSBarry Smith } 565cddf8d76SBarry Smith 566480ef9eaSBarry Smith #undef __FUNC__ 567480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw" 568480ef9eaSBarry Smith int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 569480ef9eaSBarry Smith { 570480ef9eaSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 571480ef9eaSBarry Smith int ierr; 572480ef9eaSBarry Smith Draw draw; 573480ef9eaSBarry Smith double xr,yr,xl,yl,h,w; 574480ef9eaSBarry Smith 575480ef9eaSBarry Smith PetscTruth isnull; 576480ef9eaSBarry Smith 577480ef9eaSBarry Smith PetscFunctionBegin; 57877ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw); CHKERRQ(ierr); 579480ef9eaSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); 580480ef9eaSBarry Smith if (isnull) PetscFunctionReturn(0); 581480ef9eaSBarry Smith 582480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 583480ef9eaSBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 584480ef9eaSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 585cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 586480ef9eaSBarry Smith ierr = DrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A); CHKERRQ(ierr); 587480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5883a40ed3dSBarry Smith PetscFunctionReturn(0); 589416022c9SBarry Smith } 590416022c9SBarry Smith 5915615d1e5SSatish Balay #undef __FUNC__ 592d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ" 593e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer) 594416022c9SBarry Smith { 595416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 596bcd2baecSBarry Smith ViewerType vtype; 597bcd2baecSBarry Smith int ierr; 598416022c9SBarry Smith 5993a40ed3dSBarry Smith PetscFunctionBegin; 600bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 601*3f1db9ecSBarry Smith if (PetscTypeCompare(vtype,MATLAB_VIEWER)) { 6023a40ed3dSBarry Smith ierr = ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);CHKERRQ(ierr); 603*3f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){ 6043a40ed3dSBarry Smith ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr); 605*3f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 6063a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr); 607*3f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 6083a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr); 6095cd90555SBarry Smith } else { 6105cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 61117ab2063SBarry Smith } 6123a40ed3dSBarry Smith PetscFunctionReturn(0); 61317ab2063SBarry Smith } 61419bcc07fSBarry Smith 615c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 6165615d1e5SSatish Balay #undef __FUNC__ 6175615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ" 6188f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 61917ab2063SBarry Smith { 620416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 62141c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 62243ee02c3SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0; 623416022c9SBarry Smith Scalar *aa = a->a, *ap; 62417ab2063SBarry Smith 6253a40ed3dSBarry Smith PetscFunctionBegin; 6263a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 62717ab2063SBarry Smith 62843ee02c3SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 62917ab2063SBarry Smith for ( i=1; i<m; i++ ) { 630416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 63117ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 63294a9d846SBarry Smith rmax = PetscMax(rmax,ailen[i]); 63317ab2063SBarry Smith if (fshift) { 634416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 63517ab2063SBarry Smith N = ailen[i]; 63617ab2063SBarry Smith for ( j=0; j<N; j++ ) { 63717ab2063SBarry Smith ip[j-fshift] = ip[j]; 63817ab2063SBarry Smith ap[j-fshift] = ap[j]; 63917ab2063SBarry Smith } 64017ab2063SBarry Smith } 64117ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 64217ab2063SBarry Smith } 64317ab2063SBarry Smith if (m) { 64417ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 64517ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 64617ab2063SBarry Smith } 64717ab2063SBarry Smith /* reset ilen and imax for each row */ 64817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 64917ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 65017ab2063SBarry Smith } 651416022c9SBarry Smith a->nz = ai[m] + shift; 65217ab2063SBarry Smith 65317ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 654416022c9SBarry Smith if (fshift && a->diag) { 6550452661fSBarry Smith PetscFree(a->diag); 656416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 657416022c9SBarry Smith a->diag = 0; 65817ab2063SBarry Smith } 6594e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 6604e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 6614e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 662b810aeb4SBarry Smith a->reallocs); 66394a9d846SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax); 664dd5f02e7SSatish Balay a->reallocs = 0; 6654e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 6664e220ebcSLois Curfman McInnes 66776dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 66841c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 6693a40ed3dSBarry Smith PetscFunctionReturn(0); 67017ab2063SBarry Smith } 67117ab2063SBarry Smith 6725615d1e5SSatish Balay #undef __FUNC__ 6735615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ" 6748f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A) 67517ab2063SBarry Smith { 676416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6773a40ed3dSBarry Smith 6783a40ed3dSBarry Smith PetscFunctionBegin; 679cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 6803a40ed3dSBarry Smith PetscFunctionReturn(0); 68117ab2063SBarry Smith } 682416022c9SBarry Smith 6835615d1e5SSatish Balay #undef __FUNC__ 6845615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ" 685e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A) 68617ab2063SBarry Smith { 687416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 68882bf6240SBarry Smith int ierr; 689d5d45c9bSBarry Smith 6903a40ed3dSBarry Smith PetscFunctionBegin; 69194d884c6SBarry Smith if (--A->refct > 0) PetscFunctionReturn(0); 69294d884c6SBarry Smith 69394d884c6SBarry Smith if (A->mapping) { 69494d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 69594d884c6SBarry Smith } 69694d884c6SBarry Smith if (A->bmapping) { 69794d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr); 69894d884c6SBarry Smith } 69961b13de0SBarry Smith if (A->rmap) { 70061b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 70161b13de0SBarry Smith } 70261b13de0SBarry Smith if (A->cmap) { 70361b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 70461b13de0SBarry Smith } 7053a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 706e1311b90SBarry Smith PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 70717ab2063SBarry Smith #endif 7080452661fSBarry Smith PetscFree(a->a); 7090452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 7100452661fSBarry Smith if (a->diag) PetscFree(a->diag); 7110452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 7120452661fSBarry Smith if (a->imax) PetscFree(a->imax); 7130452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 71476dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 71582bf6240SBarry Smith if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 716be6bf707SBarry Smith if (a->saved_values) PetscFree(a->saved_values); 7170452661fSBarry Smith PetscFree(a); 718eed86810SBarry Smith 719f2655603SLois Curfman McInnes PLogObjectDestroy(A); 720f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 7213a40ed3dSBarry Smith PetscFunctionReturn(0); 72217ab2063SBarry Smith } 72317ab2063SBarry Smith 7245615d1e5SSatish Balay #undef __FUNC__ 7255615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ" 7268f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A) 72717ab2063SBarry Smith { 7283a40ed3dSBarry Smith PetscFunctionBegin; 7293a40ed3dSBarry Smith PetscFunctionReturn(0); 73017ab2063SBarry Smith } 73117ab2063SBarry Smith 7325615d1e5SSatish Balay #undef __FUNC__ 7335615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ" 7348f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op) 73517ab2063SBarry Smith { 736416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7373a40ed3dSBarry Smith 7383a40ed3dSBarry Smith PetscFunctionBegin; 7396d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 7406d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 7416d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 742219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 7436d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 744c2653b3dSLois Curfman McInnes else if (op == MAT_NEW_NONZERO_LOCATION_ERROR) a->nonew = -1; 74596854ed6SLois Curfman McInnes else if (op == MAT_NEW_NONZERO_ALLOCATION_ERROR) a->nonew = -2; 7466d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 7476d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 748219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 7496d4a8577SBarry Smith op == MAT_SYMMETRIC || 7506d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 75190f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 752b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES|| 753b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) 754981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n"); 7553a40ed3dSBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) { 7563a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 7573a40ed3dSBarry Smith } else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 7586d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 7596d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 7606d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 7616d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 7623a40ed3dSBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 7633a40ed3dSBarry Smith PetscFunctionReturn(0); 76417ab2063SBarry Smith } 76517ab2063SBarry Smith 7665615d1e5SSatish Balay #undef __FUNC__ 7675615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ" 7688f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 76917ab2063SBarry Smith { 770416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7713a40ed3dSBarry Smith int i,j, n,shift = a->indexshift,ierr; 77217ab2063SBarry Smith Scalar *x, zero = 0.0; 77317ab2063SBarry Smith 7743a40ed3dSBarry Smith PetscFunctionBegin; 7753a40ed3dSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 776e1311b90SBarry Smith ierr = VecGetArray(v,&x);;CHKERRQ(ierr); 777e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr); 778a8c6a408SBarry Smith if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector"); 779416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 780416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 781416022c9SBarry Smith if (a->j[j]+shift == i) { 782416022c9SBarry Smith x[i] = a->a[j]; 78317ab2063SBarry Smith break; 78417ab2063SBarry Smith } 78517ab2063SBarry Smith } 78617ab2063SBarry Smith } 787e1311b90SBarry Smith ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr); 7883a40ed3dSBarry Smith PetscFunctionReturn(0); 78917ab2063SBarry Smith } 79017ab2063SBarry Smith 79117ab2063SBarry Smith /* -------------------------------------------------------*/ 79217ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 79317ab2063SBarry Smith /* -------------------------------------------------------*/ 7945615d1e5SSatish Balay #undef __FUNC__ 7955615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ" 79644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 79717ab2063SBarry Smith { 798416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 79917ab2063SBarry Smith Scalar *x, *y, *v, alpha; 800e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx, shift = a->indexshift; 80117ab2063SBarry Smith 8023a40ed3dSBarry Smith PetscFunctionBegin; 803e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 804e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 805cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 80617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 80717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 808416022c9SBarry Smith idx = a->j + a->i[i] + shift; 809416022c9SBarry Smith v = a->a + a->i[i] + shift; 810416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 81117ab2063SBarry Smith alpha = x[i]; 81217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 81317ab2063SBarry Smith } 814416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 815e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 816e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8173a40ed3dSBarry Smith PetscFunctionReturn(0); 81817ab2063SBarry Smith } 819d5d45c9bSBarry Smith 8205615d1e5SSatish Balay #undef __FUNC__ 8215615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ" 82244cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 82317ab2063SBarry Smith { 824416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 82517ab2063SBarry Smith Scalar *x, *y, *v, alpha; 826e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx,shift = a->indexshift; 82717ab2063SBarry Smith 8283a40ed3dSBarry Smith PetscFunctionBegin; 8292e8a6d31SBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 830e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 831e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 83217ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 83317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 834416022c9SBarry Smith idx = a->j + a->i[i] + shift; 835416022c9SBarry Smith v = a->a + a->i[i] + shift; 836416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 83717ab2063SBarry Smith alpha = x[i]; 83817ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 83917ab2063SBarry Smith } 84090f02eecSBarry Smith PLogFlops(2*a->nz); 841e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 842e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8433a40ed3dSBarry Smith PetscFunctionReturn(0); 84417ab2063SBarry Smith } 84517ab2063SBarry Smith 8465615d1e5SSatish Balay #undef __FUNC__ 8475615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ" 84844cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 84917ab2063SBarry Smith { 850416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 85117ab2063SBarry Smith Scalar *x, *y, *v, sum; 852e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 8532e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTAIJ) 854e36a17ebSSatish Balay int n, i, jrow,j; 855e36a17ebSSatish Balay #endif 85617ab2063SBarry Smith 857fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT) 858fee21e36SBarry Smith #pragma disjoint(*x,*y,*v) 859fee21e36SBarry Smith #endif 860fee21e36SBarry Smith 8613a40ed3dSBarry Smith PetscFunctionBegin; 862e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 863e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 86417ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 865416022c9SBarry Smith idx = a->j; 866d64ed03dSBarry Smith v = a->a; 867416022c9SBarry Smith ii = a->i; 868acc4d009SSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTAIJ) 86929b5ca8aSSatish Balay fortranmultaij_(&m,x,ii,idx+shift,v+shift,y); 8708d195f9aSBarry Smith #else 871d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 872d64ed03dSBarry Smith idx += shift; 87317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8749ea0dfa2SSatish Balay jrow = ii[i]; 8759ea0dfa2SSatish Balay n = ii[i+1] - jrow; 87617ab2063SBarry Smith sum = 0.0; 8779ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8789ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8799ea0dfa2SSatish Balay } 88017ab2063SBarry Smith y[i] = sum; 88117ab2063SBarry Smith } 8828d195f9aSBarry Smith #endif 883416022c9SBarry Smith PLogFlops(2*a->nz - m); 884e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 885e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 8863a40ed3dSBarry Smith PetscFunctionReturn(0); 88717ab2063SBarry Smith } 88817ab2063SBarry Smith 8895615d1e5SSatish Balay #undef __FUNC__ 8905615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ" 89144cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 89217ab2063SBarry Smith { 893416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 89417ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 895e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 8962e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTADDAIJ) 897e36a17ebSSatish Balay int n,i,jrow,j; 898e36a17ebSSatish Balay #endif 8999ea0dfa2SSatish Balay 9003a40ed3dSBarry Smith PetscFunctionBegin; 901e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 902e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 9032e8a6d31SBarry Smith if (zz != yy) { 904e1311b90SBarry Smith ierr = VecGetArray(zz,&z); CHKERRQ(ierr); 9052e8a6d31SBarry Smith } else { 9062e8a6d31SBarry Smith z = y; 9072e8a6d31SBarry Smith } 90817ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 909cddf8d76SBarry Smith idx = a->j; 910d64ed03dSBarry Smith v = a->a; 911cddf8d76SBarry Smith ii = a->i; 9122e40c62eSSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTADDAIJ) 91329b5ca8aSSatish Balay fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z); 91402ab625aSSatish Balay #else 915d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 916d64ed03dSBarry Smith idx += shift; 91717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 9189ea0dfa2SSatish Balay jrow = ii[i]; 9199ea0dfa2SSatish Balay n = ii[i+1] - jrow; 92017ab2063SBarry Smith sum = y[i]; 9219ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 9229ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 9239ea0dfa2SSatish Balay } 92417ab2063SBarry Smith z[i] = sum; 92517ab2063SBarry Smith } 92602ab625aSSatish Balay #endif 927416022c9SBarry Smith PLogFlops(2*a->nz); 928e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 929e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 9302e8a6d31SBarry Smith if (zz != yy) { 931e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr); 9322e8a6d31SBarry Smith } 9333a40ed3dSBarry Smith PetscFunctionReturn(0); 93417ab2063SBarry Smith } 93517ab2063SBarry Smith 93617ab2063SBarry Smith /* 93717ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 93817ab2063SBarry Smith */ 93917ab2063SBarry Smith 9405615d1e5SSatish Balay #undef __FUNC__ 9415615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ" 942416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 94317ab2063SBarry Smith { 944416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 945416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 94617ab2063SBarry Smith 9473a40ed3dSBarry Smith PetscFunctionBegin; 9480452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 949416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 950416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 95135b0346bSBarry Smith diag[i] = a->i[i+1]; 952416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 953416022c9SBarry Smith if (a->j[j]+shift == i) { 95417ab2063SBarry Smith diag[i] = j - shift; 95517ab2063SBarry Smith break; 95617ab2063SBarry Smith } 95717ab2063SBarry Smith } 95817ab2063SBarry Smith } 959416022c9SBarry Smith a->diag = diag; 9603a40ed3dSBarry Smith PetscFunctionReturn(0); 96117ab2063SBarry Smith } 96217ab2063SBarry Smith 9635615d1e5SSatish Balay #undef __FUNC__ 9645615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ" 96544cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 96617ab2063SBarry Smith double fshift,int its,Vec xx) 96717ab2063SBarry Smith { 968416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 969416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 970d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 97117ab2063SBarry Smith 9723a40ed3dSBarry Smith PetscFunctionBegin; 973e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 974fb2e594dSBarry Smith if (xx != bb) { 975e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 976fb2e594dSBarry Smith } else { 977fb2e594dSBarry Smith b = x; 978fb2e594dSBarry Smith } 979fb2e594dSBarry Smith 980d00d2cf4SBarry Smith if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);} 981416022c9SBarry Smith diag = a->diag; 982416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 98317ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 98417ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 98517ab2063SBarry Smith bs = b + shift; 98617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 987416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 988416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 989488ecbafSBarry Smith PLogFlops(2*n-1); 990416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 991416022c9SBarry Smith v = a->a + diag[i] + (!shift); 99217ab2063SBarry Smith sum = b[i]*d/omega; 99317ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 99417ab2063SBarry Smith x[i] = sum; 99517ab2063SBarry Smith } 996cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 997fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 9983a40ed3dSBarry Smith PetscFunctionReturn(0); 99917ab2063SBarry Smith } 100017ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 1001a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done"); 10023a40ed3dSBarry Smith } else if (flag & SOR_EISENSTAT) { 100317ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 100417ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 100517ab2063SBarry Smith 100617ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 100717ab2063SBarry Smith 100817ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 100917ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 101017ab2063SBarry Smith is the relaxation factor. 101117ab2063SBarry Smith */ 10120452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 101317ab2063SBarry Smith scale = (2.0/omega) - 1.0; 101417ab2063SBarry Smith 101517ab2063SBarry Smith /* x = (E + U)^{-1} b */ 101617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1017416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1018416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1019488ecbafSBarry Smith PLogFlops(2*n-1); 1020416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1021416022c9SBarry Smith v = a->a + diag[i] + (!shift); 102217ab2063SBarry Smith sum = b[i]; 102317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 102417ab2063SBarry Smith x[i] = omega*(sum/d); 102517ab2063SBarry Smith } 102617ab2063SBarry Smith 102717ab2063SBarry Smith /* t = b - (2*E - D)x */ 1028416022c9SBarry Smith v = a->a; 102917ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 1030488ecbafSBarry Smith PLogFlops(2*m); 1031488ecbafSBarry Smith 103217ab2063SBarry Smith 103317ab2063SBarry Smith /* t = (E + L)^{-1}t */ 1034416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 1035416022c9SBarry Smith diag = a->diag; 103617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1037416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1038416022c9SBarry Smith n = diag[i] - a->i[i]; 1039488ecbafSBarry Smith PLogFlops(2*n-1); 1040416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1041416022c9SBarry Smith v = a->a + a->i[i] + shift; 104217ab2063SBarry Smith sum = t[i]; 104317ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 104417ab2063SBarry Smith t[i] = omega*(sum/d); 104517ab2063SBarry Smith } 104617ab2063SBarry Smith 104717ab2063SBarry Smith /* x = x + t */ 104817ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 10490452661fSBarry Smith PetscFree(t); 1050cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1051fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 10523a40ed3dSBarry Smith PetscFunctionReturn(0); 105317ab2063SBarry Smith } 105417ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 105517ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 105617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1057416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1058416022c9SBarry Smith n = diag[i] - a->i[i]; 1059488ecbafSBarry Smith PLogFlops(2*n-1); 1060416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1061416022c9SBarry Smith v = a->a + a->i[i] + shift; 106217ab2063SBarry Smith sum = b[i]; 106317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 106417ab2063SBarry Smith x[i] = omega*(sum/d); 106517ab2063SBarry Smith } 106617ab2063SBarry Smith xb = x; 10673a40ed3dSBarry Smith } else xb = b; 106817ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 106917ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 107017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1071416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 107217ab2063SBarry Smith } 1073488ecbafSBarry Smith PLogFlops(m); 107417ab2063SBarry Smith } 107517ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 107617ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1077416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1078416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1079488ecbafSBarry Smith PLogFlops(2*n-1); 1080416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1081416022c9SBarry Smith v = a->a + diag[i] + (!shift); 108217ab2063SBarry Smith sum = xb[i]; 108317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 108417ab2063SBarry Smith x[i] = omega*(sum/d); 108517ab2063SBarry Smith } 108617ab2063SBarry Smith } 108717ab2063SBarry Smith its--; 108817ab2063SBarry Smith } 108917ab2063SBarry Smith while (its--) { 109017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 109117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1092416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1093416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1094488ecbafSBarry Smith PLogFlops(2*n-1); 1095416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1096416022c9SBarry Smith v = a->a + a->i[i] + shift; 109717ab2063SBarry Smith sum = b[i]; 109817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 10997e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 110017ab2063SBarry Smith } 110117ab2063SBarry Smith } 110217ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 110317ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1104416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1105416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1106488ecbafSBarry Smith PLogFlops(2*n-1); 1107416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1108416022c9SBarry Smith v = a->a + a->i[i] + shift; 110917ab2063SBarry Smith sum = b[i]; 111017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 11117e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 111217ab2063SBarry Smith } 111317ab2063SBarry Smith } 111417ab2063SBarry Smith } 1115cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1116fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 11173a40ed3dSBarry Smith PetscFunctionReturn(0); 111817ab2063SBarry Smith } 111917ab2063SBarry Smith 11205615d1e5SSatish Balay #undef __FUNC__ 11215615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ" 11228f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 112317ab2063SBarry Smith { 1124416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11254e220ebcSLois Curfman McInnes 11263a40ed3dSBarry Smith PetscFunctionBegin; 11274e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 11284e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 11294e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 11304e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 11314e220ebcSLois Curfman McInnes info->block_size = 1.0; 11324e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 11334e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 11344e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 11354e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 11364e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 11374e220ebcSLois Curfman McInnes info->memory = A->mem; 11384e220ebcSLois Curfman McInnes if (A->factor) { 11394e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 11404e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 11414e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 11424e220ebcSLois Curfman McInnes } else { 11434e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 11444e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 11454e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 11464e220ebcSLois Curfman McInnes } 11473a40ed3dSBarry Smith PetscFunctionReturn(0); 114817ab2063SBarry Smith } 114917ab2063SBarry Smith 115017ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 115117ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 115217ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 115317ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 115417ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 115517ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 115617ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 115717ab2063SBarry Smith 11585615d1e5SSatish Balay #undef __FUNC__ 11595615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ" 11608f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 116117ab2063SBarry Smith { 1162416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1163416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 116417ab2063SBarry Smith 11653a40ed3dSBarry Smith PetscFunctionBegin; 116677c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 116717ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 116817ab2063SBarry Smith if (diag) { 116917ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1170a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1171416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 1172416022c9SBarry Smith a->ilen[rows[i]] = 1; 1173416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 1174416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 11753a40ed3dSBarry Smith } else { 1176d64ed03dSBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr); 117717ab2063SBarry Smith } 117817ab2063SBarry Smith } 11793a40ed3dSBarry Smith } else { 118017ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1181a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1182416022c9SBarry Smith a->ilen[rows[i]] = 0; 118317ab2063SBarry Smith } 118417ab2063SBarry Smith } 118517ab2063SBarry Smith ISRestoreIndices(is,&rows); 118643a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11873a40ed3dSBarry Smith PetscFunctionReturn(0); 118817ab2063SBarry Smith } 118917ab2063SBarry Smith 11905615d1e5SSatish Balay #undef __FUNC__ 11915615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ" 11928f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 119317ab2063SBarry Smith { 1194416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11953a40ed3dSBarry Smith 11963a40ed3dSBarry Smith PetscFunctionBegin; 11970752156aSBarry Smith if (m) *m = a->m; 11980752156aSBarry Smith if (n) *n = a->n; 11993a40ed3dSBarry Smith PetscFunctionReturn(0); 120017ab2063SBarry Smith } 120117ab2063SBarry Smith 12025615d1e5SSatish Balay #undef __FUNC__ 12035615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ" 12048f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 120517ab2063SBarry Smith { 1206416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 12073a40ed3dSBarry Smith 12083a40ed3dSBarry Smith PetscFunctionBegin; 1209416022c9SBarry Smith *m = 0; *n = a->m; 12103a40ed3dSBarry Smith PetscFunctionReturn(0); 121117ab2063SBarry Smith } 1212026e39d0SSatish Balay 12135615d1e5SSatish Balay #undef __FUNC__ 12145615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ" 12154e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 121617ab2063SBarry Smith { 1217416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1218c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 121917ab2063SBarry Smith 12203a40ed3dSBarry Smith PetscFunctionBegin; 1221a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 122217ab2063SBarry Smith 1223416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 1224416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 122517ab2063SBarry Smith if (idx) { 1226416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 12274e093b46SBarry Smith if (*nz && shift) { 12280452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 122917ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 12304e093b46SBarry Smith } else if (*nz) { 12314e093b46SBarry Smith *idx = itmp; 123217ab2063SBarry Smith } 123317ab2063SBarry Smith else *idx = 0; 123417ab2063SBarry Smith } 12353a40ed3dSBarry Smith PetscFunctionReturn(0); 123617ab2063SBarry Smith } 123717ab2063SBarry Smith 12385615d1e5SSatish Balay #undef __FUNC__ 12395615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ" 12404e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 124117ab2063SBarry Smith { 12424e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 12433a40ed3dSBarry Smith 12443a40ed3dSBarry Smith PetscFunctionBegin; 12454e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 12463a40ed3dSBarry Smith PetscFunctionReturn(0); 124717ab2063SBarry Smith } 124817ab2063SBarry Smith 12495615d1e5SSatish Balay #undef __FUNC__ 12505615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ" 12518f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 125217ab2063SBarry Smith { 1253416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1254416022c9SBarry Smith Scalar *v = a->a; 125517ab2063SBarry Smith double sum = 0.0; 1256416022c9SBarry Smith int i, j,shift = a->indexshift; 125717ab2063SBarry Smith 12583a40ed3dSBarry Smith PetscFunctionBegin; 125917ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1260416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 12613a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 1262e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 126317ab2063SBarry Smith #else 126417ab2063SBarry Smith sum += (*v)*(*v); v++; 126517ab2063SBarry Smith #endif 126617ab2063SBarry Smith } 126717ab2063SBarry Smith *norm = sqrt(sum); 12683a40ed3dSBarry Smith } else if (type == NORM_1) { 126917ab2063SBarry Smith double *tmp; 1270416022c9SBarry Smith int *jj = a->j; 127166963ce1SSatish Balay tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp); 1272cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 127317ab2063SBarry Smith *norm = 0.0; 1274416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1275a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 127617ab2063SBarry Smith } 1277416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 127817ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 127917ab2063SBarry Smith } 12800452661fSBarry Smith PetscFree(tmp); 12813a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 128217ab2063SBarry Smith *norm = 0.0; 1283416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1284416022c9SBarry Smith v = a->a + a->i[j] + shift; 128517ab2063SBarry Smith sum = 0.0; 1286416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1287cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 128817ab2063SBarry Smith } 128917ab2063SBarry Smith if (sum > *norm) *norm = sum; 129017ab2063SBarry Smith } 12913a40ed3dSBarry Smith } else { 1292a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 129317ab2063SBarry Smith } 12943a40ed3dSBarry Smith PetscFunctionReturn(0); 129517ab2063SBarry Smith } 129617ab2063SBarry Smith 12975615d1e5SSatish Balay #undef __FUNC__ 12985615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ" 12998f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B) 130017ab2063SBarry Smith { 1301416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1302416022c9SBarry Smith Mat C; 1303416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1304416022c9SBarry Smith int shift = a->indexshift; 1305d5d45c9bSBarry Smith Scalar *array = a->a; 130617ab2063SBarry Smith 13073a40ed3dSBarry Smith PetscFunctionBegin; 1308a8c6a408SBarry Smith if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place"); 13090452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1310cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 131117ab2063SBarry Smith if (shift) { 131217ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 131317ab2063SBarry Smith } 131417ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1315416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 13160452661fSBarry Smith PetscFree(col); 131717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 131817ab2063SBarry Smith len = ai[i+1]-ai[i]; 1319416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 132017ab2063SBarry Smith array += len; aj += len; 132117ab2063SBarry Smith } 132217ab2063SBarry Smith if (shift) { 132317ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 132417ab2063SBarry Smith } 132517ab2063SBarry Smith 13266d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13276d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 132817ab2063SBarry Smith 13293638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1330416022c9SBarry Smith *B = C; 133117ab2063SBarry Smith } else { 1332f830108cSBarry Smith PetscOps *Abops; 13330a6ffc59SBarry Smith MatOps Aops; 1334f830108cSBarry Smith 1335416022c9SBarry Smith /* This isn't really an in-place transpose */ 13360452661fSBarry Smith PetscFree(a->a); 13370452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 13380452661fSBarry Smith if (a->diag) PetscFree(a->diag); 13390452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 13400452661fSBarry Smith if (a->imax) PetscFree(a->imax); 13410452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 13421073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 13430452661fSBarry Smith PetscFree(a); 1344f830108cSBarry Smith 1345488ecbafSBarry Smith 1346488ecbafSBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 1347488ecbafSBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 1348488ecbafSBarry Smith 1349f830108cSBarry Smith /* 1350f830108cSBarry Smith This is horrible, horrible code. We need to keep the 13518f0f457cSSatish Balay the bops and ops Structures, copy everything from C 13528f0f457cSSatish Balay including the function pointers.. 1353f830108cSBarry Smith */ 13548f0f457cSSatish Balay PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps)); 13558f0f457cSSatish Balay PetscMemcpy(A->bops,C->bops,sizeof(PetscOps)); 1356f830108cSBarry Smith Abops = A->bops; 1357f830108cSBarry Smith Aops = A->ops; 1358f09e8eb9SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 1359f830108cSBarry Smith A->bops = Abops; 1360f830108cSBarry Smith A->ops = Aops; 136127a8da17SBarry Smith A->qlist = 0; 1362f830108cSBarry Smith 13630452661fSBarry Smith PetscHeaderDestroy(C); 136417ab2063SBarry Smith } 13653a40ed3dSBarry Smith PetscFunctionReturn(0); 136617ab2063SBarry Smith } 136717ab2063SBarry Smith 13685615d1e5SSatish Balay #undef __FUNC__ 13695615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ" 13708f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 137117ab2063SBarry Smith { 1372416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 137317ab2063SBarry Smith Scalar *l,*r,x,*v; 1374e1311b90SBarry Smith int ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 137517ab2063SBarry Smith 13763a40ed3dSBarry Smith PetscFunctionBegin; 137717ab2063SBarry Smith if (ll) { 13783ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 13793ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 1380e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr); 1381a8c6a408SBarry Smith if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length"); 1382e1311b90SBarry Smith ierr = VecGetArray(ll,&l); CHKERRQ(ierr); 1383416022c9SBarry Smith v = a->a; 138417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 138517ab2063SBarry Smith x = l[i]; 1386416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 138717ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 138817ab2063SBarry Smith } 1389e1311b90SBarry Smith ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr); 139044cd7ae7SLois Curfman McInnes PLogFlops(nz); 139117ab2063SBarry Smith } 139217ab2063SBarry Smith if (rr) { 1393e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr); 1394a8c6a408SBarry Smith if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length"); 1395e1311b90SBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1396416022c9SBarry Smith v = a->a; jj = a->j; 139717ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 139817ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 139917ab2063SBarry Smith } 1400e1311b90SBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 140144cd7ae7SLois Curfman McInnes PLogFlops(nz); 140217ab2063SBarry Smith } 14033a40ed3dSBarry Smith PetscFunctionReturn(0); 140417ab2063SBarry Smith } 140517ab2063SBarry Smith 14065615d1e5SSatish Balay #undef __FUNC__ 14075615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ" 14086a6a5d1dSBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatGetSubMatrixCall scall,Mat *B) 140917ab2063SBarry Smith { 1410db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 1411d5db1b72SBarry Smith int *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 141299141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1413a2744918SBarry Smith register int sum,lensi; 141499141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 141599141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 141699141d43SSatish Balay Scalar *a_new,*mat_a; 1417416022c9SBarry Smith Mat C; 1418fee21e36SBarry Smith PetscTruth stride; 141917ab2063SBarry Smith 14203a40ed3dSBarry Smith PetscFunctionBegin; 1421d64ed03dSBarry Smith ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr); 1422a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted"); 1423d64ed03dSBarry Smith ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 1424a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted"); 142599141d43SSatish Balay 142617ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 142717ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 142817ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 142917ab2063SBarry Smith 1430fee21e36SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr); 1431fee21e36SBarry Smith ierr = ISStride(iscol,&stride); CHKERRQ(ierr); 1432fee21e36SBarry Smith if (stride && step == 1) { 143302834360SBarry Smith /* special case of contiguous rows */ 143457faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 143502834360SBarry Smith starts = lens + ncols; 143602834360SBarry Smith /* loop over new rows determining lens and starting points */ 143702834360SBarry Smith for (i=0; i<nrows; i++) { 1438a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1439a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 144002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1441d8ced48eSBarry Smith if (aj[k]+shift >= first) { 144202834360SBarry Smith starts[i] = k; 144302834360SBarry Smith break; 144402834360SBarry Smith } 144502834360SBarry Smith } 1446a2744918SBarry Smith sum = 0; 144702834360SBarry Smith while (k < kend) { 1448d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1449a2744918SBarry Smith sum++; 145002834360SBarry Smith } 1451a2744918SBarry Smith lens[i] = sum; 145202834360SBarry Smith } 145302834360SBarry Smith /* create submatrix */ 1454cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 145508480c60SBarry Smith int n_cols,n_rows; 145608480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1457a8c6a408SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); 1458d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 145908480c60SBarry Smith C = *B; 14603a40ed3dSBarry Smith } else { 146102834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 146208480c60SBarry Smith } 1463db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1464db02288aSLois Curfman McInnes 146502834360SBarry Smith /* loop over rows inserting into submatrix */ 1466db02288aSLois Curfman McInnes a_new = c->a; 1467db02288aSLois Curfman McInnes j_new = c->j; 1468db02288aSLois Curfman McInnes i_new = c->i; 1469db02288aSLois Curfman McInnes i_new[0] = -shift; 147002834360SBarry Smith for (i=0; i<nrows; i++) { 1471a2744918SBarry Smith ii = starts[i]; 1472a2744918SBarry Smith lensi = lens[i]; 1473a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1474a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 147502834360SBarry Smith } 1476a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1477a2744918SBarry Smith a_new += lensi; 1478a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1479a2744918SBarry Smith c->ilen[i] = lensi; 148002834360SBarry Smith } 14810452661fSBarry Smith PetscFree(lens); 14823a40ed3dSBarry Smith } else { 148302834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 14840452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 148502834360SBarry Smith ssmap = smap + shift; 148699141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1487cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 148817ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 148902834360SBarry Smith /* determine lens of each row */ 149002834360SBarry Smith for (i=0; i<nrows; i++) { 1491d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 149202834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 149302834360SBarry Smith lens[i] = 0; 149402834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1495d8ced48eSBarry Smith if (ssmap[aj[k]]) { 149602834360SBarry Smith lens[i]++; 149702834360SBarry Smith } 149802834360SBarry Smith } 149902834360SBarry Smith } 150017ab2063SBarry Smith /* Create and fill new matrix */ 1501a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 150299141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 1503a8c6a408SBarry Smith if (c->m != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size"); 150499141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1505a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros"); 150699141d43SSatish Balay } 150799141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 150808480c60SBarry Smith C = *B; 15093a40ed3dSBarry Smith } else { 151002834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 151108480c60SBarry Smith } 151299141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 151317ab2063SBarry Smith for (i=0; i<nrows; i++) { 151499141d43SSatish Balay row = irow[i]; 151599141d43SSatish Balay kstart = ai[row]+shift; 151699141d43SSatish Balay kend = kstart + a->ilen[row]; 151799141d43SSatish Balay mat_i = c->i[i]+shift; 151899141d43SSatish Balay mat_j = c->j + mat_i; 151999141d43SSatish Balay mat_a = c->a + mat_i; 152099141d43SSatish Balay mat_ilen = c->ilen + i; 152117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 152299141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 152399141d43SSatish Balay *mat_j++ = tcol - (!shift); 152499141d43SSatish Balay *mat_a++ = a->a[k]; 152599141d43SSatish Balay (*mat_ilen)++; 152699141d43SSatish Balay 152717ab2063SBarry Smith } 152817ab2063SBarry Smith } 152917ab2063SBarry Smith } 153002834360SBarry Smith /* Free work space */ 153102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 153299141d43SSatish Balay PetscFree(smap); PetscFree(lens); 153302834360SBarry Smith } 15346d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 15356d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 153617ab2063SBarry Smith 153717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1538416022c9SBarry Smith *B = C; 15393a40ed3dSBarry Smith PetscFunctionReturn(0); 154017ab2063SBarry Smith } 154117ab2063SBarry Smith 1542a871dcd8SBarry Smith /* 154363b91edcSBarry Smith note: This can only work for identity for row and col. It would 154463b91edcSBarry Smith be good to check this and otherwise generate an error. 1545a871dcd8SBarry Smith */ 15465615d1e5SSatish Balay #undef __FUNC__ 15475615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ" 15488f6be9afSLois Curfman McInnes int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1549a871dcd8SBarry Smith { 155063b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 155108480c60SBarry Smith int ierr; 155263b91edcSBarry Smith Mat outA; 155363b91edcSBarry Smith 15543a40ed3dSBarry Smith PetscFunctionBegin; 1555a8c6a408SBarry Smith if (fill != 0) SETERRQ(PETSC_ERR_SUP,0,"Only fill=0 supported"); 1556a871dcd8SBarry Smith 155763b91edcSBarry Smith outA = inA; 155863b91edcSBarry Smith inA->factor = FACTOR_LU; 155963b91edcSBarry Smith a->row = row; 156063b91edcSBarry Smith a->col = col; 156163b91edcSBarry Smith 1562f0ec6fceSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1563f0ec6fceSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 15641e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1565f0ec6fceSSatish Balay 156694a9d846SBarry Smith if (!a->solve_work) { /* this matrix may have been factored before */ 15670452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work); 156894a9d846SBarry Smith } 156963b91edcSBarry Smith 157008480c60SBarry Smith if (!a->diag) { 157108480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 157263b91edcSBarry Smith } 157363b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 15743a40ed3dSBarry Smith PetscFunctionReturn(0); 1575a871dcd8SBarry Smith } 1576a871dcd8SBarry Smith 157774cdf0dfSBarry Smith #include "pinclude/blaslapack.h" 15785615d1e5SSatish Balay #undef __FUNC__ 15795615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ" 15808f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1581f0b747eeSBarry Smith { 1582f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1583f0b747eeSBarry Smith int one = 1; 15843a40ed3dSBarry Smith 15853a40ed3dSBarry Smith PetscFunctionBegin; 1586f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1587f0b747eeSBarry Smith PLogFlops(a->nz); 15883a40ed3dSBarry Smith PetscFunctionReturn(0); 1589f0b747eeSBarry Smith } 1590f0b747eeSBarry Smith 15915615d1e5SSatish Balay #undef __FUNC__ 15925615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ" 15936a6a5d1dSBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B) 1594cddf8d76SBarry Smith { 1595cddf8d76SBarry Smith int ierr,i; 1596cddf8d76SBarry Smith 15973a40ed3dSBarry Smith PetscFunctionBegin; 1598cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 15990452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1600cddf8d76SBarry Smith } 1601cddf8d76SBarry Smith 1602cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 16036a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1604cddf8d76SBarry Smith } 16053a40ed3dSBarry Smith PetscFunctionReturn(0); 1606cddf8d76SBarry Smith } 1607cddf8d76SBarry Smith 16085615d1e5SSatish Balay #undef __FUNC__ 16095615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ" 16108f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 16115a838052SSatish Balay { 1612f830108cSBarry Smith PetscFunctionBegin; 16135a838052SSatish Balay *bs = 1; 16143a40ed3dSBarry Smith PetscFunctionReturn(0); 16155a838052SSatish Balay } 16165a838052SSatish Balay 16175615d1e5SSatish Balay #undef __FUNC__ 16185615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ" 16198f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 16204dcbc457SBarry Smith { 1621e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 162206763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 16238a047759SSatish Balay int start, end, *ai, *aj; 1624bbd702dbSSatish Balay BT table; 1625bbd702dbSSatish Balay 16263a40ed3dSBarry Smith PetscFunctionBegin; 16278a047759SSatish Balay shift = a->indexshift; 1628e4d965acSSatish Balay m = a->m; 1629e4d965acSSatish Balay ai = a->i; 16308a047759SSatish Balay aj = a->j+shift; 16318a047759SSatish Balay 1632a8c6a408SBarry Smith if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used"); 163306763907SSatish Balay 163406763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 1635bbd702dbSSatish Balay ierr = BTCreate(m,table); CHKERRQ(ierr); 163606763907SSatish Balay 1637e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1638b97fc60eSLois Curfman McInnes /* Initialize the two local arrays */ 1639e4d965acSSatish Balay isz = 0; 1640bbd702dbSSatish Balay BTMemzero(m,table); 1641e4d965acSSatish Balay 1642e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 16434dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 164477c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1645e4d965acSSatish Balay 1646dd097bc3SLois Curfman McInnes /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 1647e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 1648bbd702dbSSatish Balay if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];} 16494dcbc457SBarry Smith } 165006763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 165106763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1652e4d965acSSatish Balay 165304a348a9SBarry Smith k = 0; 165404a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap */ 165504a348a9SBarry Smith n = isz; 165606763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1657e4d965acSSatish Balay row = nidx[k]; 1658e4d965acSSatish Balay start = ai[row]; 1659e4d965acSSatish Balay end = ai[row+1]; 166004a348a9SBarry Smith for ( l = start; l<end ; l++){ 16618a047759SSatish Balay val = aj[l] + shift; 16622a8f2162SSatish Balay if (!BTLookupSet(table,val)) {nidx[isz++] = val;} 1663e4d965acSSatish Balay } 1664e4d965acSSatish Balay } 1665e4d965acSSatish Balay } 1666029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1667e4d965acSSatish Balay } 1668bbd702dbSSatish Balay BTDestroy(table); 166906763907SSatish Balay PetscFree(nidx); 16703a40ed3dSBarry Smith PetscFunctionReturn(0); 16714dcbc457SBarry Smith } 167217ab2063SBarry Smith 16730513a670SBarry Smith /* -------------------------------------------------------------- */ 16740513a670SBarry Smith #undef __FUNC__ 16750513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ" 16760513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 16770513a670SBarry Smith { 16780513a670SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 16790513a670SBarry Smith Scalar *vwork; 16800513a670SBarry Smith int i, ierr, nz, m = a->m, n = a->n, *cwork; 16810513a670SBarry Smith int *row,*col,*cnew,j,*lens; 168256cd22aeSBarry Smith IS icolp,irowp; 16830513a670SBarry Smith 16843a40ed3dSBarry Smith PetscFunctionBegin; 168556cd22aeSBarry Smith ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr); 168656cd22aeSBarry Smith ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr); 168756cd22aeSBarry Smith ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr); 168856cd22aeSBarry Smith ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr); 16890513a670SBarry Smith 16900513a670SBarry Smith /* determine lengths of permuted rows */ 16910513a670SBarry Smith lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens); 16920513a670SBarry Smith for (i=0; i<m; i++ ) { 16930513a670SBarry Smith lens[row[i]] = a->i[i+1] - a->i[i]; 16940513a670SBarry Smith } 16950513a670SBarry Smith ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr); 16960513a670SBarry Smith PetscFree(lens); 16970513a670SBarry Smith 16980513a670SBarry Smith cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew); 16990513a670SBarry Smith for (i=0; i<m; i++) { 17000513a670SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 17010513a670SBarry Smith for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];} 17020513a670SBarry Smith ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr); 17030513a670SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 17040513a670SBarry Smith } 17050513a670SBarry Smith PetscFree(cnew); 17060513a670SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 17070513a670SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 170856cd22aeSBarry Smith ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr); 170956cd22aeSBarry Smith ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr); 171056cd22aeSBarry Smith ierr = ISDestroy(irowp); CHKERRQ(ierr); 171156cd22aeSBarry Smith ierr = ISDestroy(icolp); CHKERRQ(ierr); 17123a40ed3dSBarry Smith PetscFunctionReturn(0); 17130513a670SBarry Smith } 17140513a670SBarry Smith 17155615d1e5SSatish Balay #undef __FUNC__ 17165615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ" 1717682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1718682d7d0cSBarry Smith { 1719682d7d0cSBarry Smith static int called = 0; 1720682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1721682d7d0cSBarry Smith 17223a40ed3dSBarry Smith PetscFunctionBegin; 17233a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 172476be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 172576be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n"); 172676be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n"); 172776be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_no_inode: Do not use inodes\n"); 172876be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n"); 1729682d7d0cSBarry Smith #if defined(HAVE_ESSL) 173076be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.\n"); 1731682d7d0cSBarry Smith #endif 17323a40ed3dSBarry Smith PetscFunctionReturn(0); 1733682d7d0cSBarry Smith } 17348f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1735a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1736a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1737a93ec695SBarry Smith 1738cb5b572fSBarry Smith #undef __FUNC__ 1739cb5b572fSBarry Smith #define __FUNC__ "MatCopy_SeqAIJ" 1740cb5b572fSBarry Smith int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str) 1741cb5b572fSBarry Smith { 1742be6bf707SBarry Smith int ierr; 1743cb5b572fSBarry Smith 1744cb5b572fSBarry Smith PetscFunctionBegin; 1745be6bf707SBarry Smith if (str == SAME_NONZERO_PATTERN && B->type == MATSEQAIJ) { 1746be6bf707SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1747be6bf707SBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data; 1748be6bf707SBarry Smith 1749be6bf707SBarry Smith if (a->nonew != 1) SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1750be6bf707SBarry Smith if (b->nonew != 1) SETERRQ(1,1,"Must call MatSetOption(B,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1751be6bf707SBarry Smith if (a->i[a->m]+a->indexshift != b->i[b->m]+a->indexshift) { 1752be6bf707SBarry Smith SETERRQ(1,1,"Number of nonzeros in two matrices are different"); 1753cb5b572fSBarry Smith } 1754be6bf707SBarry Smith PetscMemcpy(b->a,a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 1755cb5b572fSBarry Smith } else { 1756cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1757cb5b572fSBarry Smith } 1758cb5b572fSBarry Smith PetscFunctionReturn(0); 1759cb5b572fSBarry Smith } 1760cb5b572fSBarry Smith 1761682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 17620a6ffc59SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ, 1763cb5b572fSBarry Smith MatGetRow_SeqAIJ, 1764cb5b572fSBarry Smith MatRestoreRow_SeqAIJ, 1765cb5b572fSBarry Smith MatMult_SeqAIJ, 1766cb5b572fSBarry Smith MatMultAdd_SeqAIJ, 1767cb5b572fSBarry Smith MatMultTrans_SeqAIJ, 1768cb5b572fSBarry Smith MatMultTransAdd_SeqAIJ, 1769cb5b572fSBarry Smith MatSolve_SeqAIJ, 1770cb5b572fSBarry Smith MatSolveAdd_SeqAIJ, 1771cb5b572fSBarry Smith MatSolveTrans_SeqAIJ, 1772cb5b572fSBarry Smith MatSolveTransAdd_SeqAIJ, 1773cb5b572fSBarry Smith MatLUFactor_SeqAIJ, 1774cb5b572fSBarry Smith 0, 177517ab2063SBarry Smith MatRelax_SeqAIJ, 177617ab2063SBarry Smith MatTranspose_SeqAIJ, 1777cb5b572fSBarry Smith MatGetInfo_SeqAIJ, 1778cb5b572fSBarry Smith MatEqual_SeqAIJ, 1779cb5b572fSBarry Smith MatGetDiagonal_SeqAIJ, 1780cb5b572fSBarry Smith MatDiagonalScale_SeqAIJ, 1781cb5b572fSBarry Smith MatNorm_SeqAIJ, 1782cb5b572fSBarry Smith 0, 1783cb5b572fSBarry Smith MatAssemblyEnd_SeqAIJ, 178417ab2063SBarry Smith MatCompress_SeqAIJ, 1785cb5b572fSBarry Smith MatSetOption_SeqAIJ, 1786cb5b572fSBarry Smith MatZeroEntries_SeqAIJ, 1787cb5b572fSBarry Smith MatZeroRows_SeqAIJ, 1788cb5b572fSBarry Smith MatLUFactorSymbolic_SeqAIJ, 1789cb5b572fSBarry Smith MatLUFactorNumeric_SeqAIJ, 1790cb5b572fSBarry Smith 0, 1791cb5b572fSBarry Smith 0, 1792cb5b572fSBarry Smith MatGetSize_SeqAIJ, 1793cb5b572fSBarry Smith MatGetSize_SeqAIJ, 1794cb5b572fSBarry Smith MatGetOwnershipRange_SeqAIJ, 1795cb5b572fSBarry Smith MatILUFactorSymbolic_SeqAIJ, 1796cb5b572fSBarry Smith 0, 1797cb5b572fSBarry Smith 0, 1798cb5b572fSBarry Smith 0, 1799be6bf707SBarry Smith MatDuplicate_SeqAIJ, 1800cb5b572fSBarry Smith 0, 1801cb5b572fSBarry Smith 0, 1802cb5b572fSBarry Smith MatILUFactor_SeqAIJ, 1803cb5b572fSBarry Smith 0, 1804cb5b572fSBarry Smith 0, 1805cb5b572fSBarry Smith MatGetSubMatrices_SeqAIJ, 1806cb5b572fSBarry Smith MatIncreaseOverlap_SeqAIJ, 1807cb5b572fSBarry Smith MatGetValues_SeqAIJ, 1808cb5b572fSBarry Smith MatCopy_SeqAIJ, 1809f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1810cb5b572fSBarry Smith MatScale_SeqAIJ, 1811cb5b572fSBarry Smith 0, 1812cb5b572fSBarry Smith 0, 18136945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 18146945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 18153b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 18163b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 18173b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1818a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1819a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 18200513a670SBarry Smith MatColoringPatch_SeqAIJ, 18210513a670SBarry Smith 0, 1822cda55fadSBarry Smith MatPermute_SeqAIJ, 1823cda55fadSBarry Smith 0, 1824cda55fadSBarry Smith 0, 1825cda55fadSBarry Smith 0, 1826cda55fadSBarry Smith 0, 1827cda55fadSBarry Smith MatGetMaps_Petsc}; 182817ab2063SBarry Smith 182917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 183017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 183117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 183217ab2063SBarry Smith 1833fb2e594dSBarry Smith EXTERN_C_BEGIN 18345615d1e5SSatish Balay #undef __FUNC__ 1835bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ" 1836bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices) 1837bef8e0ddSBarry Smith { 1838bef8e0ddSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1839bef8e0ddSBarry Smith int i,nz,n; 1840bef8e0ddSBarry Smith 1841bef8e0ddSBarry Smith PetscFunctionBegin; 1842bef8e0ddSBarry Smith if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing"); 1843bef8e0ddSBarry Smith 1844bef8e0ddSBarry Smith nz = aij->maxnz; 1845bef8e0ddSBarry Smith n = aij->n; 1846bef8e0ddSBarry Smith for (i=0; i<nz; i++) { 1847bef8e0ddSBarry Smith aij->j[i] = indices[i]; 1848bef8e0ddSBarry Smith } 1849bef8e0ddSBarry Smith aij->nz = nz; 1850bef8e0ddSBarry Smith for ( i=0; i<n; i++ ) { 1851bef8e0ddSBarry Smith aij->ilen[i] = aij->imax[i]; 1852bef8e0ddSBarry Smith } 1853bef8e0ddSBarry Smith 1854bef8e0ddSBarry Smith PetscFunctionReturn(0); 1855bef8e0ddSBarry Smith } 1856fb2e594dSBarry Smith EXTERN_C_END 1857bef8e0ddSBarry Smith 1858bef8e0ddSBarry Smith #undef __FUNC__ 1859bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices" 1860bef8e0ddSBarry Smith /*@ 1861bef8e0ddSBarry Smith MatSeqAIJSetColumnIndices - Set the column indices for all the rows 1862bef8e0ddSBarry Smith in the matrix. 1863bef8e0ddSBarry Smith 1864bef8e0ddSBarry Smith Input Parameters: 1865bef8e0ddSBarry Smith + mat - the SeqAIJ matrix 1866bef8e0ddSBarry Smith - indices - the column indices 1867bef8e0ddSBarry Smith 1868bef8e0ddSBarry Smith Notes: 1869bef8e0ddSBarry Smith This can be called if you have precomputed the nonzero structure of the 1870bef8e0ddSBarry Smith matrix and want to provide it to the matrix object to improve the performance 1871bef8e0ddSBarry Smith of the MatSetValues() operation. 1872bef8e0ddSBarry Smith 1873bef8e0ddSBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 1874bef8e0ddSBarry Smith MatCreateSeqAIJ(). 1875bef8e0ddSBarry Smith 1876bef8e0ddSBarry Smith MUST be called before any calls to MatSetValues(); 1877bef8e0ddSBarry Smith 1878bef8e0ddSBarry Smith @*/ 1879bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices) 1880bef8e0ddSBarry Smith { 1881bef8e0ddSBarry Smith int ierr,(*f)(Mat,int *); 1882bef8e0ddSBarry Smith 1883bef8e0ddSBarry Smith PetscFunctionBegin; 1884bef8e0ddSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1885bef8e0ddSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f); 1886bef8e0ddSBarry Smith CHKERRQ(ierr); 1887bef8e0ddSBarry Smith if (f) { 1888bef8e0ddSBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 1889bef8e0ddSBarry Smith } else { 1890bef8e0ddSBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 1891bef8e0ddSBarry Smith } 1892bef8e0ddSBarry Smith PetscFunctionReturn(0); 1893bef8e0ddSBarry Smith } 1894bef8e0ddSBarry Smith 1895be6bf707SBarry Smith /* ----------------------------------------------------------------------------------------*/ 1896be6bf707SBarry Smith 1897fb2e594dSBarry Smith EXTERN_C_BEGIN 1898be6bf707SBarry Smith #undef __FUNC__ 1899be6bf707SBarry Smith #define __FUNC__ "MatStoreValues_SeqAIJ" 1900be6bf707SBarry Smith int MatStoreValues_SeqAIJ(Mat mat) 1901be6bf707SBarry Smith { 1902be6bf707SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1903be6bf707SBarry Smith int nz = aij->i[aij->m]+aij->indexshift; 1904be6bf707SBarry Smith 1905be6bf707SBarry Smith PetscFunctionBegin; 1906be6bf707SBarry Smith if (aij->nonew != 1) { 1907be6bf707SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1908be6bf707SBarry Smith } 1909be6bf707SBarry Smith 1910be6bf707SBarry Smith /* allocate space for values if not already there */ 1911be6bf707SBarry Smith if (!aij->saved_values) { 1912be6bf707SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 1913be6bf707SBarry Smith } 1914be6bf707SBarry Smith 1915be6bf707SBarry Smith /* copy values over */ 1916be6bf707SBarry Smith PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar)); 1917be6bf707SBarry Smith PetscFunctionReturn(0); 1918be6bf707SBarry Smith } 1919fb2e594dSBarry Smith EXTERN_C_END 1920be6bf707SBarry Smith 1921be6bf707SBarry Smith #undef __FUNC__ 1922be6bf707SBarry Smith #define __FUNC__ "MatStoreValues" 1923be6bf707SBarry Smith /*@ 1924be6bf707SBarry Smith MatStoreValues - Stashes a copy of the matrix values; this allows, for 1925be6bf707SBarry Smith example, reuse of the linear part of a Jacobian, while recomputing the 1926be6bf707SBarry Smith nonlinear portion. 1927be6bf707SBarry Smith 1928be6bf707SBarry Smith Collect on Mat 1929be6bf707SBarry Smith 1930be6bf707SBarry Smith Input Parameters: 1931be6bf707SBarry Smith . mat - the matrix (currently on AIJ matrices support this option) 1932be6bf707SBarry Smith 1933be6bf707SBarry Smith Common Usage, with SNESSolve(): 1934be6bf707SBarry Smith $ Create Jacobian matrix 1935be6bf707SBarry Smith $ Set linear terms into matrix 1936be6bf707SBarry Smith $ Apply boundary conditions to matrix, at this time matrix must have 1937be6bf707SBarry Smith $ final nonzero structure (i.e. setting the nonlinear terms and applying 1938be6bf707SBarry Smith $ boundary conditions again will not change the nonzero structure 1939be6bf707SBarry Smith $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); 1940be6bf707SBarry Smith $ ierr = MatStoreValues(mat); 1941be6bf707SBarry Smith $ Call SNESSetJacobian() with matrix 1942be6bf707SBarry Smith $ In your Jacobian routine 1943be6bf707SBarry Smith $ ierr = MatRetrieveValues(mat); 1944be6bf707SBarry Smith $ Set nonlinear terms in matrix 1945be6bf707SBarry Smith 1946be6bf707SBarry Smith Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself: 1947be6bf707SBarry Smith $ // build linear portion of Jacobian 1948be6bf707SBarry Smith $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); 1949be6bf707SBarry Smith $ ierr = MatStoreValues(mat); 1950be6bf707SBarry Smith $ loop over nonlinear iterations 1951be6bf707SBarry Smith $ ierr = MatRetrieveValues(mat); 1952be6bf707SBarry Smith $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian 1953be6bf707SBarry Smith $ // call MatAssemblyBegin/End() on matrix 1954be6bf707SBarry Smith $ Solve linear system with Jacobian 1955be6bf707SBarry Smith $ endloop 1956be6bf707SBarry Smith 1957be6bf707SBarry Smith Notes: 1958be6bf707SBarry Smith Matrix must already be assemblied before calling this routine 1959be6bf707SBarry Smith Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 1960be6bf707SBarry Smith calling this routine. 1961be6bf707SBarry Smith 1962be6bf707SBarry Smith .seealso: MatRetrieveValues() 1963be6bf707SBarry Smith 1964be6bf707SBarry Smith @*/ 1965be6bf707SBarry Smith int MatStoreValues(Mat mat) 1966be6bf707SBarry Smith { 1967be6bf707SBarry Smith int ierr,(*f)(Mat); 1968be6bf707SBarry Smith 1969be6bf707SBarry Smith PetscFunctionBegin; 1970be6bf707SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1971be6bf707SBarry Smith if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1972be6bf707SBarry Smith if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1973be6bf707SBarry Smith 1974be6bf707SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void **)&f); 1975be6bf707SBarry Smith CHKERRQ(ierr); 1976be6bf707SBarry Smith if (f) { 1977be6bf707SBarry Smith ierr = (*f)(mat);CHKERRQ(ierr); 1978be6bf707SBarry Smith } else { 1979be6bf707SBarry Smith SETERRQ(1,1,"Wrong type of matrix to store values"); 1980be6bf707SBarry Smith } 1981be6bf707SBarry Smith PetscFunctionReturn(0); 1982be6bf707SBarry Smith } 1983be6bf707SBarry Smith 1984fb2e594dSBarry Smith EXTERN_C_BEGIN 1985be6bf707SBarry Smith #undef __FUNC__ 1986be6bf707SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqAIJ" 1987be6bf707SBarry Smith int MatRetrieveValues_SeqAIJ(Mat mat) 1988be6bf707SBarry Smith { 1989be6bf707SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1990be6bf707SBarry Smith int nz = aij->i[aij->m]+aij->indexshift; 1991be6bf707SBarry Smith 1992be6bf707SBarry Smith PetscFunctionBegin; 1993be6bf707SBarry Smith if (aij->nonew != 1) { 1994be6bf707SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1995be6bf707SBarry Smith } 1996be6bf707SBarry Smith if (!aij->saved_values) { 1997be6bf707SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 1998be6bf707SBarry Smith } 1999be6bf707SBarry Smith 2000be6bf707SBarry Smith /* copy values over */ 2001be6bf707SBarry Smith PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar)); 2002be6bf707SBarry Smith PetscFunctionReturn(0); 2003be6bf707SBarry Smith } 2004fb2e594dSBarry Smith EXTERN_C_END 2005be6bf707SBarry Smith 2006be6bf707SBarry Smith #undef __FUNC__ 2007be6bf707SBarry Smith #define __FUNC__ "MatRetrieveValues" 2008be6bf707SBarry Smith /*@ 2009be6bf707SBarry Smith MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 2010be6bf707SBarry Smith example, reuse of the linear part of a Jacobian, while recomputing the 2011be6bf707SBarry Smith nonlinear portion. 2012be6bf707SBarry Smith 2013be6bf707SBarry Smith Collect on Mat 2014be6bf707SBarry Smith 2015be6bf707SBarry Smith Input Parameters: 2016be6bf707SBarry Smith . mat - the matrix (currently on AIJ matrices support this option) 2017be6bf707SBarry Smith 2018be6bf707SBarry Smith .seealso: MatStoreValues() 2019be6bf707SBarry Smith 2020be6bf707SBarry Smith @*/ 2021be6bf707SBarry Smith int MatRetrieveValues(Mat mat) 2022be6bf707SBarry Smith { 2023be6bf707SBarry Smith int ierr,(*f)(Mat); 2024be6bf707SBarry Smith 2025be6bf707SBarry Smith PetscFunctionBegin; 2026be6bf707SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 2027be6bf707SBarry Smith if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2028be6bf707SBarry Smith if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2029be6bf707SBarry Smith 2030be6bf707SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void **)&f); 2031be6bf707SBarry Smith CHKERRQ(ierr); 2032be6bf707SBarry Smith if (f) { 2033be6bf707SBarry Smith ierr = (*f)(mat);CHKERRQ(ierr); 2034be6bf707SBarry Smith } else { 2035be6bf707SBarry Smith SETERRQ(1,1,"Wrong type of matrix to retrieve values"); 2036be6bf707SBarry Smith } 2037be6bf707SBarry Smith PetscFunctionReturn(0); 2038be6bf707SBarry Smith } 2039be6bf707SBarry Smith 2040be6bf707SBarry Smith /* --------------------------------------------------------------------------------*/ 2041cb5b572fSBarry Smith 2042bef8e0ddSBarry Smith #undef __FUNC__ 20435615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ" 204417ab2063SBarry Smith /*@C 2045682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 20460d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 20476e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 20482bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 20492bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 205017ab2063SBarry Smith 2051db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2052db81eaa0SLois Curfman McInnes 205317ab2063SBarry Smith Input Parameters: 2054db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 205517ab2063SBarry Smith . m - number of rows 205617ab2063SBarry Smith . n - number of columns 205717ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 2058db81eaa0SLois Curfman McInnes - nzz - array containing the number of nonzeros in the various rows 20592bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 206017ab2063SBarry Smith 206117ab2063SBarry Smith Output Parameter: 2062416022c9SBarry Smith . A - the matrix 206317ab2063SBarry Smith 2064b259b22eSLois Curfman McInnes Notes: 206517ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 206617ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 20670002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 206844cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 206917ab2063SBarry Smith 207017ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2071a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 20723d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 20736da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 207417ab2063SBarry Smith 2075682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 20764fca80b9SLois Curfman McInnes improve numerical efficiency of matrix-vector products and solves. We 2077682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 20786c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 20796c7ebb05SLois Curfman McInnes 20806c7ebb05SLois Curfman McInnes Options Database Keys: 2081db81eaa0SLois Curfman McInnes + -mat_aij_no_inode - Do not use inodes 2082db81eaa0SLois Curfman McInnes . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5) 2083db81eaa0SLois Curfman McInnes - -mat_aij_oneindex - Internally use indexing starting at 1 2084db81eaa0SLois Curfman McInnes rather than 0. Note that when calling MatSetValues(), 2085db81eaa0SLois Curfman McInnes the user still MUST index entries starting at 0! 208617ab2063SBarry Smith 2087bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices() 208817ab2063SBarry Smith @*/ 2089416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 209017ab2063SBarry Smith { 2091416022c9SBarry Smith Mat B; 2092416022c9SBarry Smith Mat_SeqAIJ *b; 20936945ee14SBarry Smith int i, len, ierr, flg,size; 20946945ee14SBarry Smith 20953a40ed3dSBarry Smith PetscFunctionBegin; 20966945ee14SBarry Smith MPI_Comm_size(comm,&size); 2097a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1"); 2098d5d45c9bSBarry Smith 2099416022c9SBarry Smith *A = 0; 2100*3f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,"Mat",comm,MatDestroy,MatView); 2101416022c9SBarry Smith PLogObjectCreate(B); 21020452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 210344cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 21040a6ffc59SBarry Smith PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)); 2105e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqAIJ; 2106e1311b90SBarry Smith B->ops->view = MatView_SeqAIJ; 2107416022c9SBarry Smith B->factor = 0; 2108416022c9SBarry Smith B->lupivotthreshold = 1.0; 210990f02eecSBarry Smith B->mapping = 0; 2110e1311b90SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg);CHKERRQ(ierr); 21117a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 2112e1311b90SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*)&b->ilu_preserve_row_sums);CHKERRQ(ierr); 2113416022c9SBarry Smith b->row = 0; 2114416022c9SBarry Smith b->col = 0; 211582bf6240SBarry Smith b->icol = 0; 2116416022c9SBarry Smith b->indexshift = 0; 2117b810aeb4SBarry Smith b->reallocs = 0; 211869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 211969957df2SSatish Balay if (flg) b->indexshift = -1; 212017ab2063SBarry Smith 212144cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 212244cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 2123a5ae1ecdSBarry Smith 21244d197716SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 21254d197716SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 2126a5ae1ecdSBarry Smith 21270452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 2128b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 21297b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 21307b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 2131416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 213217ab2063SBarry Smith nz = nz*m; 21333a40ed3dSBarry Smith } else { 213417ab2063SBarry Smith nz = 0; 2135416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 213617ab2063SBarry Smith } 213717ab2063SBarry Smith 213817ab2063SBarry Smith /* allocate the matrix space */ 2139416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 21400452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 2141416022c9SBarry Smith b->j = (int *) (b->a + nz); 2142cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 2143416022c9SBarry Smith b->i = b->j + nz; 2144416022c9SBarry Smith b->singlemalloc = 1; 214517ab2063SBarry Smith 2146416022c9SBarry Smith b->i[0] = -b->indexshift; 214717ab2063SBarry Smith for (i=1; i<m+1; i++) { 2148416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 214917ab2063SBarry Smith } 215017ab2063SBarry Smith 2151416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 21520452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 2153f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 2154416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 215517ab2063SBarry Smith 2156416022c9SBarry Smith b->nz = 0; 2157416022c9SBarry Smith b->maxnz = nz; 2158416022c9SBarry Smith b->sorted = 0; 2159416022c9SBarry Smith b->roworiented = 1; 2160416022c9SBarry Smith b->nonew = 0; 2161416022c9SBarry Smith b->diag = 0; 2162416022c9SBarry Smith b->solve_work = 0; 216371bd300dSLois Curfman McInnes b->spptr = 0; 2164754ec7b1SSatish Balay b->inode.node_count = 0; 2165754ec7b1SSatish Balay b->inode.size = 0; 21666c7ebb05SLois Curfman McInnes b->inode.limit = 5; 21676c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 2168be6bf707SBarry Smith b->saved_values = 0; 21694e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 217017ab2063SBarry Smith 2171416022c9SBarry Smith *A = B; 21724e220ebcSLois Curfman McInnes 21734b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 21744b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 217569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 217669957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 21774b14c69eSBarry Smith #endif 217869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 217969957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 218069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 218169957df2SSatish Balay if (flg) { 2182a8c6a408SBarry Smith if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml"); 2183416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 218417ab2063SBarry Smith } 218569957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 218669957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 2187bef8e0ddSBarry Smith 2188bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C", 2189bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 2190bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 2191be6bf707SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 2192be6bf707SBarry Smith "MatStoreValues_SeqAIJ", 2193be6bf707SBarry Smith (void*)MatStoreValues_SeqAIJ);CHKERRQ(ierr); 2194be6bf707SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 2195be6bf707SBarry Smith "MatRetrieveValues_SeqAIJ", 2196be6bf707SBarry Smith (void*)MatRetrieveValues_SeqAIJ);CHKERRQ(ierr); 21973a40ed3dSBarry Smith PetscFunctionReturn(0); 219817ab2063SBarry Smith } 219917ab2063SBarry Smith 22005615d1e5SSatish Balay #undef __FUNC__ 2201be6bf707SBarry Smith #define __FUNC__ "MatDuplicate_SeqAIJ" 2202be6bf707SBarry Smith int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 220317ab2063SBarry Smith { 2204416022c9SBarry Smith Mat C; 2205416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 2206bef8e0ddSBarry Smith int i,len, m = a->m,shift = a->indexshift,ierr; 220717ab2063SBarry Smith 22083a40ed3dSBarry Smith PetscFunctionBegin; 22094043dd9cSLois Curfman McInnes *B = 0; 2210*3f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,"Mat",A->comm,MatDestroy,MatView); 2211416022c9SBarry Smith PLogObjectCreate(C); 22120452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 2213f830108cSBarry Smith PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 2214e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqAIJ; 2215e1311b90SBarry Smith C->ops->view = MatView_SeqAIJ; 2216416022c9SBarry Smith C->factor = A->factor; 2217416022c9SBarry Smith c->row = 0; 2218416022c9SBarry Smith c->col = 0; 221982bf6240SBarry Smith c->icol = 0; 2220416022c9SBarry Smith c->indexshift = shift; 2221c456f294SBarry Smith C->assembled = PETSC_TRUE; 222217ab2063SBarry Smith 222344cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 222444cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 222544cd7ae7SLois Curfman McInnes C->M = a->m; 222644cd7ae7SLois Curfman McInnes C->N = a->n; 222717ab2063SBarry Smith 22280452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 22290452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 223017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 2231416022c9SBarry Smith c->imax[i] = a->imax[i]; 2232416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 223317ab2063SBarry Smith } 223417ab2063SBarry Smith 223517ab2063SBarry Smith /* allocate the matrix space */ 2236416022c9SBarry Smith c->singlemalloc = 1; 2237416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 22380452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 2239416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 2240416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 2241416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 224217ab2063SBarry Smith if (m > 0) { 2243416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 2244be6bf707SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 2245416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 2246be6bf707SBarry Smith } else { 2247be6bf707SBarry Smith PetscMemzero(c->a,(a->i[m]+shift)*sizeof(Scalar)); 224817ab2063SBarry Smith } 224908480c60SBarry Smith } 225017ab2063SBarry Smith 2251f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 2252416022c9SBarry Smith c->sorted = a->sorted; 2253416022c9SBarry Smith c->roworiented = a->roworiented; 2254416022c9SBarry Smith c->nonew = a->nonew; 22557a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 2256be6bf707SBarry Smith c->saved_values = 0; 225717ab2063SBarry Smith 2258416022c9SBarry Smith if (a->diag) { 22590452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 2260416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 226117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 2262416022c9SBarry Smith c->diag[i] = a->diag[i]; 226317ab2063SBarry Smith } 22643a40ed3dSBarry Smith } else c->diag = 0; 22656c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 22666c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 2267754ec7b1SSatish Balay if (a->inode.size){ 2268daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 2269754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 2270daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 2271754ec7b1SSatish Balay } else { 2272754ec7b1SSatish Balay c->inode.size = 0; 2273754ec7b1SSatish Balay c->inode.node_count = 0; 2274754ec7b1SSatish Balay } 2275416022c9SBarry Smith c->nz = a->nz; 2276416022c9SBarry Smith c->maxnz = a->maxnz; 2277416022c9SBarry Smith c->solve_work = 0; 227876dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 2279754ec7b1SSatish Balay 2280416022c9SBarry Smith *B = C; 2281bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)C,"MatSeqAIJSetColumnIndices_C", 2282bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 2283bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 22843a40ed3dSBarry Smith PetscFunctionReturn(0); 228517ab2063SBarry Smith } 228617ab2063SBarry Smith 22875615d1e5SSatish Balay #undef __FUNC__ 22885615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ" 228919bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 229017ab2063SBarry Smith { 2291416022c9SBarry Smith Mat_SeqAIJ *a; 2292416022c9SBarry Smith Mat B; 229317699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 2294bcd2baecSBarry Smith MPI_Comm comm; 229517ab2063SBarry Smith 22963a40ed3dSBarry Smith PetscFunctionBegin; 229719bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 229817699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 2299a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor"); 230090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 23010752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 2302a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file"); 230317ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 230417ab2063SBarry Smith 2305d64ed03dSBarry Smith if (nz < 0) { 2306a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ"); 2307d64ed03dSBarry Smith } 2308d64ed03dSBarry Smith 230917ab2063SBarry Smith /* read in row lengths */ 23100452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 23110752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 231217ab2063SBarry Smith 231317ab2063SBarry Smith /* create our matrix */ 2314416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 2315416022c9SBarry Smith B = *A; 2316416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 2317416022c9SBarry Smith shift = a->indexshift; 231817ab2063SBarry Smith 231917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 23200752156aSBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr); 232117ab2063SBarry Smith if (shift) { 232217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 2323416022c9SBarry Smith a->j[i] += 1; 232417ab2063SBarry Smith } 232517ab2063SBarry Smith } 232617ab2063SBarry Smith 232717ab2063SBarry Smith /* read in nonzero values */ 23280752156aSBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr); 232917ab2063SBarry Smith 233017ab2063SBarry Smith /* set matrix "i" values */ 2331416022c9SBarry Smith a->i[0] = -shift; 233217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 2333416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 2334416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 233517ab2063SBarry Smith } 23360452661fSBarry Smith PetscFree(rowlengths); 233717ab2063SBarry Smith 23386d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 23396d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 23403a40ed3dSBarry Smith PetscFunctionReturn(0); 234117ab2063SBarry Smith } 234217ab2063SBarry Smith 23435615d1e5SSatish Balay #undef __FUNC__ 23445615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ" 23458f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 23467264ac53SSatish Balay { 23477264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 23487264ac53SSatish Balay 23493a40ed3dSBarry Smith PetscFunctionBegin; 2350a8c6a408SBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 23517264ac53SSatish Balay 23527264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 23537264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 2354bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 23553a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2356bcd2baecSBarry Smith } 23577264ac53SSatish Balay 23587264ac53SSatish Balay /* if the a->i are the same */ 23598108c231SLois Curfman McInnes if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) { 23603a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 23617264ac53SSatish Balay } 23627264ac53SSatish Balay 23637264ac53SSatish Balay /* if a->j are the same */ 2364bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 23653a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2366bcd2baecSBarry Smith } 2367bcd2baecSBarry Smith 2368bcd2baecSBarry Smith /* if a->a are the same */ 236919bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 23703a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 23717264ac53SSatish Balay } 237277c4ece6SBarry Smith *flg = PETSC_TRUE; 23733a40ed3dSBarry Smith PetscFunctionReturn(0); 23747264ac53SSatish Balay 23757264ac53SSatish Balay } 2376