1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*4787f768SSatish Balay static char vcid[] = "$Id: aij.c,v 1.297 1999/02/11 19:51:52 bsmith Exp balay $"; 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; 273f1db9ecSBarry 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]; 1595ef9f2a5SBarry Smith if (row < 0) continue; 1603a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 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 */ 1675ef9f2a5SBarry Smith if (in[l] < 0) continue; 1683a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 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) { 1735ef9f2a5SBarry Smith value = v[l + k*n]; 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); 3390ef38995SBarry Smith if (flg1 || flg2) ViewerASCIIPrintf(viewer," not using I-node routines\n"); 3400ef38995SBarry Smith else ViewerASCIIPrintf(viewer," using I-node routines: found %d nodes, limit used is %d\n",a->inode.node_count,a->inode.limit); 3413a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 342d00d2cf4SBarry Smith int nofinalvalue = 0; 343d00d2cf4SBarry Smith if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != a->n-!shift)) { 344d00d2cf4SBarry Smith nofinalvalue = 1; 345d00d2cf4SBarry Smith } 346416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3474e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 348d00d2cf4SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",a->nz+nofinalvalue); 34917ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 35017ab2063SBarry Smith 35117ab2063SBarry Smith for (i=0; i<m; i++) { 352416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 3533a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 354e20fef11SSatish 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])); 35517ab2063SBarry Smith #else 3567a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 35717ab2063SBarry Smith #endif 35817ab2063SBarry Smith } 35917ab2063SBarry Smith } 360d00d2cf4SBarry Smith if (nofinalvalue) { 361d00d2cf4SBarry Smith fprintf(fd,"%d %d %18.16e\n", m, a->n, 0.0); 362d00d2cf4SBarry Smith } 363e24b481bSBarry Smith if (outputname) fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 364e24b481bSBarry Smith else fprintf(fd,"];\n M = spconvert(zzz);\n"); 3653a40ed3dSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 36644cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 36744cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 36844cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 3693a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 370e20fef11SSatish Balay if (PetscImaginary(a->a[j]) > 0.0 && PetscReal(a->a[j]) != 0.0) 371e20fef11SSatish Balay fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j])); 372e20fef11SSatish Balay else if (PetscImaginary(a->a[j]) < 0.0 && PetscReal(a->a[j]) != 0.0) 373e20fef11SSatish Balay fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j])); 374e20fef11SSatish Balay else if (PetscReal(a->a[j]) != 0.0) 375e20fef11SSatish Balay fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j])); 37644cd7ae7SLois Curfman McInnes #else 37744cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 37844cd7ae7SLois Curfman McInnes #endif 37944cd7ae7SLois Curfman McInnes } 38044cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 38144cd7ae7SLois Curfman McInnes } 38202594712SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_SYMMODU) { 383496be53dSLois Curfman McInnes int nzd=0, fshift=1, *sptr; 3842e44a96cSLois Curfman McInnes sptr = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(sptr); 385496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 386496be53dSLois Curfman McInnes sptr[i] = nzd+1; 387496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 388496be53dSLois Curfman McInnes if (a->j[j] >= i) { 3893a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 390e20fef11SSatish Balay if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) nzd++; 391496be53dSLois Curfman McInnes #else 392496be53dSLois Curfman McInnes if (a->a[j] != 0.0) nzd++; 393496be53dSLois Curfman McInnes #endif 394496be53dSLois Curfman McInnes } 395496be53dSLois Curfman McInnes } 396496be53dSLois Curfman McInnes } 3972e44a96cSLois Curfman McInnes sptr[m] = nzd+1; 398496be53dSLois Curfman McInnes fprintf(fd," %d %d\n\n",m,nzd); 3992e44a96cSLois Curfman McInnes for ( i=0; i<m+1; i+=6 ) { 4002e44a96cSLois 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]); 4012e44a96cSLois 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]); 4022e44a96cSLois 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]); 4032e44a96cSLois Curfman McInnes else if (i+1<m) fprintf(fd," %d %d %d\n",sptr[i],sptr[i+1],sptr[i+2]); 4042e44a96cSLois Curfman McInnes else if (i<m) fprintf(fd," %d %d\n",sptr[i],sptr[i+1]); 4057272d637SLois Curfman McInnes else fprintf(fd," %d\n",sptr[i]); 406496be53dSLois Curfman McInnes } 407496be53dSLois Curfman McInnes fprintf(fd,"\n"); 408496be53dSLois Curfman McInnes PetscFree(sptr); 409496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 410496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 411496be53dSLois Curfman McInnes if (a->j[j] >= i) fprintf(fd," %d ",a->j[j]+fshift); 412496be53dSLois Curfman McInnes } 413496be53dSLois Curfman McInnes fprintf(fd,"\n"); 414496be53dSLois Curfman McInnes } 415496be53dSLois Curfman McInnes fprintf(fd,"\n"); 416496be53dSLois Curfman McInnes for ( i=0; i<m; i++ ) { 417496be53dSLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 418496be53dSLois Curfman McInnes if (a->j[j] >= i) { 4193a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 420e20fef11SSatish Balay if (PetscImaginary(a->a[j]) != 0.0 || PetscReal(a->a[j]) != 0.0) 421e20fef11SSatish Balay fprintf(fd," %18.16e %18.16e ",PetscReal(a->a[j]),PetscImaginary(a->a[j])); 422496be53dSLois Curfman McInnes #else 423496be53dSLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %18.16e ",a->a[j]); 424496be53dSLois Curfman McInnes #endif 425496be53dSLois Curfman McInnes } 426496be53dSLois Curfman McInnes } 427496be53dSLois Curfman McInnes fprintf(fd,"\n"); 428496be53dSLois Curfman McInnes } 42902594712SBarry Smith } else if (format == VIEWER_FORMAT_ASCII_DENSE) { 43002594712SBarry Smith int cnt = 0,jcnt; 43102594712SBarry Smith Scalar value; 43202594712SBarry Smith 43302594712SBarry Smith for ( i=0; i<m; i++ ) { 43402594712SBarry Smith jcnt = 0; 43502594712SBarry Smith for ( j=0; j<a->n; j++ ) { 436e24b481bSBarry Smith if ( jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) { 43702594712SBarry Smith value = a->a[cnt++]; 438e24b481bSBarry Smith jcnt++; 43902594712SBarry Smith } else { 44002594712SBarry Smith value = 0.0; 44102594712SBarry Smith } 44202594712SBarry Smith #if defined(USE_PETSC_COMPLEX) 44302594712SBarry Smith fprintf(fd," %7.5e+%7.5e i ",PetscReal(value),PetscImaginary(value)); 44402594712SBarry Smith #else 44502594712SBarry Smith fprintf(fd," %7.5e ",value); 44602594712SBarry Smith #endif 44702594712SBarry Smith } 44802594712SBarry Smith fprintf(fd,"\n"); 44902594712SBarry Smith } 4503a40ed3dSBarry Smith } else { 45117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 45217ab2063SBarry Smith fprintf(fd,"row %d:",i); 453416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 4543a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 455e20fef11SSatish Balay if (PetscImaginary(a->a[j]) > 0.0) { 456e20fef11SSatish Balay fprintf(fd," %d %g + %g i",a->j[j]+shift,PetscReal(a->a[j]),PetscImaginary(a->a[j])); 457e20fef11SSatish Balay } else if (PetscImaginary(a->a[j]) < 0.0) { 458e20fef11SSatish Balay fprintf(fd," %d %g - %g i",a->j[j]+shift,PetscReal(a->a[j]),-PetscImaginary(a->a[j])); 4593a40ed3dSBarry Smith } else { 460e20fef11SSatish Balay fprintf(fd," %d %g ",a->j[j]+shift,PetscReal(a->a[j])); 46117ab2063SBarry Smith } 46217ab2063SBarry Smith #else 463416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 46417ab2063SBarry Smith #endif 46517ab2063SBarry Smith } 46617ab2063SBarry Smith fprintf(fd,"\n"); 46717ab2063SBarry Smith } 46817ab2063SBarry Smith } 46917ab2063SBarry Smith fflush(fd); 4703a40ed3dSBarry Smith PetscFunctionReturn(0); 471416022c9SBarry Smith } 472416022c9SBarry Smith 4735615d1e5SSatish Balay #undef __FUNC__ 474480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw_Zoom" 475480ef9eaSBarry Smith int MatView_SeqAIJ_Draw_Zoom(Draw draw,void *Aa) 476416022c9SBarry Smith { 477480ef9eaSBarry Smith Mat A = (Mat) Aa; 478416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 47977ed5343SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,color,rank; 4800513a670SBarry Smith int format; 481480ef9eaSBarry Smith double xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0; 482480ef9eaSBarry Smith Viewer viewer; 48377ed5343SBarry Smith MPI_Comm comm; 484cddf8d76SBarry Smith 4853a40ed3dSBarry Smith PetscFunctionBegin; 48677ed5343SBarry Smith /* 48777ed5343SBarry Smith This is nasty. If this is called from an originally parallel matrix 48877ed5343SBarry Smith then all processes call this, but only the first has the matrix so the 48977ed5343SBarry Smith rest should return immediately. 49077ed5343SBarry Smith */ 49177ed5343SBarry Smith ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr); 49277ed5343SBarry Smith MPI_Comm_rank(comm,&rank); 49377ed5343SBarry Smith if (rank) PetscFunctionReturn(0); 49477ed5343SBarry Smith 495480ef9eaSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*) &viewer);CHKERRQ(ierr); 4960513a670SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 49719bcc07fSBarry Smith 498480ef9eaSBarry Smith ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr); CHKERRQ(ierr); 499416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 5000513a670SBarry Smith 5010513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 5020513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 503cddf8d76SBarry Smith color = DRAW_BLUE; 504416022c9SBarry Smith for ( i=0; i<m; i++ ) { 505cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 506416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 507cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5083a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 509e20fef11SSatish Balay if (PetscReal(a->a[j]) >= 0.) continue; 510cddf8d76SBarry Smith #else 511cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 512cddf8d76SBarry Smith #endif 513480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 514cddf8d76SBarry Smith } 515cddf8d76SBarry Smith } 516cddf8d76SBarry Smith color = DRAW_CYAN; 517cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 518cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 519cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 520cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 521cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 522480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 523cddf8d76SBarry Smith } 524cddf8d76SBarry Smith } 525cddf8d76SBarry Smith color = DRAW_RED; 526cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 527cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 528cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 529cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5303a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 531e20fef11SSatish Balay if (PetscReal(a->a[j]) <= 0.) continue; 532cddf8d76SBarry Smith #else 533cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 534cddf8d76SBarry Smith #endif 535480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 536416022c9SBarry Smith } 537416022c9SBarry Smith } 5380513a670SBarry Smith } else { 5390513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 5400513a670SBarry Smith /* first determine max of all nonzero values */ 5410513a670SBarry Smith int nz = a->nz,count; 5420513a670SBarry Smith Draw popup; 543480ef9eaSBarry Smith double scale; 5440513a670SBarry Smith 5450513a670SBarry Smith for ( i=0; i<nz; i++ ) { 5460513a670SBarry Smith if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 5470513a670SBarry Smith } 548480ef9eaSBarry Smith scale = (245.0 - DRAW_BASIC_COLORS)/maxv; 549522c5e43SBarry Smith ierr = DrawGetPopup(draw,&popup); CHKERRQ(ierr); 5500513a670SBarry Smith ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr); 5510513a670SBarry Smith count = 0; 5520513a670SBarry Smith for ( i=0; i<m; i++ ) { 5530513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 5540513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 5550513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 5566d6b6c30SSatish Balay color = DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(a->a[count])); 557480ef9eaSBarry Smith ierr = DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 5580513a670SBarry Smith count++; 5590513a670SBarry Smith } 5600513a670SBarry Smith } 5610513a670SBarry Smith } 562480ef9eaSBarry Smith PetscFunctionReturn(0); 563480ef9eaSBarry Smith } 564cddf8d76SBarry Smith 565480ef9eaSBarry Smith #undef __FUNC__ 566480ef9eaSBarry Smith #define __FUNC__ "MatView_SeqAIJ_Draw" 567480ef9eaSBarry Smith int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 568480ef9eaSBarry Smith { 569480ef9eaSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 570480ef9eaSBarry Smith int ierr; 571480ef9eaSBarry Smith Draw draw; 572480ef9eaSBarry Smith double xr,yr,xl,yl,h,w; 573480ef9eaSBarry Smith 574480ef9eaSBarry Smith PetscTruth isnull; 575480ef9eaSBarry Smith 576480ef9eaSBarry Smith PetscFunctionBegin; 57777ed5343SBarry Smith ierr = ViewerDrawGetDraw(viewer,0,&draw); CHKERRQ(ierr); 578480ef9eaSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); 579480ef9eaSBarry Smith if (isnull) PetscFunctionReturn(0); 580480ef9eaSBarry Smith 581480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 582480ef9eaSBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 583480ef9eaSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 584cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 585480ef9eaSBarry Smith ierr = DrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A); CHKERRQ(ierr); 586480ef9eaSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 5873a40ed3dSBarry Smith PetscFunctionReturn(0); 588416022c9SBarry Smith } 589416022c9SBarry Smith 5905615d1e5SSatish Balay #undef __FUNC__ 591d4bb536fSBarry Smith #define __FUNC__ "MatView_SeqAIJ" 592e1311b90SBarry Smith int MatView_SeqAIJ(Mat A,Viewer viewer) 593416022c9SBarry Smith { 594416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 595bcd2baecSBarry Smith ViewerType vtype; 596bcd2baecSBarry Smith int ierr; 597416022c9SBarry Smith 5983a40ed3dSBarry Smith PetscFunctionBegin; 599bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 6007b2a1423SBarry Smith if (PetscTypeCompare(vtype,SOCKET_VIEWER)) { 6017b2a1423SBarry Smith ierr = ViewerSocketPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j);CHKERRQ(ierr); 6023f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,ASCII_VIEWER)){ 6033a40ed3dSBarry Smith ierr = MatView_SeqAIJ_ASCII(A,viewer); CHKERRQ(ierr); 6043f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,BINARY_VIEWER)) { 6053a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Binary(A,viewer); CHKERRQ(ierr); 6063f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,DRAW_VIEWER)) { 6073a40ed3dSBarry Smith ierr = MatView_SeqAIJ_Draw(A,viewer); CHKERRQ(ierr); 6085cd90555SBarry Smith } else { 6095cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported by PETSc object"); 61017ab2063SBarry Smith } 6113a40ed3dSBarry Smith PetscFunctionReturn(0); 61217ab2063SBarry Smith } 61319bcc07fSBarry Smith 614c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 6155615d1e5SSatish Balay #undef __FUNC__ 6165615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ" 6178f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 61817ab2063SBarry Smith { 619416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 62041c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 62143ee02c3SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift,rmax = 0; 622416022c9SBarry Smith Scalar *aa = a->a, *ap; 62317ab2063SBarry Smith 6243a40ed3dSBarry Smith PetscFunctionBegin; 6253a40ed3dSBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); 62617ab2063SBarry Smith 62743ee02c3SBarry Smith if (m) rmax = ailen[0]; /* determine row with most nonzeros */ 62817ab2063SBarry Smith for ( i=1; i<m; i++ ) { 629416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 63017ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 63194a9d846SBarry Smith rmax = PetscMax(rmax,ailen[i]); 63217ab2063SBarry Smith if (fshift) { 633416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 63417ab2063SBarry Smith N = ailen[i]; 63517ab2063SBarry Smith for ( j=0; j<N; j++ ) { 63617ab2063SBarry Smith ip[j-fshift] = ip[j]; 63717ab2063SBarry Smith ap[j-fshift] = ap[j]; 63817ab2063SBarry Smith } 63917ab2063SBarry Smith } 64017ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 64117ab2063SBarry Smith } 64217ab2063SBarry Smith if (m) { 64317ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 64417ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 64517ab2063SBarry Smith } 64617ab2063SBarry Smith /* reset ilen and imax for each row */ 64717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 64817ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 64917ab2063SBarry Smith } 650416022c9SBarry Smith a->nz = ai[m] + shift; 65117ab2063SBarry Smith 65217ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 653416022c9SBarry Smith if (fshift && a->diag) { 6540452661fSBarry Smith PetscFree(a->diag); 655416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 656416022c9SBarry Smith a->diag = 0; 65717ab2063SBarry Smith } 6584e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 6594e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 6604e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 661b810aeb4SBarry Smith a->reallocs); 66294a9d846SBarry Smith PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Most nonzeros in any row is %d\n",rmax); 663dd5f02e7SSatish Balay a->reallocs = 0; 6644e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 6654e220ebcSLois Curfman McInnes 66676dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 66741c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 6683a40ed3dSBarry Smith PetscFunctionReturn(0); 66917ab2063SBarry Smith } 67017ab2063SBarry Smith 6715615d1e5SSatish Balay #undef __FUNC__ 6725615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ" 6738f6be9afSLois Curfman McInnes int MatZeroEntries_SeqAIJ(Mat A) 67417ab2063SBarry Smith { 675416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6763a40ed3dSBarry Smith 6773a40ed3dSBarry Smith PetscFunctionBegin; 678cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 6793a40ed3dSBarry Smith PetscFunctionReturn(0); 68017ab2063SBarry Smith } 681416022c9SBarry Smith 6825615d1e5SSatish Balay #undef __FUNC__ 6835615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ" 684e1311b90SBarry Smith int MatDestroy_SeqAIJ(Mat A) 68517ab2063SBarry Smith { 686416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 68782bf6240SBarry Smith int ierr; 688d5d45c9bSBarry Smith 6893a40ed3dSBarry Smith PetscFunctionBegin; 69094d884c6SBarry Smith if (--A->refct > 0) PetscFunctionReturn(0); 69194d884c6SBarry Smith 69294d884c6SBarry Smith if (A->mapping) { 69394d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 69494d884c6SBarry Smith } 69594d884c6SBarry Smith if (A->bmapping) { 69694d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->bmapping); CHKERRQ(ierr); 69794d884c6SBarry Smith } 69861b13de0SBarry Smith if (A->rmap) { 69961b13de0SBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 70061b13de0SBarry Smith } 70161b13de0SBarry Smith if (A->cmap) { 70261b13de0SBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 70361b13de0SBarry Smith } 7043a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 705e1311b90SBarry Smith PLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 70617ab2063SBarry Smith #endif 7070452661fSBarry Smith PetscFree(a->a); 7080452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 7090452661fSBarry Smith if (a->diag) PetscFree(a->diag); 7100452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 7110452661fSBarry Smith if (a->imax) PetscFree(a->imax); 7120452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 71376dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 71482bf6240SBarry Smith if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} 715be6bf707SBarry Smith if (a->saved_values) PetscFree(a->saved_values); 7160452661fSBarry Smith PetscFree(a); 717eed86810SBarry Smith 718f2655603SLois Curfman McInnes PLogObjectDestroy(A); 719f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 7203a40ed3dSBarry Smith PetscFunctionReturn(0); 72117ab2063SBarry Smith } 72217ab2063SBarry Smith 7235615d1e5SSatish Balay #undef __FUNC__ 7245615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ" 7258f6be9afSLois Curfman McInnes int MatCompress_SeqAIJ(Mat A) 72617ab2063SBarry Smith { 7273a40ed3dSBarry Smith PetscFunctionBegin; 7283a40ed3dSBarry Smith PetscFunctionReturn(0); 72917ab2063SBarry Smith } 73017ab2063SBarry Smith 7315615d1e5SSatish Balay #undef __FUNC__ 7325615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ" 7338f6be9afSLois Curfman McInnes int MatSetOption_SeqAIJ(Mat A,MatOption op) 73417ab2063SBarry Smith { 735416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7363a40ed3dSBarry Smith 7373a40ed3dSBarry Smith PetscFunctionBegin; 7386d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 7396d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 7406d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 741219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 7426d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 743*4787f768SSatish Balay else if (op == MAT_NEW_NONZERO_LOCATION_ERR) a->nonew = -1; 744*4787f768SSatish Balay else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR) a->nonew = -2; 7456d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 7466d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 747219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 7486d4a8577SBarry Smith op == MAT_SYMMETRIC || 7496d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 75090f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 751b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES|| 752b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) 753981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqAIJ:Option ignored\n"); 7543a40ed3dSBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) { 7553a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS"); 7563a40ed3dSBarry Smith } else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 7576d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 7586d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 7596d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 7606d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 7613a40ed3dSBarry Smith else SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 7623a40ed3dSBarry Smith PetscFunctionReturn(0); 76317ab2063SBarry Smith } 76417ab2063SBarry Smith 7655615d1e5SSatish Balay #undef __FUNC__ 7665615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ" 7678f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 76817ab2063SBarry Smith { 769416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 7703a40ed3dSBarry Smith int i,j, n,shift = a->indexshift,ierr; 77117ab2063SBarry Smith Scalar *x, zero = 0.0; 77217ab2063SBarry Smith 7733a40ed3dSBarry Smith PetscFunctionBegin; 7743a40ed3dSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 775e1311b90SBarry Smith ierr = VecGetArray(v,&x);;CHKERRQ(ierr); 776e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);;CHKERRQ(ierr); 777a8c6a408SBarry Smith if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector"); 778416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 779416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 780416022c9SBarry Smith if (a->j[j]+shift == i) { 781416022c9SBarry Smith x[i] = a->a[j]; 78217ab2063SBarry Smith break; 78317ab2063SBarry Smith } 78417ab2063SBarry Smith } 78517ab2063SBarry Smith } 786e1311b90SBarry Smith ierr = VecRestoreArray(v,&x);;CHKERRQ(ierr); 7873a40ed3dSBarry Smith PetscFunctionReturn(0); 78817ab2063SBarry Smith } 78917ab2063SBarry Smith 79017ab2063SBarry Smith /* -------------------------------------------------------*/ 79117ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 79217ab2063SBarry Smith /* -------------------------------------------------------*/ 7935615d1e5SSatish Balay #undef __FUNC__ 7945615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ" 79544cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 79617ab2063SBarry Smith { 797416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 79817ab2063SBarry Smith Scalar *x, *y, *v, alpha; 799e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx, shift = a->indexshift; 80017ab2063SBarry Smith 8013a40ed3dSBarry Smith PetscFunctionBegin; 802e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 803e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 804cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 80517ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 80617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 807416022c9SBarry Smith idx = a->j + a->i[i] + shift; 808416022c9SBarry Smith v = a->a + a->i[i] + shift; 809416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 81017ab2063SBarry Smith alpha = x[i]; 81117ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 81217ab2063SBarry Smith } 813416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 814e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 815e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8163a40ed3dSBarry Smith PetscFunctionReturn(0); 81717ab2063SBarry Smith } 818d5d45c9bSBarry Smith 8195615d1e5SSatish Balay #undef __FUNC__ 8205615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ" 82144cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 82217ab2063SBarry Smith { 823416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 82417ab2063SBarry Smith Scalar *x, *y, *v, alpha; 825e1311b90SBarry Smith int ierr,m = a->m, n, i, *idx,shift = a->indexshift; 82617ab2063SBarry Smith 8273a40ed3dSBarry Smith PetscFunctionBegin; 8282e8a6d31SBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 829e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 830e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 83117ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 83217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 833416022c9SBarry Smith idx = a->j + a->i[i] + shift; 834416022c9SBarry Smith v = a->a + a->i[i] + shift; 835416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 83617ab2063SBarry Smith alpha = x[i]; 83717ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 83817ab2063SBarry Smith } 83990f02eecSBarry Smith PLogFlops(2*a->nz); 840e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 841e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8423a40ed3dSBarry Smith PetscFunctionReturn(0); 84317ab2063SBarry Smith } 84417ab2063SBarry Smith 8455615d1e5SSatish Balay #undef __FUNC__ 8465615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ" 84744cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 84817ab2063SBarry Smith { 849416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 85017ab2063SBarry Smith Scalar *x, *y, *v, sum; 851e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 8522e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTAIJ) 853e36a17ebSSatish Balay int n, i, jrow,j; 854e36a17ebSSatish Balay #endif 85517ab2063SBarry Smith 856fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT) 857fee21e36SBarry Smith #pragma disjoint(*x,*y,*v) 858fee21e36SBarry Smith #endif 859fee21e36SBarry Smith 8603a40ed3dSBarry Smith PetscFunctionBegin; 861e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 862e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 86317ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 864416022c9SBarry Smith idx = a->j; 865d64ed03dSBarry Smith v = a->a; 866416022c9SBarry Smith ii = a->i; 867acc4d009SSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTAIJ) 86829b5ca8aSSatish Balay fortranmultaij_(&m,x,ii,idx+shift,v+shift,y); 8698d195f9aSBarry Smith #else 870d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 871d64ed03dSBarry Smith idx += shift; 87217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 8739ea0dfa2SSatish Balay jrow = ii[i]; 8749ea0dfa2SSatish Balay n = ii[i+1] - jrow; 87517ab2063SBarry Smith sum = 0.0; 8769ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 8779ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 8789ea0dfa2SSatish Balay } 87917ab2063SBarry Smith y[i] = sum; 88017ab2063SBarry Smith } 8818d195f9aSBarry Smith #endif 882416022c9SBarry Smith PLogFlops(2*a->nz - m); 883e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 884e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 8853a40ed3dSBarry Smith PetscFunctionReturn(0); 88617ab2063SBarry Smith } 88717ab2063SBarry Smith 8885615d1e5SSatish Balay #undef __FUNC__ 8895615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ" 89044cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 89117ab2063SBarry Smith { 892416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 89317ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 894e1311b90SBarry Smith int ierr,m = a->m, *idx, shift = a->indexshift,*ii; 8952e40c62eSSatish Balay #if !defined(USE_FORTRAN_KERNEL_MULTADDAIJ) 896e36a17ebSSatish Balay int n,i,jrow,j; 897e36a17ebSSatish Balay #endif 8989ea0dfa2SSatish Balay 8993a40ed3dSBarry Smith PetscFunctionBegin; 900e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 901e1311b90SBarry Smith ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 9022e8a6d31SBarry Smith if (zz != yy) { 903e1311b90SBarry Smith ierr = VecGetArray(zz,&z); CHKERRQ(ierr); 9042e8a6d31SBarry Smith } else { 9052e8a6d31SBarry Smith z = y; 9062e8a6d31SBarry Smith } 90717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 908cddf8d76SBarry Smith idx = a->j; 909d64ed03dSBarry Smith v = a->a; 910cddf8d76SBarry Smith ii = a->i; 9112e40c62eSSatish Balay #if defined(USE_FORTRAN_KERNEL_MULTADDAIJ) 91229b5ca8aSSatish Balay fortranmultaddaij_(&m,x,ii,idx+shift,v+shift,y,z); 91302ab625aSSatish Balay #else 914d64ed03dSBarry Smith v += shift; /* shift for Fortran start by 1 indexing */ 915d64ed03dSBarry Smith idx += shift; 91617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 9179ea0dfa2SSatish Balay jrow = ii[i]; 9189ea0dfa2SSatish Balay n = ii[i+1] - jrow; 91917ab2063SBarry Smith sum = y[i]; 9209ea0dfa2SSatish Balay for ( j=0; j<n; j++) { 9219ea0dfa2SSatish Balay sum += v[jrow]*x[idx[jrow]]; jrow++; 9229ea0dfa2SSatish Balay } 92317ab2063SBarry Smith z[i] = sum; 92417ab2063SBarry Smith } 92502ab625aSSatish Balay #endif 926416022c9SBarry Smith PLogFlops(2*a->nz); 927e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 928e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y); CHKERRQ(ierr); 9292e8a6d31SBarry Smith if (zz != yy) { 930e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z); CHKERRQ(ierr); 9312e8a6d31SBarry Smith } 9323a40ed3dSBarry Smith PetscFunctionReturn(0); 93317ab2063SBarry Smith } 93417ab2063SBarry Smith 93517ab2063SBarry Smith /* 93617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 93717ab2063SBarry Smith */ 9385615d1e5SSatish Balay #undef __FUNC__ 9395615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ" 940416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 94117ab2063SBarry Smith { 942416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 943416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 94417ab2063SBarry Smith 9453a40ed3dSBarry Smith PetscFunctionBegin; 9460452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 947416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 948416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 94935b0346bSBarry Smith diag[i] = a->i[i+1]; 950416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 951416022c9SBarry Smith if (a->j[j]+shift == i) { 95217ab2063SBarry Smith diag[i] = j - shift; 95317ab2063SBarry Smith break; 95417ab2063SBarry Smith } 95517ab2063SBarry Smith } 95617ab2063SBarry Smith } 957416022c9SBarry Smith a->diag = diag; 9583a40ed3dSBarry Smith PetscFunctionReturn(0); 95917ab2063SBarry Smith } 96017ab2063SBarry Smith 961be5855fcSBarry Smith /* 962be5855fcSBarry Smith Checks for missing diagonals 963be5855fcSBarry Smith */ 964be5855fcSBarry Smith #undef __FUNC__ 965be5855fcSBarry Smith #define __FUNC__ "MatMissingDiag_SeqAIJ" 966be5855fcSBarry Smith int MatMissingDiag_SeqAIJ(Mat A) 967be5855fcSBarry Smith { 968be5855fcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 969be5855fcSBarry Smith int *diag = a->diag, *jj = a->j,i,shift = a->indexshift; 970be5855fcSBarry Smith 971be5855fcSBarry Smith PetscFunctionBegin; 972be5855fcSBarry Smith for ( i=0; i<a->m; i++ ) { 973be5855fcSBarry Smith if (jj[diag[i]+shift] != i-shift) { 974be5855fcSBarry Smith SETERRQ1(1,1,"Matrix is missing diagonal number %d",i); 975be5855fcSBarry Smith } 976be5855fcSBarry Smith } 977be5855fcSBarry Smith PetscFunctionReturn(0); 978be5855fcSBarry Smith } 979be5855fcSBarry Smith 9805615d1e5SSatish Balay #undef __FUNC__ 9815615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ" 982be5855fcSBarry Smith int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag,double fshift,int its,Vec xx) 98317ab2063SBarry Smith { 984416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 985416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 986d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 98717ab2063SBarry Smith 9883a40ed3dSBarry Smith PetscFunctionBegin; 989e1311b90SBarry Smith ierr = VecGetArray(xx,&x); CHKERRQ(ierr); 990fb2e594dSBarry Smith if (xx != bb) { 991e1311b90SBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 992fb2e594dSBarry Smith } else { 993fb2e594dSBarry Smith b = x; 994fb2e594dSBarry Smith } 995fb2e594dSBarry Smith 996d00d2cf4SBarry Smith if (!a->diag) {ierr = MatMarkDiag_SeqAIJ(A);CHKERRQ(ierr);} 997416022c9SBarry Smith diag = a->diag; 998416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 99917ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 100017ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 100117ab2063SBarry Smith bs = b + shift; 100217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1003416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1004416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1005488ecbafSBarry Smith PLogFlops(2*n-1); 1006416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1007416022c9SBarry Smith v = a->a + diag[i] + (!shift); 100817ab2063SBarry Smith sum = b[i]*d/omega; 100917ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 101017ab2063SBarry Smith x[i] = sum; 101117ab2063SBarry Smith } 1012cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1013fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 10143a40ed3dSBarry Smith PetscFunctionReturn(0); 101517ab2063SBarry Smith } 101617ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 1017a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"SOR_APPLY_LOWER is not done"); 10183a40ed3dSBarry Smith } else if (flag & SOR_EISENSTAT) { 101917ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 102017ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 102117ab2063SBarry Smith 102217ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 102317ab2063SBarry Smith 102417ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 102517ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 102617ab2063SBarry Smith is the relaxation factor. 102717ab2063SBarry Smith */ 10280452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 102917ab2063SBarry Smith scale = (2.0/omega) - 1.0; 103017ab2063SBarry Smith 103117ab2063SBarry Smith /* x = (E + U)^{-1} b */ 103217ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1033416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1034416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1035488ecbafSBarry Smith PLogFlops(2*n-1); 1036416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1037416022c9SBarry Smith v = a->a + diag[i] + (!shift); 103817ab2063SBarry Smith sum = b[i]; 103917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 104017ab2063SBarry Smith x[i] = omega*(sum/d); 104117ab2063SBarry Smith } 104217ab2063SBarry Smith 104317ab2063SBarry Smith /* t = b - (2*E - D)x */ 1044416022c9SBarry Smith v = a->a; 104517ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 1046488ecbafSBarry Smith PLogFlops(2*m); 1047488ecbafSBarry Smith 104817ab2063SBarry Smith 104917ab2063SBarry Smith /* t = (E + L)^{-1}t */ 1050416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 1051416022c9SBarry Smith diag = a->diag; 105217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1053416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1054416022c9SBarry Smith n = diag[i] - a->i[i]; 1055488ecbafSBarry Smith PLogFlops(2*n-1); 1056416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1057416022c9SBarry Smith v = a->a + a->i[i] + shift; 105817ab2063SBarry Smith sum = t[i]; 105917ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 106017ab2063SBarry Smith t[i] = omega*(sum/d); 106117ab2063SBarry Smith } 106217ab2063SBarry Smith 106317ab2063SBarry Smith /* x = x + t */ 106417ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 10650452661fSBarry Smith PetscFree(t); 1066cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1067fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 10683a40ed3dSBarry Smith PetscFunctionReturn(0); 106917ab2063SBarry Smith } 107017ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 107117ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 107217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1073416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1074416022c9SBarry Smith n = diag[i] - a->i[i]; 1075488ecbafSBarry Smith PLogFlops(2*n-1); 1076416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1077416022c9SBarry Smith v = a->a + a->i[i] + shift; 107817ab2063SBarry Smith sum = b[i]; 107917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 108017ab2063SBarry Smith x[i] = omega*(sum/d); 108117ab2063SBarry Smith } 108217ab2063SBarry Smith xb = x; 10833a40ed3dSBarry Smith } else xb = b; 108417ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 108517ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 108617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1087416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 108817ab2063SBarry Smith } 1089488ecbafSBarry Smith PLogFlops(m); 109017ab2063SBarry Smith } 109117ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 109217ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1093416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1094416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 1095488ecbafSBarry Smith PLogFlops(2*n-1); 1096416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 1097416022c9SBarry Smith v = a->a + diag[i] + (!shift); 109817ab2063SBarry Smith sum = xb[i]; 109917ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 110017ab2063SBarry Smith x[i] = omega*(sum/d); 110117ab2063SBarry Smith } 110217ab2063SBarry Smith } 110317ab2063SBarry Smith its--; 110417ab2063SBarry Smith } 110517ab2063SBarry Smith while (its--) { 110617ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 110717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1108416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 1109416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1110488ecbafSBarry Smith PLogFlops(2*n-1); 1111416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1112416022c9SBarry Smith v = a->a + a->i[i] + shift; 111317ab2063SBarry Smith sum = b[i]; 111417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 11157e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 111617ab2063SBarry Smith } 111717ab2063SBarry Smith } 111817ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 111917ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 1120416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 1121416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 1122488ecbafSBarry Smith PLogFlops(2*n-1); 1123416022c9SBarry Smith idx = a->j + a->i[i] + shift; 1124416022c9SBarry Smith v = a->a + a->i[i] + shift; 112517ab2063SBarry Smith sum = b[i]; 112617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 11277e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 112817ab2063SBarry Smith } 112917ab2063SBarry Smith } 113017ab2063SBarry Smith } 1131cb5b572fSBarry Smith ierr = VecRestoreArray(xx,&x); CHKERRQ(ierr); 1132fb2e594dSBarry Smith if (bb != xx) {ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);} 11333a40ed3dSBarry Smith PetscFunctionReturn(0); 113417ab2063SBarry Smith } 113517ab2063SBarry Smith 11365615d1e5SSatish Balay #undef __FUNC__ 11375615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ" 11388f6be9afSLois Curfman McInnes int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 113917ab2063SBarry Smith { 1140416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11414e220ebcSLois Curfman McInnes 11423a40ed3dSBarry Smith PetscFunctionBegin; 11434e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 11444e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 11454e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 11464e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 11474e220ebcSLois Curfman McInnes info->block_size = 1.0; 11484e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 11494e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 11504e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 11514e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 11524e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 11534e220ebcSLois Curfman McInnes info->memory = A->mem; 11544e220ebcSLois Curfman McInnes if (A->factor) { 11554e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 11564e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 11574e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 11584e220ebcSLois Curfman McInnes } else { 11594e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 11604e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 11614e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 11624e220ebcSLois Curfman McInnes } 11633a40ed3dSBarry Smith PetscFunctionReturn(0); 116417ab2063SBarry Smith } 116517ab2063SBarry Smith 116617ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 116717ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 116817ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 116917ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 117017ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 117117ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 117217ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 117317ab2063SBarry Smith 11745615d1e5SSatish Balay #undef __FUNC__ 11755615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ" 11768f6be9afSLois Curfman McInnes int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 117717ab2063SBarry Smith { 1178416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1179416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 118017ab2063SBarry Smith 11813a40ed3dSBarry Smith PetscFunctionBegin; 118277c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 118317ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 118417ab2063SBarry Smith if (diag) { 118517ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1186a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 11877ae801bdSBarry Smith if (a->ilen[rows[i]] > 0) { 1188416022c9SBarry Smith a->ilen[rows[i]] = 1; 1189416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 1190416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 11917ae801bdSBarry Smith } else { /* in case row was completely empty */ 1192d64ed03dSBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES);CHKERRQ(ierr); 119317ab2063SBarry Smith } 119417ab2063SBarry Smith } 11953a40ed3dSBarry Smith } else { 119617ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1197a8c6a408SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"row out of range"); 1198416022c9SBarry Smith a->ilen[rows[i]] = 0; 119917ab2063SBarry Smith } 120017ab2063SBarry Smith } 12017ae801bdSBarry Smith ierr = ISRestoreIndices(is,&rows);CHKERRQ(ierr); 120243a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 12033a40ed3dSBarry Smith PetscFunctionReturn(0); 120417ab2063SBarry Smith } 120517ab2063SBarry Smith 12065615d1e5SSatish Balay #undef __FUNC__ 12075615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ" 12088f6be9afSLois Curfman McInnes int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 120917ab2063SBarry Smith { 1210416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 12113a40ed3dSBarry Smith 12123a40ed3dSBarry Smith PetscFunctionBegin; 12130752156aSBarry Smith if (m) *m = a->m; 12140752156aSBarry Smith if (n) *n = a->n; 12153a40ed3dSBarry Smith PetscFunctionReturn(0); 121617ab2063SBarry Smith } 121717ab2063SBarry Smith 12185615d1e5SSatish Balay #undef __FUNC__ 12195615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ" 12208f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 122117ab2063SBarry Smith { 1222416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 12233a40ed3dSBarry Smith 12243a40ed3dSBarry Smith PetscFunctionBegin; 1225416022c9SBarry Smith *m = 0; *n = a->m; 12263a40ed3dSBarry Smith PetscFunctionReturn(0); 122717ab2063SBarry Smith } 1228026e39d0SSatish Balay 12295615d1e5SSatish Balay #undef __FUNC__ 12305615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ" 12314e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 123217ab2063SBarry Smith { 1233416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1234c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 123517ab2063SBarry Smith 12363a40ed3dSBarry Smith PetscFunctionBegin; 1237a8c6a408SBarry Smith if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range"); 123817ab2063SBarry Smith 1239416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 1240416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 124117ab2063SBarry Smith if (idx) { 1242416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 12434e093b46SBarry Smith if (*nz && shift) { 12440452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 124517ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 12464e093b46SBarry Smith } else if (*nz) { 12474e093b46SBarry Smith *idx = itmp; 124817ab2063SBarry Smith } 124917ab2063SBarry Smith else *idx = 0; 125017ab2063SBarry Smith } 12513a40ed3dSBarry Smith PetscFunctionReturn(0); 125217ab2063SBarry Smith } 125317ab2063SBarry Smith 12545615d1e5SSatish Balay #undef __FUNC__ 12555615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ" 12564e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 125717ab2063SBarry Smith { 12584e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 12593a40ed3dSBarry Smith 12603a40ed3dSBarry Smith PetscFunctionBegin; 12614e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 12623a40ed3dSBarry Smith PetscFunctionReturn(0); 126317ab2063SBarry Smith } 126417ab2063SBarry Smith 12655615d1e5SSatish Balay #undef __FUNC__ 12665615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ" 12678f6be9afSLois Curfman McInnes int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 126817ab2063SBarry Smith { 1269416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1270416022c9SBarry Smith Scalar *v = a->a; 127117ab2063SBarry Smith double sum = 0.0; 1272416022c9SBarry Smith int i, j,shift = a->indexshift; 127317ab2063SBarry Smith 12743a40ed3dSBarry Smith PetscFunctionBegin; 127517ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1276416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 12773a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 1278e20fef11SSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 127917ab2063SBarry Smith #else 128017ab2063SBarry Smith sum += (*v)*(*v); v++; 128117ab2063SBarry Smith #endif 128217ab2063SBarry Smith } 128317ab2063SBarry Smith *norm = sqrt(sum); 12843a40ed3dSBarry Smith } else if (type == NORM_1) { 128517ab2063SBarry Smith double *tmp; 1286416022c9SBarry Smith int *jj = a->j; 128766963ce1SSatish Balay tmp = (double *) PetscMalloc( (a->n+1)*sizeof(double) ); CHKPTRQ(tmp); 1288cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 128917ab2063SBarry Smith *norm = 0.0; 1290416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1291a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 129217ab2063SBarry Smith } 1293416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 129417ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 129517ab2063SBarry Smith } 12960452661fSBarry Smith PetscFree(tmp); 12973a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 129817ab2063SBarry Smith *norm = 0.0; 1299416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1300416022c9SBarry Smith v = a->a + a->i[j] + shift; 130117ab2063SBarry Smith sum = 0.0; 1302416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1303cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 130417ab2063SBarry Smith } 130517ab2063SBarry Smith if (sum > *norm) *norm = sum; 130617ab2063SBarry Smith } 13073a40ed3dSBarry Smith } else { 1308a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No support for two norm"); 130917ab2063SBarry Smith } 13103a40ed3dSBarry Smith PetscFunctionReturn(0); 131117ab2063SBarry Smith } 131217ab2063SBarry Smith 13135615d1e5SSatish Balay #undef __FUNC__ 13145615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ" 13158f6be9afSLois Curfman McInnes int MatTranspose_SeqAIJ(Mat A,Mat *B) 131617ab2063SBarry Smith { 1317416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1318416022c9SBarry Smith Mat C; 1319416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1320416022c9SBarry Smith int shift = a->indexshift; 1321d5d45c9bSBarry Smith Scalar *array = a->a; 132217ab2063SBarry Smith 13233a40ed3dSBarry Smith PetscFunctionBegin; 1324a8c6a408SBarry Smith if (B == PETSC_NULL && m != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Square matrix only for in-place"); 13250452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1326cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 132717ab2063SBarry Smith if (shift) { 132817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 132917ab2063SBarry Smith } 133017ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1331416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 13320452661fSBarry Smith PetscFree(col); 133317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 133417ab2063SBarry Smith len = ai[i+1]-ai[i]; 1335416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 133617ab2063SBarry Smith array += len; aj += len; 133717ab2063SBarry Smith } 133817ab2063SBarry Smith if (shift) { 133917ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 134017ab2063SBarry Smith } 134117ab2063SBarry Smith 13426d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13436d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 134417ab2063SBarry Smith 13453638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1346416022c9SBarry Smith *B = C; 134717ab2063SBarry Smith } else { 1348f830108cSBarry Smith PetscOps *Abops; 13490a6ffc59SBarry Smith MatOps Aops; 1350f830108cSBarry Smith 1351416022c9SBarry Smith /* This isn't really an in-place transpose */ 13520452661fSBarry Smith PetscFree(a->a); 13530452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 13540452661fSBarry Smith if (a->diag) PetscFree(a->diag); 13550452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 13560452661fSBarry Smith if (a->imax) PetscFree(a->imax); 13570452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 13581073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 13590452661fSBarry Smith PetscFree(a); 1360f830108cSBarry Smith 1361488ecbafSBarry Smith 1362488ecbafSBarry Smith ierr = MapDestroy(A->rmap);CHKERRQ(ierr); 1363488ecbafSBarry Smith ierr = MapDestroy(A->cmap);CHKERRQ(ierr); 1364488ecbafSBarry Smith 1365f830108cSBarry Smith /* 1366f830108cSBarry Smith This is horrible, horrible code. We need to keep the 13678f0f457cSSatish Balay the bops and ops Structures, copy everything from C 13688f0f457cSSatish Balay including the function pointers.. 1369f830108cSBarry Smith */ 13708f0f457cSSatish Balay PetscMemcpy(A->ops,C->ops,sizeof(struct _MatOps)); 13718f0f457cSSatish Balay PetscMemcpy(A->bops,C->bops,sizeof(PetscOps)); 1372f830108cSBarry Smith Abops = A->bops; 1373f830108cSBarry Smith Aops = A->ops; 1374f09e8eb9SSatish Balay PetscMemcpy(A,C,sizeof(struct _p_Mat)); 1375f830108cSBarry Smith A->bops = Abops; 1376f830108cSBarry Smith A->ops = Aops; 137727a8da17SBarry Smith A->qlist = 0; 1378f830108cSBarry Smith 13790452661fSBarry Smith PetscHeaderDestroy(C); 138017ab2063SBarry Smith } 13813a40ed3dSBarry Smith PetscFunctionReturn(0); 138217ab2063SBarry Smith } 138317ab2063SBarry Smith 13845615d1e5SSatish Balay #undef __FUNC__ 13855615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ" 13868f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 138717ab2063SBarry Smith { 1388416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 138917ab2063SBarry Smith Scalar *l,*r,x,*v; 1390e1311b90SBarry Smith int ierr,i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 139117ab2063SBarry Smith 13923a40ed3dSBarry Smith PetscFunctionBegin; 139317ab2063SBarry Smith if (ll) { 13943ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 13953ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 1396e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr); 1397a8c6a408SBarry Smith if (m != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length"); 1398e1311b90SBarry Smith ierr = VecGetArray(ll,&l); CHKERRQ(ierr); 1399416022c9SBarry Smith v = a->a; 140017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 140117ab2063SBarry Smith x = l[i]; 1402416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 140317ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 140417ab2063SBarry Smith } 1405e1311b90SBarry Smith ierr = VecRestoreArray(ll,&l); CHKERRQ(ierr); 140644cd7ae7SLois Curfman McInnes PLogFlops(nz); 140717ab2063SBarry Smith } 140817ab2063SBarry Smith if (rr) { 1409e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&n);CHKERRQ(ierr); 1410a8c6a408SBarry Smith if (n != a->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length"); 1411e1311b90SBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1412416022c9SBarry Smith v = a->a; jj = a->j; 141317ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 141417ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 141517ab2063SBarry Smith } 1416e1311b90SBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 141744cd7ae7SLois Curfman McInnes PLogFlops(nz); 141817ab2063SBarry Smith } 14193a40ed3dSBarry Smith PetscFunctionReturn(0); 142017ab2063SBarry Smith } 142117ab2063SBarry Smith 14225615d1e5SSatish Balay #undef __FUNC__ 14235615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ" 14247b2a1423SBarry Smith int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,int csize,MatReuse scall,Mat *B) 142517ab2063SBarry Smith { 1426db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 1427d5db1b72SBarry Smith int *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 142899141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1429a2744918SBarry Smith register int sum,lensi; 143099141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 143199141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 143299141d43SSatish Balay Scalar *a_new,*mat_a; 1433416022c9SBarry Smith Mat C; 1434fee21e36SBarry Smith PetscTruth stride; 143517ab2063SBarry Smith 14363a40ed3dSBarry Smith PetscFunctionBegin; 1437d64ed03dSBarry Smith ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr); 1438a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted"); 1439d64ed03dSBarry Smith ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 1440a8c6a408SBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted"); 144199141d43SSatish Balay 144217ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 144317ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 144417ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 144517ab2063SBarry Smith 1446fee21e36SBarry Smith ierr = ISStrideGetInfo(iscol,&first,&step); CHKERRQ(ierr); 1447fee21e36SBarry Smith ierr = ISStride(iscol,&stride); CHKERRQ(ierr); 1448fee21e36SBarry Smith if (stride && step == 1) { 144902834360SBarry Smith /* special case of contiguous rows */ 145057faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 145102834360SBarry Smith starts = lens + ncols; 145202834360SBarry Smith /* loop over new rows determining lens and starting points */ 145302834360SBarry Smith for (i=0; i<nrows; i++) { 1454a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1455a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 145602834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1457d8ced48eSBarry Smith if (aj[k]+shift >= first) { 145802834360SBarry Smith starts[i] = k; 145902834360SBarry Smith break; 146002834360SBarry Smith } 146102834360SBarry Smith } 1462a2744918SBarry Smith sum = 0; 146302834360SBarry Smith while (k < kend) { 1464d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1465a2744918SBarry Smith sum++; 146602834360SBarry Smith } 1467a2744918SBarry Smith lens[i] = sum; 146802834360SBarry Smith } 146902834360SBarry Smith /* create submatrix */ 1470cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 147108480c60SBarry Smith int n_cols,n_rows; 147208480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1473a8c6a408SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); 1474d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 147508480c60SBarry Smith C = *B; 14763a40ed3dSBarry Smith } else { 147702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 147808480c60SBarry Smith } 1479db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1480db02288aSLois Curfman McInnes 148102834360SBarry Smith /* loop over rows inserting into submatrix */ 1482db02288aSLois Curfman McInnes a_new = c->a; 1483db02288aSLois Curfman McInnes j_new = c->j; 1484db02288aSLois Curfman McInnes i_new = c->i; 1485db02288aSLois Curfman McInnes i_new[0] = -shift; 148602834360SBarry Smith for (i=0; i<nrows; i++) { 1487a2744918SBarry Smith ii = starts[i]; 1488a2744918SBarry Smith lensi = lens[i]; 1489a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1490a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 149102834360SBarry Smith } 1492a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1493a2744918SBarry Smith a_new += lensi; 1494a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1495a2744918SBarry Smith c->ilen[i] = lensi; 149602834360SBarry Smith } 14970452661fSBarry Smith PetscFree(lens); 14983a40ed3dSBarry Smith } else { 149902834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 15000452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 150102834360SBarry Smith ssmap = smap + shift; 150299141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1503cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 150417ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 150502834360SBarry Smith /* determine lens of each row */ 150602834360SBarry Smith for (i=0; i<nrows; i++) { 1507d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 150802834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 150902834360SBarry Smith lens[i] = 0; 151002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1511d8ced48eSBarry Smith if (ssmap[aj[k]]) { 151202834360SBarry Smith lens[i]++; 151302834360SBarry Smith } 151402834360SBarry Smith } 151502834360SBarry Smith } 151617ab2063SBarry Smith /* Create and fill new matrix */ 1517a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 151899141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 1519a8c6a408SBarry Smith if (c->m != nrows || c->n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size"); 152099141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1521a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros"); 152299141d43SSatish Balay } 152399141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 152408480c60SBarry Smith C = *B; 15253a40ed3dSBarry Smith } else { 152602834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 152708480c60SBarry Smith } 152899141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 152917ab2063SBarry Smith for (i=0; i<nrows; i++) { 153099141d43SSatish Balay row = irow[i]; 153199141d43SSatish Balay kstart = ai[row]+shift; 153299141d43SSatish Balay kend = kstart + a->ilen[row]; 153399141d43SSatish Balay mat_i = c->i[i]+shift; 153499141d43SSatish Balay mat_j = c->j + mat_i; 153599141d43SSatish Balay mat_a = c->a + mat_i; 153699141d43SSatish Balay mat_ilen = c->ilen + i; 153717ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 153899141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 153999141d43SSatish Balay *mat_j++ = tcol - (!shift); 154099141d43SSatish Balay *mat_a++ = a->a[k]; 154199141d43SSatish Balay (*mat_ilen)++; 154299141d43SSatish Balay 154317ab2063SBarry Smith } 154417ab2063SBarry Smith } 154517ab2063SBarry Smith } 154602834360SBarry Smith /* Free work space */ 154702834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 154899141d43SSatish Balay PetscFree(smap); PetscFree(lens); 154902834360SBarry Smith } 15506d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 15516d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 155217ab2063SBarry Smith 155317ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1554416022c9SBarry Smith *B = C; 15553a40ed3dSBarry Smith PetscFunctionReturn(0); 155617ab2063SBarry Smith } 155717ab2063SBarry Smith 1558a871dcd8SBarry Smith /* 155963b91edcSBarry Smith note: This can only work for identity for row and col. It would 156063b91edcSBarry Smith be good to check this and otherwise generate an error. 1561a871dcd8SBarry Smith */ 15625615d1e5SSatish Balay #undef __FUNC__ 15635615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ" 15645ef9f2a5SBarry Smith int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatILUInfo *info) 1565a871dcd8SBarry Smith { 156663b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 156708480c60SBarry Smith int ierr; 156863b91edcSBarry Smith Mat outA; 156963b91edcSBarry Smith 15703a40ed3dSBarry Smith PetscFunctionBegin; 15715ef9f2a5SBarry Smith if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels=0 supported"); 1572a871dcd8SBarry Smith 157363b91edcSBarry Smith outA = inA; 157463b91edcSBarry Smith inA->factor = FACTOR_LU; 157563b91edcSBarry Smith a->row = row; 157663b91edcSBarry Smith a->col = col; 157763b91edcSBarry Smith 1578f0ec6fceSSatish Balay /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */ 1579f0ec6fceSSatish Balay ierr = ISInvertPermutation(col,&(a->icol)); CHKERRQ(ierr); 15801e4dd532SSatish Balay PLogObjectParent(inA,a->icol); 1581f0ec6fceSSatish Balay 158294a9d846SBarry Smith if (!a->solve_work) { /* this matrix may have been factored before */ 15830452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar));CHKPTRQ(a->solve_work); 158494a9d846SBarry Smith } 158563b91edcSBarry Smith 158608480c60SBarry Smith if (!a->diag) { 158708480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 158863b91edcSBarry Smith } 158963b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 15903a40ed3dSBarry Smith PetscFunctionReturn(0); 1591a871dcd8SBarry Smith } 1592a871dcd8SBarry Smith 159374cdf0dfSBarry Smith #include "pinclude/blaslapack.h" 15945615d1e5SSatish Balay #undef __FUNC__ 15955615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ" 15968f6be9afSLois Curfman McInnes int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1597f0b747eeSBarry Smith { 1598f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1599f0b747eeSBarry Smith int one = 1; 16003a40ed3dSBarry Smith 16013a40ed3dSBarry Smith PetscFunctionBegin; 1602f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1603f0b747eeSBarry Smith PLogFlops(a->nz); 16043a40ed3dSBarry Smith PetscFunctionReturn(0); 1605f0b747eeSBarry Smith } 1606f0b747eeSBarry Smith 16075615d1e5SSatish Balay #undef __FUNC__ 16085615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ" 16097b2a1423SBarry Smith int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B) 1610cddf8d76SBarry Smith { 1611cddf8d76SBarry Smith int ierr,i; 1612cddf8d76SBarry Smith 16133a40ed3dSBarry Smith PetscFunctionBegin; 1614cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 16150452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1616cddf8d76SBarry Smith } 1617cddf8d76SBarry Smith 1618cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 16196a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1620cddf8d76SBarry Smith } 16213a40ed3dSBarry Smith PetscFunctionReturn(0); 1622cddf8d76SBarry Smith } 1623cddf8d76SBarry Smith 16245615d1e5SSatish Balay #undef __FUNC__ 16255615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ" 16268f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 16275a838052SSatish Balay { 1628f830108cSBarry Smith PetscFunctionBegin; 16295a838052SSatish Balay *bs = 1; 16303a40ed3dSBarry Smith PetscFunctionReturn(0); 16315a838052SSatish Balay } 16325a838052SSatish Balay 16335615d1e5SSatish Balay #undef __FUNC__ 16345615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ" 16358f6be9afSLois Curfman McInnes int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 16364dcbc457SBarry Smith { 1637e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 163806763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 16398a047759SSatish Balay int start, end, *ai, *aj; 1640bbd702dbSSatish Balay BT table; 1641bbd702dbSSatish Balay 16423a40ed3dSBarry Smith PetscFunctionBegin; 16438a047759SSatish Balay shift = a->indexshift; 1644e4d965acSSatish Balay m = a->m; 1645e4d965acSSatish Balay ai = a->i; 16468a047759SSatish Balay aj = a->j+shift; 16478a047759SSatish Balay 1648a8c6a408SBarry Smith if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"illegal overlap value used"); 164906763907SSatish Balay 165006763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 1651bbd702dbSSatish Balay ierr = BTCreate(m,table); CHKERRQ(ierr); 165206763907SSatish Balay 1653e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1654b97fc60eSLois Curfman McInnes /* Initialize the two local arrays */ 1655e4d965acSSatish Balay isz = 0; 1656bbd702dbSSatish Balay BTMemzero(m,table); 1657e4d965acSSatish Balay 1658e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 16594dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 166077c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1661e4d965acSSatish Balay 1662dd097bc3SLois Curfman McInnes /* Enter these into the temp arrays. I.e., mark table[row], enter row into new index */ 1663e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 1664bbd702dbSSatish Balay if(!BTLookupSet(table, idx[j])) { nidx[isz++] = idx[j];} 16654dcbc457SBarry Smith } 166606763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 166706763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1668e4d965acSSatish Balay 166904a348a9SBarry Smith k = 0; 167004a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap */ 167104a348a9SBarry Smith n = isz; 167206763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1673e4d965acSSatish Balay row = nidx[k]; 1674e4d965acSSatish Balay start = ai[row]; 1675e4d965acSSatish Balay end = ai[row+1]; 167604a348a9SBarry Smith for ( l = start; l<end ; l++){ 16778a047759SSatish Balay val = aj[l] + shift; 16782a8f2162SSatish Balay if (!BTLookupSet(table,val)) {nidx[isz++] = val;} 1679e4d965acSSatish Balay } 1680e4d965acSSatish Balay } 1681e4d965acSSatish Balay } 1682029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1683e4d965acSSatish Balay } 1684bbd702dbSSatish Balay BTDestroy(table); 168506763907SSatish Balay PetscFree(nidx); 16863a40ed3dSBarry Smith PetscFunctionReturn(0); 16874dcbc457SBarry Smith } 168817ab2063SBarry Smith 16890513a670SBarry Smith /* -------------------------------------------------------------- */ 16900513a670SBarry Smith #undef __FUNC__ 16910513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ" 16920513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 16930513a670SBarry Smith { 16940513a670SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 16950513a670SBarry Smith Scalar *vwork; 16960513a670SBarry Smith int i, ierr, nz, m = a->m, n = a->n, *cwork; 16970513a670SBarry Smith int *row,*col,*cnew,j,*lens; 169856cd22aeSBarry Smith IS icolp,irowp; 16990513a670SBarry Smith 17003a40ed3dSBarry Smith PetscFunctionBegin; 170156cd22aeSBarry Smith ierr = ISInvertPermutation(rowp,&irowp); CHKERRQ(ierr); 170256cd22aeSBarry Smith ierr = ISGetIndices(irowp,&row); CHKERRQ(ierr); 170356cd22aeSBarry Smith ierr = ISInvertPermutation(colp,&icolp); CHKERRQ(ierr); 170456cd22aeSBarry Smith ierr = ISGetIndices(icolp,&col); CHKERRQ(ierr); 17050513a670SBarry Smith 17060513a670SBarry Smith /* determine lengths of permuted rows */ 17070513a670SBarry Smith lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens); 17080513a670SBarry Smith for (i=0; i<m; i++ ) { 17090513a670SBarry Smith lens[row[i]] = a->i[i+1] - a->i[i]; 17100513a670SBarry Smith } 17110513a670SBarry Smith ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr); 17120513a670SBarry Smith PetscFree(lens); 17130513a670SBarry Smith 17140513a670SBarry Smith cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew); 17150513a670SBarry Smith for (i=0; i<m; i++) { 17160513a670SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 17170513a670SBarry Smith for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];} 17180513a670SBarry Smith ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr); 17190513a670SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 17200513a670SBarry Smith } 17210513a670SBarry Smith PetscFree(cnew); 17220513a670SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 17230513a670SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 172456cd22aeSBarry Smith ierr = ISRestoreIndices(irowp,&row); CHKERRQ(ierr); 172556cd22aeSBarry Smith ierr = ISRestoreIndices(icolp,&col); CHKERRQ(ierr); 172656cd22aeSBarry Smith ierr = ISDestroy(irowp); CHKERRQ(ierr); 172756cd22aeSBarry Smith ierr = ISDestroy(icolp); CHKERRQ(ierr); 17283a40ed3dSBarry Smith PetscFunctionReturn(0); 17290513a670SBarry Smith } 17300513a670SBarry Smith 17315615d1e5SSatish Balay #undef __FUNC__ 17325615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ" 1733682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1734682d7d0cSBarry Smith { 1735682d7d0cSBarry Smith static int called = 0; 1736682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1737682d7d0cSBarry Smith 17383a40ed3dSBarry Smith PetscFunctionBegin; 17393a40ed3dSBarry Smith if (called) {PetscFunctionReturn(0);} else called = 1; 174076be9ce4SBarry Smith (*PetscHelpPrintf)(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 174176be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_lu_pivotthreshold <threshold>: Set pivoting threshold\n"); 174276be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_oneindex: internal indices begin at 1 instead of the default 0.\n"); 174376be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_no_inode: Do not use inodes\n"); 174476be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_inode_limit <limit>: Set inode limit (max limit=5)\n"); 1745682d7d0cSBarry Smith #if defined(HAVE_ESSL) 174676be9ce4SBarry Smith (*PetscHelpPrintf)(comm," -mat_aij_essl: Use IBM sparse LU factorization and solve.\n"); 1747682d7d0cSBarry Smith #endif 17483a40ed3dSBarry Smith PetscFunctionReturn(0); 1749682d7d0cSBarry Smith } 17508f6be9afSLois Curfman McInnes extern int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1751a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1752a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1753a93ec695SBarry Smith 1754cb5b572fSBarry Smith #undef __FUNC__ 1755cb5b572fSBarry Smith #define __FUNC__ "MatCopy_SeqAIJ" 1756cb5b572fSBarry Smith int MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str) 1757cb5b572fSBarry Smith { 1758be6bf707SBarry Smith int ierr; 1759cb5b572fSBarry Smith 1760cb5b572fSBarry Smith PetscFunctionBegin; 1761be6bf707SBarry Smith if (str == SAME_NONZERO_PATTERN && B->type == MATSEQAIJ) { 1762be6bf707SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1763be6bf707SBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *) B->data; 1764be6bf707SBarry Smith 1765be6bf707SBarry Smith if (a->nonew != 1) SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1766be6bf707SBarry Smith if (b->nonew != 1) SETERRQ(1,1,"Must call MatSetOption(B,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1767be6bf707SBarry Smith if (a->i[a->m]+a->indexshift != b->i[b->m]+a->indexshift) { 1768be6bf707SBarry Smith SETERRQ(1,1,"Number of nonzeros in two matrices are different"); 1769cb5b572fSBarry Smith } 1770be6bf707SBarry Smith PetscMemcpy(b->a,a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 1771cb5b572fSBarry Smith } else { 1772cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1773cb5b572fSBarry Smith } 1774cb5b572fSBarry Smith PetscFunctionReturn(0); 1775cb5b572fSBarry Smith } 1776cb5b572fSBarry Smith 1777682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 17780a6ffc59SBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqAIJ, 1779cb5b572fSBarry Smith MatGetRow_SeqAIJ, 1780cb5b572fSBarry Smith MatRestoreRow_SeqAIJ, 1781cb5b572fSBarry Smith MatMult_SeqAIJ, 1782cb5b572fSBarry Smith MatMultAdd_SeqAIJ, 1783cb5b572fSBarry Smith MatMultTrans_SeqAIJ, 1784cb5b572fSBarry Smith MatMultTransAdd_SeqAIJ, 1785cb5b572fSBarry Smith MatSolve_SeqAIJ, 1786cb5b572fSBarry Smith MatSolveAdd_SeqAIJ, 1787cb5b572fSBarry Smith MatSolveTrans_SeqAIJ, 1788cb5b572fSBarry Smith MatSolveTransAdd_SeqAIJ, 1789cb5b572fSBarry Smith MatLUFactor_SeqAIJ, 1790cb5b572fSBarry Smith 0, 179117ab2063SBarry Smith MatRelax_SeqAIJ, 179217ab2063SBarry Smith MatTranspose_SeqAIJ, 1793cb5b572fSBarry Smith MatGetInfo_SeqAIJ, 1794cb5b572fSBarry Smith MatEqual_SeqAIJ, 1795cb5b572fSBarry Smith MatGetDiagonal_SeqAIJ, 1796cb5b572fSBarry Smith MatDiagonalScale_SeqAIJ, 1797cb5b572fSBarry Smith MatNorm_SeqAIJ, 1798cb5b572fSBarry Smith 0, 1799cb5b572fSBarry Smith MatAssemblyEnd_SeqAIJ, 180017ab2063SBarry Smith MatCompress_SeqAIJ, 1801cb5b572fSBarry Smith MatSetOption_SeqAIJ, 1802cb5b572fSBarry Smith MatZeroEntries_SeqAIJ, 1803cb5b572fSBarry Smith MatZeroRows_SeqAIJ, 1804cb5b572fSBarry Smith MatLUFactorSymbolic_SeqAIJ, 1805cb5b572fSBarry Smith MatLUFactorNumeric_SeqAIJ, 1806cb5b572fSBarry Smith 0, 1807cb5b572fSBarry Smith 0, 1808cb5b572fSBarry Smith MatGetSize_SeqAIJ, 1809cb5b572fSBarry Smith MatGetSize_SeqAIJ, 1810cb5b572fSBarry Smith MatGetOwnershipRange_SeqAIJ, 1811cb5b572fSBarry Smith MatILUFactorSymbolic_SeqAIJ, 1812cb5b572fSBarry Smith 0, 1813cb5b572fSBarry Smith 0, 1814cb5b572fSBarry Smith 0, 1815be6bf707SBarry Smith MatDuplicate_SeqAIJ, 1816cb5b572fSBarry Smith 0, 1817cb5b572fSBarry Smith 0, 1818cb5b572fSBarry Smith MatILUFactor_SeqAIJ, 1819cb5b572fSBarry Smith 0, 1820cb5b572fSBarry Smith 0, 1821cb5b572fSBarry Smith MatGetSubMatrices_SeqAIJ, 1822cb5b572fSBarry Smith MatIncreaseOverlap_SeqAIJ, 1823cb5b572fSBarry Smith MatGetValues_SeqAIJ, 1824cb5b572fSBarry Smith MatCopy_SeqAIJ, 1825f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1826cb5b572fSBarry Smith MatScale_SeqAIJ, 1827cb5b572fSBarry Smith 0, 1828cb5b572fSBarry Smith 0, 18296945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 18306945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 18313b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 18323b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 18333b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1834a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1835a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 18360513a670SBarry Smith MatColoringPatch_SeqAIJ, 18370513a670SBarry Smith 0, 1838cda55fadSBarry Smith MatPermute_SeqAIJ, 1839cda55fadSBarry Smith 0, 1840cda55fadSBarry Smith 0, 1841cda55fadSBarry Smith 0, 1842cda55fadSBarry Smith 0, 1843cda55fadSBarry Smith MatGetMaps_Petsc}; 184417ab2063SBarry Smith 184517ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 184617ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 184717ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 184817ab2063SBarry Smith 1849fb2e594dSBarry Smith EXTERN_C_BEGIN 18505615d1e5SSatish Balay #undef __FUNC__ 1851bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices_SeqAIJ" 1852bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices_SeqAIJ(Mat mat,int *indices) 1853bef8e0ddSBarry Smith { 1854bef8e0ddSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1855bef8e0ddSBarry Smith int i,nz,n; 1856bef8e0ddSBarry Smith 1857bef8e0ddSBarry Smith PetscFunctionBegin; 1858bef8e0ddSBarry Smith if (aij->indexshift) SETERRQ(1,1,"No support with 1 based indexing"); 1859bef8e0ddSBarry Smith 1860bef8e0ddSBarry Smith nz = aij->maxnz; 1861bef8e0ddSBarry Smith n = aij->n; 1862bef8e0ddSBarry Smith for (i=0; i<nz; i++) { 1863bef8e0ddSBarry Smith aij->j[i] = indices[i]; 1864bef8e0ddSBarry Smith } 1865bef8e0ddSBarry Smith aij->nz = nz; 1866bef8e0ddSBarry Smith for ( i=0; i<n; i++ ) { 1867bef8e0ddSBarry Smith aij->ilen[i] = aij->imax[i]; 1868bef8e0ddSBarry Smith } 1869bef8e0ddSBarry Smith 1870bef8e0ddSBarry Smith PetscFunctionReturn(0); 1871bef8e0ddSBarry Smith } 1872fb2e594dSBarry Smith EXTERN_C_END 1873bef8e0ddSBarry Smith 1874bef8e0ddSBarry Smith #undef __FUNC__ 1875bef8e0ddSBarry Smith #define __FUNC__ "MatSeqAIJSetColumnIndices" 1876bef8e0ddSBarry Smith /*@ 1877bef8e0ddSBarry Smith MatSeqAIJSetColumnIndices - Set the column indices for all the rows 1878bef8e0ddSBarry Smith in the matrix. 1879bef8e0ddSBarry Smith 1880bef8e0ddSBarry Smith Input Parameters: 1881bef8e0ddSBarry Smith + mat - the SeqAIJ matrix 1882bef8e0ddSBarry Smith - indices - the column indices 1883bef8e0ddSBarry Smith 1884bef8e0ddSBarry Smith Notes: 1885bef8e0ddSBarry Smith This can be called if you have precomputed the nonzero structure of the 1886bef8e0ddSBarry Smith matrix and want to provide it to the matrix object to improve the performance 1887bef8e0ddSBarry Smith of the MatSetValues() operation. 1888bef8e0ddSBarry Smith 1889bef8e0ddSBarry Smith You MUST have set the correct numbers of nonzeros per row in the call to 1890bef8e0ddSBarry Smith MatCreateSeqAIJ(). 1891bef8e0ddSBarry Smith 1892bef8e0ddSBarry Smith MUST be called before any calls to MatSetValues(); 1893bef8e0ddSBarry Smith 1894bef8e0ddSBarry Smith @*/ 1895bef8e0ddSBarry Smith int MatSeqAIJSetColumnIndices(Mat mat,int *indices) 1896bef8e0ddSBarry Smith { 1897bef8e0ddSBarry Smith int ierr,(*f)(Mat,int *); 1898bef8e0ddSBarry Smith 1899bef8e0ddSBarry Smith PetscFunctionBegin; 1900bef8e0ddSBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1901bef8e0ddSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void **)&f); 1902bef8e0ddSBarry Smith CHKERRQ(ierr); 1903bef8e0ddSBarry Smith if (f) { 1904bef8e0ddSBarry Smith ierr = (*f)(mat,indices);CHKERRQ(ierr); 1905bef8e0ddSBarry Smith } else { 1906bef8e0ddSBarry Smith SETERRQ(1,1,"Wrong type of matrix to set column indices"); 1907bef8e0ddSBarry Smith } 1908bef8e0ddSBarry Smith PetscFunctionReturn(0); 1909bef8e0ddSBarry Smith } 1910bef8e0ddSBarry Smith 1911be6bf707SBarry Smith /* ----------------------------------------------------------------------------------------*/ 1912be6bf707SBarry Smith 1913fb2e594dSBarry Smith EXTERN_C_BEGIN 1914be6bf707SBarry Smith #undef __FUNC__ 1915be6bf707SBarry Smith #define __FUNC__ "MatStoreValues_SeqAIJ" 1916be6bf707SBarry Smith int MatStoreValues_SeqAIJ(Mat mat) 1917be6bf707SBarry Smith { 1918be6bf707SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 1919be6bf707SBarry Smith int nz = aij->i[aij->m]+aij->indexshift; 1920be6bf707SBarry Smith 1921be6bf707SBarry Smith PetscFunctionBegin; 1922be6bf707SBarry Smith if (aij->nonew != 1) { 1923be6bf707SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 1924be6bf707SBarry Smith } 1925be6bf707SBarry Smith 1926be6bf707SBarry Smith /* allocate space for values if not already there */ 1927be6bf707SBarry Smith if (!aij->saved_values) { 1928be6bf707SBarry Smith aij->saved_values = (Scalar *) PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values); 1929be6bf707SBarry Smith } 1930be6bf707SBarry Smith 1931be6bf707SBarry Smith /* copy values over */ 1932be6bf707SBarry Smith PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar)); 1933be6bf707SBarry Smith PetscFunctionReturn(0); 1934be6bf707SBarry Smith } 1935fb2e594dSBarry Smith EXTERN_C_END 1936be6bf707SBarry Smith 1937be6bf707SBarry Smith #undef __FUNC__ 1938be6bf707SBarry Smith #define __FUNC__ "MatStoreValues" 1939be6bf707SBarry Smith /*@ 1940be6bf707SBarry Smith MatStoreValues - Stashes a copy of the matrix values; this allows, for 1941be6bf707SBarry Smith example, reuse of the linear part of a Jacobian, while recomputing the 1942be6bf707SBarry Smith nonlinear portion. 1943be6bf707SBarry Smith 1944be6bf707SBarry Smith Collect on Mat 1945be6bf707SBarry Smith 1946be6bf707SBarry Smith Input Parameters: 1947be6bf707SBarry Smith . mat - the matrix (currently on AIJ matrices support this option) 1948be6bf707SBarry Smith 1949be6bf707SBarry Smith Common Usage, with SNESSolve(): 1950be6bf707SBarry Smith $ Create Jacobian matrix 1951be6bf707SBarry Smith $ Set linear terms into matrix 1952be6bf707SBarry Smith $ Apply boundary conditions to matrix, at this time matrix must have 1953be6bf707SBarry Smith $ final nonzero structure (i.e. setting the nonlinear terms and applying 1954be6bf707SBarry Smith $ boundary conditions again will not change the nonzero structure 1955be6bf707SBarry Smith $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); 1956be6bf707SBarry Smith $ ierr = MatStoreValues(mat); 1957be6bf707SBarry Smith $ Call SNESSetJacobian() with matrix 1958be6bf707SBarry Smith $ In your Jacobian routine 1959be6bf707SBarry Smith $ ierr = MatRetrieveValues(mat); 1960be6bf707SBarry Smith $ Set nonlinear terms in matrix 1961be6bf707SBarry Smith 1962be6bf707SBarry Smith Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself: 1963be6bf707SBarry Smith $ // build linear portion of Jacobian 1964be6bf707SBarry Smith $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); 1965be6bf707SBarry Smith $ ierr = MatStoreValues(mat); 1966be6bf707SBarry Smith $ loop over nonlinear iterations 1967be6bf707SBarry Smith $ ierr = MatRetrieveValues(mat); 1968be6bf707SBarry Smith $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian 1969be6bf707SBarry Smith $ // call MatAssemblyBegin/End() on matrix 1970be6bf707SBarry Smith $ Solve linear system with Jacobian 1971be6bf707SBarry Smith $ endloop 1972be6bf707SBarry Smith 1973be6bf707SBarry Smith Notes: 1974be6bf707SBarry Smith Matrix must already be assemblied before calling this routine 1975be6bf707SBarry Smith Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before 1976be6bf707SBarry Smith calling this routine. 1977be6bf707SBarry Smith 1978be6bf707SBarry Smith .seealso: MatRetrieveValues() 1979be6bf707SBarry Smith 1980be6bf707SBarry Smith @*/ 1981be6bf707SBarry Smith int MatStoreValues(Mat mat) 1982be6bf707SBarry Smith { 1983be6bf707SBarry Smith int ierr,(*f)(Mat); 1984be6bf707SBarry Smith 1985be6bf707SBarry Smith PetscFunctionBegin; 1986be6bf707SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 1987be6bf707SBarry Smith if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1988be6bf707SBarry Smith if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1989be6bf707SBarry Smith 1990be6bf707SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void **)&f); 1991be6bf707SBarry Smith CHKERRQ(ierr); 1992be6bf707SBarry Smith if (f) { 1993be6bf707SBarry Smith ierr = (*f)(mat);CHKERRQ(ierr); 1994be6bf707SBarry Smith } else { 1995be6bf707SBarry Smith SETERRQ(1,1,"Wrong type of matrix to store values"); 1996be6bf707SBarry Smith } 1997be6bf707SBarry Smith PetscFunctionReturn(0); 1998be6bf707SBarry Smith } 1999be6bf707SBarry Smith 2000fb2e594dSBarry Smith EXTERN_C_BEGIN 2001be6bf707SBarry Smith #undef __FUNC__ 2002be6bf707SBarry Smith #define __FUNC__ "MatRetrieveValues_SeqAIJ" 2003be6bf707SBarry Smith int MatRetrieveValues_SeqAIJ(Mat mat) 2004be6bf707SBarry Smith { 2005be6bf707SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; 2006be6bf707SBarry Smith int nz = aij->i[aij->m]+aij->indexshift; 2007be6bf707SBarry Smith 2008be6bf707SBarry Smith PetscFunctionBegin; 2009be6bf707SBarry Smith if (aij->nonew != 1) { 2010be6bf707SBarry Smith SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); 2011be6bf707SBarry Smith } 2012be6bf707SBarry Smith if (!aij->saved_values) { 2013be6bf707SBarry Smith SETERRQ(1,1,"Must call MatStoreValues(A);first"); 2014be6bf707SBarry Smith } 2015be6bf707SBarry Smith 2016be6bf707SBarry Smith /* copy values over */ 2017be6bf707SBarry Smith PetscMemcpy(aij->a, aij->saved_values,nz*sizeof(Scalar)); 2018be6bf707SBarry Smith PetscFunctionReturn(0); 2019be6bf707SBarry Smith } 2020fb2e594dSBarry Smith EXTERN_C_END 2021be6bf707SBarry Smith 2022be6bf707SBarry Smith #undef __FUNC__ 2023be6bf707SBarry Smith #define __FUNC__ "MatRetrieveValues" 2024be6bf707SBarry Smith /*@ 2025be6bf707SBarry Smith MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for 2026be6bf707SBarry Smith example, reuse of the linear part of a Jacobian, while recomputing the 2027be6bf707SBarry Smith nonlinear portion. 2028be6bf707SBarry Smith 2029be6bf707SBarry Smith Collect on Mat 2030be6bf707SBarry Smith 2031be6bf707SBarry Smith Input Parameters: 2032be6bf707SBarry Smith . mat - the matrix (currently on AIJ matrices support this option) 2033be6bf707SBarry Smith 2034be6bf707SBarry Smith .seealso: MatStoreValues() 2035be6bf707SBarry Smith 2036be6bf707SBarry Smith @*/ 2037be6bf707SBarry Smith int MatRetrieveValues(Mat mat) 2038be6bf707SBarry Smith { 2039be6bf707SBarry Smith int ierr,(*f)(Mat); 2040be6bf707SBarry Smith 2041be6bf707SBarry Smith PetscFunctionBegin; 2042be6bf707SBarry Smith PetscValidHeaderSpecific(mat,MAT_COOKIE); 2043be6bf707SBarry Smith if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2044be6bf707SBarry Smith if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2045be6bf707SBarry Smith 2046be6bf707SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void **)&f); 2047be6bf707SBarry Smith CHKERRQ(ierr); 2048be6bf707SBarry Smith if (f) { 2049be6bf707SBarry Smith ierr = (*f)(mat);CHKERRQ(ierr); 2050be6bf707SBarry Smith } else { 2051be6bf707SBarry Smith SETERRQ(1,1,"Wrong type of matrix to retrieve values"); 2052be6bf707SBarry Smith } 2053be6bf707SBarry Smith PetscFunctionReturn(0); 2054be6bf707SBarry Smith } 2055be6bf707SBarry Smith 2056be6bf707SBarry Smith /* --------------------------------------------------------------------------------*/ 2057cb5b572fSBarry Smith 2058bef8e0ddSBarry Smith #undef __FUNC__ 20595615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ" 206017ab2063SBarry Smith /*@C 2061682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 20620d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 20636e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 20642bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 20652bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 206617ab2063SBarry Smith 2067db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2068db81eaa0SLois Curfman McInnes 206917ab2063SBarry Smith Input Parameters: 2070db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 207117ab2063SBarry Smith . m - number of rows 207217ab2063SBarry Smith . n - number of columns 207317ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 2074db81eaa0SLois Curfman McInnes - nzz - array containing the number of nonzeros in the various rows 20752bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 207617ab2063SBarry Smith 207717ab2063SBarry Smith Output Parameter: 2078416022c9SBarry Smith . A - the matrix 207917ab2063SBarry Smith 2080b259b22eSLois Curfman McInnes Notes: 208117ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 208217ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 20830002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 208444cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 208517ab2063SBarry Smith 208617ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 2087a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 20883d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 20896da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 209017ab2063SBarry Smith 2091682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 20924fca80b9SLois Curfman McInnes improve numerical efficiency of matrix-vector products and solves. We 2093682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 20946c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 20956c7ebb05SLois Curfman McInnes 20966c7ebb05SLois Curfman McInnes Options Database Keys: 2097db81eaa0SLois Curfman McInnes + -mat_aij_no_inode - Do not use inodes 2098db81eaa0SLois Curfman McInnes . -mat_aij_inode_limit <limit> - Sets inode limit (max limit=5) 2099db81eaa0SLois Curfman McInnes - -mat_aij_oneindex - Internally use indexing starting at 1 2100db81eaa0SLois Curfman McInnes rather than 0. Note that when calling MatSetValues(), 2101db81eaa0SLois Curfman McInnes the user still MUST index entries starting at 0! 210217ab2063SBarry Smith 2103027ccd11SLois Curfman McInnes Level: intermediate 2104027ccd11SLois Curfman McInnes 2105bef8e0ddSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices() 210617ab2063SBarry Smith @*/ 2107416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 210817ab2063SBarry Smith { 2109416022c9SBarry Smith Mat B; 2110416022c9SBarry Smith Mat_SeqAIJ *b; 21116945ee14SBarry Smith int i, len, ierr, flg,size; 21126945ee14SBarry Smith 21133a40ed3dSBarry Smith PetscFunctionBegin; 21146945ee14SBarry Smith MPI_Comm_size(comm,&size); 2115a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Comm must be of size 1"); 2116d5d45c9bSBarry Smith 2117416022c9SBarry Smith *A = 0; 21183f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,"Mat",comm,MatDestroy,MatView); 2119416022c9SBarry Smith PLogObjectCreate(B); 21200452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 212144cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 21220a6ffc59SBarry Smith PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps)); 2123e1311b90SBarry Smith B->ops->destroy = MatDestroy_SeqAIJ; 2124e1311b90SBarry Smith B->ops->view = MatView_SeqAIJ; 2125416022c9SBarry Smith B->factor = 0; 2126416022c9SBarry Smith B->lupivotthreshold = 1.0; 212790f02eecSBarry Smith B->mapping = 0; 2128e1311b90SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold,&flg);CHKERRQ(ierr); 21297a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 2130e1311b90SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums",(int*)&b->ilu_preserve_row_sums);CHKERRQ(ierr); 2131416022c9SBarry Smith b->row = 0; 2132416022c9SBarry Smith b->col = 0; 213382bf6240SBarry Smith b->icol = 0; 2134416022c9SBarry Smith b->indexshift = 0; 2135b810aeb4SBarry Smith b->reallocs = 0; 213669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 213769957df2SSatish Balay if (flg) b->indexshift = -1; 213817ab2063SBarry Smith 213944cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 214044cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 2141a5ae1ecdSBarry Smith 21424d197716SBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 21434d197716SBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 2144a5ae1ecdSBarry Smith 21450452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 2146b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 21477b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 21487b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 2149416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 215017ab2063SBarry Smith nz = nz*m; 21513a40ed3dSBarry Smith } else { 215217ab2063SBarry Smith nz = 0; 2153416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 215417ab2063SBarry Smith } 215517ab2063SBarry Smith 215617ab2063SBarry Smith /* allocate the matrix space */ 2157416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 21580452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 2159416022c9SBarry Smith b->j = (int *) (b->a + nz); 2160cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 2161416022c9SBarry Smith b->i = b->j + nz; 2162416022c9SBarry Smith b->singlemalloc = 1; 216317ab2063SBarry Smith 2164416022c9SBarry Smith b->i[0] = -b->indexshift; 216517ab2063SBarry Smith for (i=1; i<m+1; i++) { 2166416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 216717ab2063SBarry Smith } 216817ab2063SBarry Smith 2169416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 21700452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 2171f09e8eb9SSatish Balay PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 2172416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 217317ab2063SBarry Smith 2174416022c9SBarry Smith b->nz = 0; 2175416022c9SBarry Smith b->maxnz = nz; 2176416022c9SBarry Smith b->sorted = 0; 2177416022c9SBarry Smith b->roworiented = 1; 2178416022c9SBarry Smith b->nonew = 0; 2179416022c9SBarry Smith b->diag = 0; 2180416022c9SBarry Smith b->solve_work = 0; 218171bd300dSLois Curfman McInnes b->spptr = 0; 2182754ec7b1SSatish Balay b->inode.node_count = 0; 2183754ec7b1SSatish Balay b->inode.size = 0; 21846c7ebb05SLois Curfman McInnes b->inode.limit = 5; 21856c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 2186be6bf707SBarry Smith b->saved_values = 0; 21874e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 218817ab2063SBarry Smith 2189416022c9SBarry Smith *A = B; 21904e220ebcSLois Curfman McInnes 21914b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 21924b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 219369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 219469957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 21954b14c69eSBarry Smith #endif 219669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 219769957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 219869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 219969957df2SSatish Balay if (flg) { 2200a8c6a408SBarry Smith if (!b->indexshift) SETERRQ( PETSC_ERR_LIB,0,"need -mat_aij_oneindex with -mat_aij_dxml"); 2201416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 220217ab2063SBarry Smith } 220369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 220469957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 2205bef8e0ddSBarry Smith 2206bef8e0ddSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqAIJSetColumnIndices_C", 2207bef8e0ddSBarry Smith "MatSeqAIJSetColumnIndices_SeqAIJ", 2208bef8e0ddSBarry Smith (void*)MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); 2209be6bf707SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C", 2210be6bf707SBarry Smith "MatStoreValues_SeqAIJ", 2211be6bf707SBarry Smith (void*)MatStoreValues_SeqAIJ);CHKERRQ(ierr); 2212be6bf707SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C", 2213be6bf707SBarry Smith "MatRetrieveValues_SeqAIJ", 2214be6bf707SBarry Smith (void*)MatRetrieveValues_SeqAIJ);CHKERRQ(ierr); 22153a40ed3dSBarry Smith PetscFunctionReturn(0); 221617ab2063SBarry Smith } 221717ab2063SBarry Smith 22185615d1e5SSatish Balay #undef __FUNC__ 2219be6bf707SBarry Smith #define __FUNC__ "MatDuplicate_SeqAIJ" 2220be6bf707SBarry Smith int MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) 222117ab2063SBarry Smith { 2222416022c9SBarry Smith Mat C; 2223416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 2224bef8e0ddSBarry Smith int i,len, m = a->m,shift = a->indexshift,ierr; 222517ab2063SBarry Smith 22263a40ed3dSBarry Smith PetscFunctionBegin; 22274043dd9cSLois Curfman McInnes *B = 0; 22283f1db9ecSBarry Smith PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQAIJ,"Mat",A->comm,MatDestroy,MatView); 2229416022c9SBarry Smith PLogObjectCreate(C); 22300452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 2231f830108cSBarry Smith PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps)); 2232e1311b90SBarry Smith C->ops->destroy = MatDestroy_SeqAIJ; 2233e1311b90SBarry Smith C->ops->view = MatView_SeqAIJ; 2234416022c9SBarry Smith C->factor = A->factor; 2235416022c9SBarry Smith c->row = 0; 2236416022c9SBarry Smith c->col = 0; 223782bf6240SBarry Smith c->icol = 0; 2238416022c9SBarry Smith c->indexshift = shift; 2239c456f294SBarry Smith C->assembled = PETSC_TRUE; 224017ab2063SBarry Smith 224144cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 224244cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 224344cd7ae7SLois Curfman McInnes C->M = a->m; 224444cd7ae7SLois Curfman McInnes C->N = a->n; 224517ab2063SBarry Smith 22460452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 22470452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 224817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 2249416022c9SBarry Smith c->imax[i] = a->imax[i]; 2250416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 225117ab2063SBarry Smith } 225217ab2063SBarry Smith 225317ab2063SBarry Smith /* allocate the matrix space */ 2254416022c9SBarry Smith c->singlemalloc = 1; 2255416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 22560452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 2257416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 2258416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 2259416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 226017ab2063SBarry Smith if (m > 0) { 2261416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 2262be6bf707SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 2263416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 2264be6bf707SBarry Smith } else { 2265be6bf707SBarry Smith PetscMemzero(c->a,(a->i[m]+shift)*sizeof(Scalar)); 226617ab2063SBarry Smith } 226708480c60SBarry Smith } 226817ab2063SBarry Smith 2269f09e8eb9SSatish Balay PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqAIJ)); 2270416022c9SBarry Smith c->sorted = a->sorted; 2271416022c9SBarry Smith c->roworiented = a->roworiented; 2272416022c9SBarry Smith c->nonew = a->nonew; 22737a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 2274be6bf707SBarry Smith c->saved_values = 0; 227517ab2063SBarry Smith 2276416022c9SBarry Smith if (a->diag) { 22770452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 2278416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 227917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 2280416022c9SBarry Smith c->diag[i] = a->diag[i]; 228117ab2063SBarry Smith } 22823a40ed3dSBarry Smith } else c->diag = 0; 22836c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 22846c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 2285754ec7b1SSatish Balay if (a->inode.size){ 2286daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 2287754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 2288daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 2289754ec7b1SSatish Balay } else { 2290754ec7b1SSatish Balay c->inode.size = 0; 2291754ec7b1SSatish Balay c->inode.node_count = 0; 2292754ec7b1SSatish Balay } 2293416022c9SBarry Smith c->nz = a->nz; 2294416022c9SBarry Smith c->maxnz = a->maxnz; 2295416022c9SBarry Smith c->solve_work = 0; 229676dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 2297754ec7b1SSatish Balay 2298416022c9SBarry Smith *B = C; 22997b2a1423SBarry Smith ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); 23003a40ed3dSBarry Smith PetscFunctionReturn(0); 230117ab2063SBarry Smith } 230217ab2063SBarry Smith 23035615d1e5SSatish Balay #undef __FUNC__ 23045615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ" 230519bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 230617ab2063SBarry Smith { 2307416022c9SBarry Smith Mat_SeqAIJ *a; 2308416022c9SBarry Smith Mat B; 230917699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 2310bcd2baecSBarry Smith MPI_Comm comm; 231117ab2063SBarry Smith 23123a40ed3dSBarry Smith PetscFunctionBegin; 231319bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 231417699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 2315a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,0,"view must have one processor"); 231690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 23170752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 2318a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not matrix object in file"); 231917ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 232017ab2063SBarry Smith 2321d64ed03dSBarry Smith if (nz < 0) { 2322a8c6a408SBarry Smith SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format on disk, cannot load as SeqAIJ"); 2323d64ed03dSBarry Smith } 2324d64ed03dSBarry Smith 232517ab2063SBarry Smith /* read in row lengths */ 23260452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 23270752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 232817ab2063SBarry Smith 232917ab2063SBarry Smith /* create our matrix */ 2330416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 2331416022c9SBarry Smith B = *A; 2332416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 2333416022c9SBarry Smith shift = a->indexshift; 233417ab2063SBarry Smith 233517ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 23360752156aSBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT); CHKERRQ(ierr); 233717ab2063SBarry Smith if (shift) { 233817ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 2339416022c9SBarry Smith a->j[i] += 1; 234017ab2063SBarry Smith } 234117ab2063SBarry Smith } 234217ab2063SBarry Smith 234317ab2063SBarry Smith /* read in nonzero values */ 23440752156aSBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR); CHKERRQ(ierr); 234517ab2063SBarry Smith 234617ab2063SBarry Smith /* set matrix "i" values */ 2347416022c9SBarry Smith a->i[0] = -shift; 234817ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 2349416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 2350416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 235117ab2063SBarry Smith } 23520452661fSBarry Smith PetscFree(rowlengths); 235317ab2063SBarry Smith 23546d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 23556d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 23563a40ed3dSBarry Smith PetscFunctionReturn(0); 235717ab2063SBarry Smith } 235817ab2063SBarry Smith 23595615d1e5SSatish Balay #undef __FUNC__ 23605615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ" 23618f6be9afSLois Curfman McInnes int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 23627264ac53SSatish Balay { 23637264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 23647264ac53SSatish Balay 23653a40ed3dSBarry Smith PetscFunctionBegin; 2366a8c6a408SBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 23677264ac53SSatish Balay 23687264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 23697264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 2370bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 23713a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2372bcd2baecSBarry Smith } 23737264ac53SSatish Balay 23747264ac53SSatish Balay /* if the a->i are the same */ 23758108c231SLois Curfman McInnes if (PetscMemcmp(a->i,b->i,(a->m+1)*sizeof(int))) { 23763a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 23777264ac53SSatish Balay } 23787264ac53SSatish Balay 23797264ac53SSatish Balay /* if a->j are the same */ 2380bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 23813a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 2382bcd2baecSBarry Smith } 2383bcd2baecSBarry Smith 2384bcd2baecSBarry Smith /* if a->a are the same */ 238519bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 23863a40ed3dSBarry Smith *flg = PETSC_FALSE; PetscFunctionReturn(0); 23877264ac53SSatish Balay } 238877c4ece6SBarry Smith *flg = PETSC_TRUE; 23893a40ed3dSBarry Smith PetscFunctionReturn(0); 23907264ac53SSatish Balay 23917264ac53SSatish Balay } 2392