16d84be18SBarry Smith 26945ee14SBarry Smith 317ab2063SBarry Smith #ifndef lint 4*0513a670SBarry Smith static char vcid[] = "$Id: aij.c,v 1.203 1997/01/06 20:24:23 balay Exp bsmith $"; 517ab2063SBarry Smith #endif 617ab2063SBarry Smith 7d5d45c9bSBarry Smith /* 85a838052SSatish Balay B Defines the basic matrix operations for the AIJ (compressed row) 9d5d45c9bSBarry Smith matrix storage format. 10d5d45c9bSBarry Smith */ 1170f55243SBarry Smith #include "src/mat/impls/aij/seq/aij.h" 12f5eb4b81SSatish Balay #include "src/vec/vecimpl.h" 13f5eb4b81SSatish Balay #include "src/inline/spops.h" 14e4d965acSSatish Balay #include "petsc.h" 15f5eb4b81SSatish Balay #include "src/inline/bitarray.h" 1617ab2063SBarry Smith 17a2ce50c7SBarry Smith /* 18a2ce50c7SBarry Smith Basic AIJ format ILU based on drop tolerance 19a2ce50c7SBarry Smith */ 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 int ierr = 1; 26a2ce50c7SBarry Smith 27e3372554SBarry Smith SETERRQ(ierr,0,"Not implemented"); 28a2ce50c7SBarry Smith } 29a2ce50c7SBarry Smith 30bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); 3117ab2063SBarry Smith 325615d1e5SSatish Balay #undef __FUNC__ 335615d1e5SSatish Balay #define __FUNC__ "MatGetRowIJ_SeqAIJ" 3431625ec5SSatish Balay static int MatGetRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja, 356945ee14SBarry Smith PetscTruth *done) 3617ab2063SBarry Smith { 37416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 386945ee14SBarry Smith int ierr,i,ishift; 3917ab2063SBarry Smith 4031625ec5SSatish Balay *m = A->m; 416945ee14SBarry Smith if (!ia) return 0; 426945ee14SBarry Smith ishift = a->indexshift; 436945ee14SBarry Smith if (symmetric) { 4431625ec5SSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 456945ee14SBarry Smith } else if (oshift == 0 && ishift == -1) { 4631625ec5SSatish Balay int nz = a->i[a->m]; 473b2fbd54SBarry Smith /* malloc space and subtract 1 from i and j indices */ 4831625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 493b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 503b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] - 1; 5131625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] - 1; 526945ee14SBarry Smith } else if (oshift == 1 && ishift == 0) { 5331625ec5SSatish Balay int nz = a->i[a->m] + 1; 543b2fbd54SBarry Smith /* malloc space and add 1 to i and j indices */ 5531625ec5SSatish Balay *ia = (int *) PetscMalloc( (a->m+1)*sizeof(int) ); CHKPTRQ(*ia); 563b2fbd54SBarry Smith *ja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(*ja); 573b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) (*ja)[i] = a->j[i] + 1; 5831625ec5SSatish Balay for ( i=0; i<a->m+1; i++ ) (*ia)[i] = a->i[i] + 1; 596945ee14SBarry Smith } else { 606945ee14SBarry Smith *ia = a->i; *ja = a->j; 61a2ce50c7SBarry Smith } 62a2ce50c7SBarry Smith 63a2744918SBarry Smith return 0; 64a2744918SBarry Smith } 65a2744918SBarry Smith 665615d1e5SSatish Balay #undef __FUNC__ 675615d1e5SSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqAIJ" 683b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia,int **ja, 696945ee14SBarry Smith PetscTruth *done) 706945ee14SBarry Smith { 716945ee14SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 723b2fbd54SBarry Smith int ishift = a->indexshift; 736945ee14SBarry Smith 746945ee14SBarry Smith if (!ia) return 0; 753b2fbd54SBarry Smith if (symmetric || (oshift == 0 && ishift == -1) || (oshift == 1 && ishift == 0)) { 766945ee14SBarry Smith PetscFree(*ia); 776945ee14SBarry Smith PetscFree(*ja); 78bcd2baecSBarry Smith } 7917ab2063SBarry Smith return 0; 8017ab2063SBarry Smith } 8117ab2063SBarry Smith 825615d1e5SSatish Balay #undef __FUNC__ 835615d1e5SSatish Balay #define __FUNC__ "MatGetColumnIJ_SeqAIJ" 8443a90d84SBarry Smith int MatGetColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja, 853b2fbd54SBarry Smith PetscTruth *done) 863b2fbd54SBarry Smith { 873b2fbd54SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 88a93ec695SBarry Smith int ierr,i,ishift = a->indexshift,*collengths,*cia,*cja,n = A->n,m = A->m; 89a93ec695SBarry Smith int nz = a->i[m]+ishift,row,*jj,mr,col; 903b2fbd54SBarry Smith 913b2fbd54SBarry Smith *nn = A->n; 923b2fbd54SBarry Smith if (!ia) return 0; 933b2fbd54SBarry Smith if (symmetric) { 94179192dfSSatish Balay ierr = MatToSymmetricIJ_SeqAIJ(a->m,a->i,a->j,ishift,oshift,ia,ja); CHKERRQ(ierr); 953b2fbd54SBarry Smith } else { 9661d2ded1SBarry Smith collengths = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(collengths); 973b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 983b2fbd54SBarry Smith cia = (int *) PetscMalloc( (n+1)*sizeof(int) ); CHKPTRQ(cia); 99a93ec695SBarry Smith cja = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cja); 1003b2fbd54SBarry Smith jj = a->j; 1013b2fbd54SBarry Smith for ( i=0; i<nz; i++ ) { 1023b2fbd54SBarry Smith collengths[jj[i] + ishift]++; 1033b2fbd54SBarry Smith } 1043b2fbd54SBarry Smith cia[0] = oshift; 1053b2fbd54SBarry Smith for ( i=0; i<n; i++) { 1063b2fbd54SBarry Smith cia[i+1] = cia[i] + collengths[i]; 1073b2fbd54SBarry Smith } 1083b2fbd54SBarry Smith PetscMemzero(collengths,n*sizeof(int)); 1093b2fbd54SBarry Smith jj = a->j; 110a93ec695SBarry Smith for ( row=0; row<m; row++ ) { 111a93ec695SBarry Smith mr = a->i[row+1] - a->i[row]; 112a93ec695SBarry Smith for ( i=0; i<mr; i++ ) { 1133b2fbd54SBarry Smith col = *jj++ + ishift; 1143b2fbd54SBarry Smith cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1153b2fbd54SBarry Smith } 1163b2fbd54SBarry Smith } 1173b2fbd54SBarry Smith PetscFree(collengths); 1183b2fbd54SBarry Smith *ia = cia; *ja = cja; 1193b2fbd54SBarry Smith } 1203b2fbd54SBarry Smith 1213b2fbd54SBarry Smith return 0; 1223b2fbd54SBarry Smith } 1233b2fbd54SBarry Smith 1245615d1e5SSatish Balay #undef __FUNC__ 1255615d1e5SSatish Balay #define __FUNC__ "MatRestoreColumnIJ_SeqAIJ" 12643a90d84SBarry Smith int MatRestoreColumnIJ_SeqAIJ(Mat A,int oshift,PetscTruth symmetric,int *n,int **ia, 1273b2fbd54SBarry Smith int **ja,PetscTruth *done) 1283b2fbd54SBarry Smith { 1293b2fbd54SBarry Smith if (!ia) return 0; 1303b2fbd54SBarry Smith 1313b2fbd54SBarry Smith PetscFree(*ia); 1323b2fbd54SBarry Smith PetscFree(*ja); 1333b2fbd54SBarry Smith 1343b2fbd54SBarry Smith return 0; 1353b2fbd54SBarry Smith } 1363b2fbd54SBarry Smith 137227d817aSBarry Smith #define CHUNKSIZE 15 13817ab2063SBarry Smith 1395615d1e5SSatish Balay #undef __FUNC__ 1405615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqAIJ" 14161d2ded1SBarry Smith int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 14217ab2063SBarry Smith { 143416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 144416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 1454b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 146d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 147416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 14817ab2063SBarry Smith 14917ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 150416022c9SBarry Smith row = im[k]; 1513b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 152e3372554SBarry Smith if (row < 0) SETERRQ(1,0,"Negative row"); 153e3372554SBarry Smith if (row >= a->m) SETERRQ(1,0,"Row too large"); 1543b2fbd54SBarry Smith #endif 15517ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 15617ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 157416022c9SBarry Smith low = 0; 15817ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 1593b2fbd54SBarry Smith #if defined(PETSC_BOPT_g) 160e3372554SBarry Smith if (in[l] < 0) SETERRQ(1,0,"Negative column"); 161e3372554SBarry Smith if (in[l] >= a->n) SETERRQ(1,0,"Column too large"); 1623b2fbd54SBarry Smith #endif 1634b0e389bSBarry Smith col = in[l] - shift; 1644b0e389bSBarry Smith if (roworiented) { 1654b0e389bSBarry Smith value = *v++; 1664b0e389bSBarry Smith } 1674b0e389bSBarry Smith else { 1684b0e389bSBarry Smith value = v[k + l*m]; 1694b0e389bSBarry Smith } 170416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 171416022c9SBarry Smith while (high-low > 5) { 172416022c9SBarry Smith t = (low+high)/2; 173416022c9SBarry Smith if (rp[t] > col) high = t; 174416022c9SBarry Smith else low = t; 17517ab2063SBarry Smith } 176416022c9SBarry Smith for ( i=low; i<high; i++ ) { 17717ab2063SBarry Smith if (rp[i] > col) break; 17817ab2063SBarry Smith if (rp[i] == col) { 179416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 18017ab2063SBarry Smith else ap[i] = value; 18117ab2063SBarry Smith goto noinsert; 18217ab2063SBarry Smith } 18317ab2063SBarry Smith } 18417ab2063SBarry Smith if (nonew) goto noinsert; 18517ab2063SBarry Smith if (nrow >= rmax) { 18617ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 187416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 18817ab2063SBarry Smith Scalar *new_a; 18917ab2063SBarry Smith 19017ab2063SBarry Smith /* malloc new storage space */ 191416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1920452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 19317ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 19417ab2063SBarry Smith new_i = new_j + new_nz; 19517ab2063SBarry Smith 19617ab2063SBarry Smith /* copy over old data into new slots */ 19717ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 198416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 199416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 200416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 201416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 20217ab2063SBarry Smith len*sizeof(int)); 203416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 204416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 20517ab2063SBarry Smith len*sizeof(Scalar)); 20617ab2063SBarry Smith /* free up old matrix storage */ 2070452661fSBarry Smith PetscFree(a->a); 2080452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 209416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 210416022c9SBarry Smith a->singlemalloc = 1; 21117ab2063SBarry Smith 21217ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 213416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 214416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 215416022c9SBarry Smith a->maxnz += CHUNKSIZE; 216b810aeb4SBarry Smith a->reallocs++; 21717ab2063SBarry Smith } 218416022c9SBarry Smith N = nrow++ - 1; a->nz++; 219416022c9SBarry Smith /* shift up all the later entries in this row */ 220416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 22117ab2063SBarry Smith rp[ii+1] = rp[ii]; 22217ab2063SBarry Smith ap[ii+1] = ap[ii]; 22317ab2063SBarry Smith } 22417ab2063SBarry Smith rp[i] = col; 22517ab2063SBarry Smith ap[i] = value; 22617ab2063SBarry Smith noinsert:; 227416022c9SBarry Smith low = i + 1; 22817ab2063SBarry Smith } 22917ab2063SBarry Smith ailen[row] = nrow; 23017ab2063SBarry Smith } 23117ab2063SBarry Smith return 0; 23217ab2063SBarry Smith } 23317ab2063SBarry Smith 2345615d1e5SSatish Balay #undef __FUNC__ 2355615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqAIJ" 2367eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 2377eb43aa7SLois Curfman McInnes { 2387eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 239b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 2407eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 2417eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 2427eb43aa7SLois Curfman McInnes 2437eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 2447eb43aa7SLois Curfman McInnes row = im[k]; 245e3372554SBarry Smith if (row < 0) SETERRQ(1,0,"Negative row"); 246e3372554SBarry Smith if (row >= a->m) SETERRQ(1,0,"Row too large"); 2477eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 2487eb43aa7SLois Curfman McInnes nrow = ailen[row]; 2497eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 250e3372554SBarry Smith if (in[l] < 0) SETERRQ(1,0,"Negative column"); 251e3372554SBarry Smith if (in[l] >= a->n) SETERRQ(1,0,"Column too large"); 2527eb43aa7SLois Curfman McInnes col = in[l] - shift; 2537eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 2547eb43aa7SLois Curfman McInnes while (high-low > 5) { 2557eb43aa7SLois Curfman McInnes t = (low+high)/2; 2567eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 2577eb43aa7SLois Curfman McInnes else low = t; 2587eb43aa7SLois Curfman McInnes } 2597eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 2607eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 2617eb43aa7SLois Curfman McInnes if (rp[i] == col) { 262b49de8d1SLois Curfman McInnes *v++ = ap[i]; 2637eb43aa7SLois Curfman McInnes goto finished; 2647eb43aa7SLois Curfman McInnes } 2657eb43aa7SLois Curfman McInnes } 266b49de8d1SLois Curfman McInnes *v++ = zero; 2677eb43aa7SLois Curfman McInnes finished:; 2687eb43aa7SLois Curfman McInnes } 2697eb43aa7SLois Curfman McInnes } 2707eb43aa7SLois Curfman McInnes return 0; 2717eb43aa7SLois Curfman McInnes } 2727eb43aa7SLois Curfman McInnes 27317ab2063SBarry Smith #include "draw.h" 27417ab2063SBarry Smith #include "pinclude/pviewer.h" 27577c4ece6SBarry Smith #include "sys.h" 27617ab2063SBarry Smith 2775615d1e5SSatish Balay #undef __FUNC__ 2785615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Binary" 279416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 28017ab2063SBarry Smith { 281416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 282416022c9SBarry Smith int i, fd, *col_lens, ierr; 28317ab2063SBarry Smith 28490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 2850452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 286416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 287416022c9SBarry Smith col_lens[1] = a->m; 288416022c9SBarry Smith col_lens[2] = a->n; 289416022c9SBarry Smith col_lens[3] = a->nz; 290416022c9SBarry Smith 291416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 292416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 293416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 29417ab2063SBarry Smith } 29577c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr); 2960452661fSBarry Smith PetscFree(col_lens); 297416022c9SBarry Smith 298416022c9SBarry Smith /* store column indices (zero start index) */ 299416022c9SBarry Smith if (a->indexshift) { 300416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 30117ab2063SBarry Smith } 30277c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->j,a->nz,BINARY_INT,0); CHKERRQ(ierr); 303416022c9SBarry Smith if (a->indexshift) { 304416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 30517ab2063SBarry Smith } 306416022c9SBarry Smith 307416022c9SBarry Smith /* store nonzero values */ 30877c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,a->a,a->nz,BINARY_SCALAR,0); CHKERRQ(ierr); 30917ab2063SBarry Smith return 0; 31017ab2063SBarry Smith } 311416022c9SBarry Smith 3125615d1e5SSatish Balay #undef __FUNC__ 3135615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_ASCII" 314416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 315416022c9SBarry Smith { 316416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 317496e697eSBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift, format, flg1,flg2; 31817ab2063SBarry Smith FILE *fd; 31917ab2063SBarry Smith char *outputname; 32017ab2063SBarry Smith 32190ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 322416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 32390ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 324a93ec695SBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO) { 32595e01e2fSLois Curfman McInnes return 0; 32695e01e2fSLois Curfman McInnes } 327a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_INFO_LONG) { 328496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg1); CHKERRQ(ierr); 329496e697eSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg2); CHKERRQ(ierr); 330496e697eSBarry Smith if (flg1 || flg2) fprintf(fd," not using I-node routines\n"); 33195e01e2fSLois Curfman McInnes else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 33295e01e2fSLois Curfman McInnes a->inode.node_count,a->inode.limit); 33317ab2063SBarry Smith } 334a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_MATLAB) { 335416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 3364e220ebcSLois Curfman McInnes fprintf(fd,"%% Nonzeros = %d \n",a->nz); 3374e220ebcSLois Curfman McInnes fprintf(fd,"zzz = zeros(%d,3);\n",a->nz); 33817ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 33917ab2063SBarry Smith 34017ab2063SBarry Smith for (i=0; i<m; i++) { 341416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 34217ab2063SBarry Smith #if defined(PETSC_COMPLEX) 3436945ee14SBarry Smith fprintf(fd,"%d %d %18.16e + %18.16e i \n",i+1,a->j[j]+!shift,real(a->a[j]), 344416022c9SBarry Smith imag(a->a[j])); 34517ab2063SBarry Smith #else 3467a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 34717ab2063SBarry Smith #endif 34817ab2063SBarry Smith } 34917ab2063SBarry Smith } 35017ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 35117ab2063SBarry Smith } 352a93ec695SBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 35344cd7ae7SLois Curfman McInnes for ( i=0; i<m; i++ ) { 35444cd7ae7SLois Curfman McInnes fprintf(fd,"row %d:",i); 35544cd7ae7SLois Curfman McInnes for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 35644cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX) 35744cd7ae7SLois Curfman McInnes if (imag(a->a[j]) != 0.0 && real(a->a[j]) != 0.0) 35844cd7ae7SLois Curfman McInnes fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 35944cd7ae7SLois Curfman McInnes else if (real(a->a[j]) != 0.0) 36044cd7ae7SLois Curfman McInnes fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 36144cd7ae7SLois Curfman McInnes #else 36244cd7ae7SLois Curfman McInnes if (a->a[j] != 0.0) fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 36344cd7ae7SLois Curfman McInnes #endif 36444cd7ae7SLois Curfman McInnes } 36544cd7ae7SLois Curfman McInnes fprintf(fd,"\n"); 36644cd7ae7SLois Curfman McInnes } 36744cd7ae7SLois Curfman McInnes } 36817ab2063SBarry Smith else { 36917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 37017ab2063SBarry Smith fprintf(fd,"row %d:",i); 371416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 37217ab2063SBarry Smith #if defined(PETSC_COMPLEX) 373416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 374416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 37517ab2063SBarry Smith } 37617ab2063SBarry Smith else { 377416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 37817ab2063SBarry Smith } 37917ab2063SBarry Smith #else 380416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 38117ab2063SBarry Smith #endif 38217ab2063SBarry Smith } 38317ab2063SBarry Smith fprintf(fd,"\n"); 38417ab2063SBarry Smith } 38517ab2063SBarry Smith } 38617ab2063SBarry Smith fflush(fd); 387416022c9SBarry Smith return 0; 388416022c9SBarry Smith } 389416022c9SBarry Smith 3905615d1e5SSatish Balay #undef __FUNC__ 3915615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ_Draw" 392416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 393416022c9SBarry Smith { 394416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 395cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 396*0513a670SBarry Smith int format; 397*0513a670SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r,maxv; 398bcd2baecSBarry Smith Draw draw; 399cddf8d76SBarry Smith DrawButton button; 40019bcc07fSBarry Smith PetscTruth isnull; 401cddf8d76SBarry Smith 402*0513a670SBarry Smith ierr = ViewerDrawGetDraw(viewer,&draw); CHKERRQ(ierr); 403*0513a670SBarry Smith ierr = DrawClear(draw); CHKERRQ(ierr); 404*0513a670SBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 40519bcc07fSBarry Smith ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0; 40619bcc07fSBarry Smith 407416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 408416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 409416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 410416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 411*0513a670SBarry Smith 412*0513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 413*0513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 414cddf8d76SBarry Smith color = DRAW_BLUE; 415416022c9SBarry Smith for ( i=0; i<m; i++ ) { 416cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 417416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 418cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 419cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 420cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 421cddf8d76SBarry Smith #else 422cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 423cddf8d76SBarry Smith #endif 424cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 425cddf8d76SBarry Smith } 426cddf8d76SBarry Smith } 427cddf8d76SBarry Smith color = DRAW_CYAN; 428cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 429cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 430cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 431cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 432cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 433cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 434cddf8d76SBarry Smith } 435cddf8d76SBarry Smith } 436cddf8d76SBarry Smith color = DRAW_RED; 437cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 438cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 439cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 440cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 441cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 442cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 443cddf8d76SBarry Smith #else 444cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 445cddf8d76SBarry Smith #endif 446cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 447416022c9SBarry Smith } 448416022c9SBarry Smith } 449*0513a670SBarry Smith } else { 450*0513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 451*0513a670SBarry Smith /* first determine max of all nonzero values */ 452*0513a670SBarry Smith int nz = a->nz,count; 453*0513a670SBarry Smith Draw popup; 454*0513a670SBarry Smith 455*0513a670SBarry Smith maxv = 0.0; 456*0513a670SBarry Smith for ( i=0; i<nz; i++ ) { 457*0513a670SBarry Smith if (PetscAbsScalar(a->a[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); 458*0513a670SBarry Smith } 459*0513a670SBarry Smith ierr = DrawCreatePopUp(draw,&popup); CHKERRQ(ierr); 460*0513a670SBarry Smith ierr = DrawScalePopup(popup,0.0,maxv); CHKERRQ(ierr); 461*0513a670SBarry Smith count = 0; 462*0513a670SBarry Smith for ( i=0; i<m; i++ ) { 463*0513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 464*0513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 465*0513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 466*0513a670SBarry Smith color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv); 467*0513a670SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 468*0513a670SBarry Smith count++; 469*0513a670SBarry Smith } 470*0513a670SBarry Smith } 471*0513a670SBarry Smith } 472416022c9SBarry Smith DrawFlush(draw); 473cddf8d76SBarry Smith DrawGetPause(draw,&pause); 474cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 475cddf8d76SBarry Smith 476cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 4776945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 478cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 479cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 480cddf8d76SBarry Smith DrawClear(draw); 481cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 482cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 483cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 484cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 485cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 486cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 487cddf8d76SBarry Smith w *= scale; h *= scale; 488cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 489*0513a670SBarry Smith if (format != VIEWER_FORMAT_DRAW_CONTOUR) { 490*0513a670SBarry Smith /* Blue for negative, Cyan for zero and Red for positive */ 491cddf8d76SBarry Smith color = DRAW_BLUE; 492cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 493cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 494cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 495cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 496cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 497cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 498cddf8d76SBarry Smith #else 499cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 500cddf8d76SBarry Smith #endif 501cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 502cddf8d76SBarry Smith } 503cddf8d76SBarry Smith } 504cddf8d76SBarry Smith color = DRAW_CYAN; 505cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 506cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 507cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 508cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 509cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 510cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 511cddf8d76SBarry Smith } 512cddf8d76SBarry Smith } 513cddf8d76SBarry Smith color = DRAW_RED; 514cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 515cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 516cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 517cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 518cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 519cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 520cddf8d76SBarry Smith #else 521cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 522cddf8d76SBarry Smith #endif 523cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 524cddf8d76SBarry Smith } 525cddf8d76SBarry Smith } 526*0513a670SBarry Smith } else { 527*0513a670SBarry Smith /* use contour shading to indicate magnitude of values */ 528*0513a670SBarry Smith int count = 0; 529*0513a670SBarry Smith for ( i=0; i<m; i++ ) { 530*0513a670SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 531*0513a670SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 532*0513a670SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 533*0513a670SBarry Smith color = 32 + (int) ((200.0 - 32.0)*PetscAbsScalar(a->a[count])/maxv); 534*0513a670SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 535*0513a670SBarry Smith count++; 536*0513a670SBarry Smith } 537*0513a670SBarry Smith } 538*0513a670SBarry Smith } 539*0513a670SBarry Smith 5406945ee14SBarry Smith ierr = DrawCheckResizedWindow(draw); 541cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 542cddf8d76SBarry Smith } 543416022c9SBarry Smith return 0; 544416022c9SBarry Smith } 545416022c9SBarry Smith 5465615d1e5SSatish Balay #undef __FUNC__ 5475615d1e5SSatish Balay #define __FUNC__ "MatView_SeqAIJ" 548416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 549416022c9SBarry Smith { 550416022c9SBarry Smith Mat A = (Mat) obj; 551416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 552bcd2baecSBarry Smith ViewerType vtype; 553bcd2baecSBarry Smith int ierr; 554416022c9SBarry Smith 555bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 556bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 557416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 558416022c9SBarry Smith } 559bcd2baecSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){ 560416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 561416022c9SBarry Smith } 562bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 563416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 564416022c9SBarry Smith } 565bcd2baecSBarry Smith else if (vtype == DRAW_VIEWER) { 566bcd2baecSBarry Smith return MatView_SeqAIJ_Draw(A,viewer); 56717ab2063SBarry Smith } 56817ab2063SBarry Smith return 0; 56917ab2063SBarry Smith } 57019bcc07fSBarry Smith 571c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 5725615d1e5SSatish Balay #undef __FUNC__ 5735615d1e5SSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqAIJ" 574416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 57517ab2063SBarry Smith { 576416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 57741c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 578416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 579416022c9SBarry Smith Scalar *aa = a->a, *ap; 58017ab2063SBarry Smith 5816d4a8577SBarry Smith if (mode == MAT_FLUSH_ASSEMBLY) return 0; 58217ab2063SBarry Smith 58317ab2063SBarry Smith for ( i=1; i<m; i++ ) { 584416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 58517ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 58617ab2063SBarry Smith if (fshift) { 587416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 58817ab2063SBarry Smith N = ailen[i]; 58917ab2063SBarry Smith for ( j=0; j<N; j++ ) { 59017ab2063SBarry Smith ip[j-fshift] = ip[j]; 59117ab2063SBarry Smith ap[j-fshift] = ap[j]; 59217ab2063SBarry Smith } 59317ab2063SBarry Smith } 59417ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 59517ab2063SBarry Smith } 59617ab2063SBarry Smith if (m) { 59717ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 59817ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 59917ab2063SBarry Smith } 60017ab2063SBarry Smith /* reset ilen and imax for each row */ 60117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 60217ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 60317ab2063SBarry Smith } 604416022c9SBarry Smith a->nz = ai[m] + shift; 60517ab2063SBarry Smith 60617ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 607416022c9SBarry Smith if (fshift && a->diag) { 6080452661fSBarry Smith PetscFree(a->diag); 609416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 610416022c9SBarry Smith a->diag = 0; 61117ab2063SBarry Smith } 6124e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Matrix size: %d X %d; storage space: %d unneeded, %d used\n", 6134e220ebcSLois Curfman McInnes m,a->n,fshift,a->nz); 6144e220ebcSLois Curfman McInnes PLogInfo(A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues is %d\n", 615b810aeb4SBarry Smith a->reallocs); 616dd5f02e7SSatish Balay a->reallocs = 0; 6174e220ebcSLois Curfman McInnes A->info.nz_unneeded = (double)fshift; 6184e220ebcSLois Curfman McInnes 61976dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 62041c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 62117ab2063SBarry Smith return 0; 62217ab2063SBarry Smith } 62317ab2063SBarry Smith 6245615d1e5SSatish Balay #undef __FUNC__ 6255615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqAIJ" 626416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 62717ab2063SBarry Smith { 628416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 629cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 63017ab2063SBarry Smith return 0; 63117ab2063SBarry Smith } 632416022c9SBarry Smith 6335615d1e5SSatish Balay #undef __FUNC__ 6345615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqAIJ" 63517ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 63617ab2063SBarry Smith { 637416022c9SBarry Smith Mat A = (Mat) obj; 638416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 63990f02eecSBarry Smith int ierr; 640d5d45c9bSBarry Smith 64117ab2063SBarry Smith #if defined(PETSC_LOG) 642416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 64317ab2063SBarry Smith #endif 6440452661fSBarry Smith PetscFree(a->a); 6450452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 6460452661fSBarry Smith if (a->diag) PetscFree(a->diag); 6470452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 6480452661fSBarry Smith if (a->imax) PetscFree(a->imax); 6490452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 65076dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 6510452661fSBarry Smith PetscFree(a); 65290f02eecSBarry Smith if (A->mapping) { 65390f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(A->mapping); CHKERRQ(ierr); 65490f02eecSBarry Smith } 655f2655603SLois Curfman McInnes PLogObjectDestroy(A); 656f2655603SLois Curfman McInnes PetscHeaderDestroy(A); 65717ab2063SBarry Smith return 0; 65817ab2063SBarry Smith } 65917ab2063SBarry Smith 6605615d1e5SSatish Balay #undef __FUNC__ 6615615d1e5SSatish Balay #define __FUNC__ "MatCompress_SeqAIJ" 662416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 66317ab2063SBarry Smith { 66417ab2063SBarry Smith return 0; 66517ab2063SBarry Smith } 66617ab2063SBarry Smith 6675615d1e5SSatish Balay #undef __FUNC__ 6685615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqAIJ" 669416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 67017ab2063SBarry Smith { 671416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6726d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) a->roworiented = 1; 6736d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) a->roworiented = 0; 6746d4a8577SBarry Smith else if (op == MAT_COLUMNS_SORTED) a->sorted = 1; 675219d9a1aSLois Curfman McInnes else if (op == MAT_COLUMNS_UNSORTED) a->sorted = 0; 6766d4a8577SBarry Smith else if (op == MAT_NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 6776d4a8577SBarry Smith else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 6786d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 679219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 6806d4a8577SBarry Smith op == MAT_SYMMETRIC || 6816d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 68290f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 68390f02eecSBarry Smith op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES) 68494a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 6856d4a8577SBarry Smith else if (op == MAT_NO_NEW_DIAGONALS) 686e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");} 6876d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_1) a->inode.limit = 1; 6886d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_2) a->inode.limit = 2; 6896d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_3) a->inode.limit = 3; 6906d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_4) a->inode.limit = 4; 6916d4a8577SBarry Smith else if (op == MAT_INODE_LIMIT_5) a->inode.limit = 5; 692e2f28af5SBarry Smith else 693e3372554SBarry Smith {SETERRQ(PETSC_ERR_SUP,0,"unknown option");} 69417ab2063SBarry Smith return 0; 69517ab2063SBarry Smith } 69617ab2063SBarry Smith 6975615d1e5SSatish Balay #undef __FUNC__ 6985615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqAIJ" 699416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 70017ab2063SBarry Smith { 701416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 702416022c9SBarry Smith int i,j, n,shift = a->indexshift; 70317ab2063SBarry Smith Scalar *x, zero = 0.0; 70417ab2063SBarry Smith 70517ab2063SBarry Smith VecSet(&zero,v); 70690f02eecSBarry Smith VecGetArray_Fast(v,x); VecGetLocalSize(v,&n); 707e3372554SBarry Smith if (n != a->m) SETERRQ(1,0,"Nonconforming matrix and vector"); 708416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 709416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 710416022c9SBarry Smith if (a->j[j]+shift == i) { 711416022c9SBarry Smith x[i] = a->a[j]; 71217ab2063SBarry Smith break; 71317ab2063SBarry Smith } 71417ab2063SBarry Smith } 71517ab2063SBarry Smith } 71617ab2063SBarry Smith return 0; 71717ab2063SBarry Smith } 71817ab2063SBarry Smith 71917ab2063SBarry Smith /* -------------------------------------------------------*/ 72017ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 72117ab2063SBarry Smith /* -------------------------------------------------------*/ 7225615d1e5SSatish Balay #undef __FUNC__ 7235615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqAIJ" 72444cd7ae7SLois Curfman McInnes int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 72517ab2063SBarry Smith { 726416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 72717ab2063SBarry Smith Scalar *x, *y, *v, alpha; 728416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 72917ab2063SBarry Smith 73090f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 731cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 73217ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 73317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 734416022c9SBarry Smith idx = a->j + a->i[i] + shift; 735416022c9SBarry Smith v = a->a + a->i[i] + shift; 736416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 73717ab2063SBarry Smith alpha = x[i]; 73817ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 73917ab2063SBarry Smith } 740416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 74117ab2063SBarry Smith return 0; 74217ab2063SBarry Smith } 743d5d45c9bSBarry Smith 7445615d1e5SSatish Balay #undef __FUNC__ 7455615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqAIJ" 74644cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 74717ab2063SBarry Smith { 748416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 74917ab2063SBarry Smith Scalar *x, *y, *v, alpha; 750416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 75117ab2063SBarry Smith 75290f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 75317ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 75417ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 75517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 756416022c9SBarry Smith idx = a->j + a->i[i] + shift; 757416022c9SBarry Smith v = a->a + a->i[i] + shift; 758416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 75917ab2063SBarry Smith alpha = x[i]; 76017ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 76117ab2063SBarry Smith } 76290f02eecSBarry Smith PLogFlops(2*a->nz); 76317ab2063SBarry Smith return 0; 76417ab2063SBarry Smith } 76517ab2063SBarry Smith 7665615d1e5SSatish Balay #undef __FUNC__ 7675615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqAIJ" 76844cd7ae7SLois Curfman McInnes int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 76917ab2063SBarry Smith { 770416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 77117ab2063SBarry Smith Scalar *x, *y, *v, sum; 772416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 77317ab2063SBarry Smith 77490f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); 77517ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 776416022c9SBarry Smith idx = a->j; 777416022c9SBarry Smith v = a->a; 778416022c9SBarry Smith ii = a->i; 77917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 780416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 78117ab2063SBarry Smith sum = 0.0; 78217ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 78317ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 784416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 78517ab2063SBarry Smith y[i] = sum; 78617ab2063SBarry Smith } 787416022c9SBarry Smith PLogFlops(2*a->nz - m); 78817ab2063SBarry Smith return 0; 78917ab2063SBarry Smith } 79017ab2063SBarry Smith 7915615d1e5SSatish Balay #undef __FUNC__ 7925615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqAIJ" 79344cd7ae7SLois Curfman McInnes int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 79417ab2063SBarry Smith { 795416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 79617ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 797cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 79817ab2063SBarry Smith 79990f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(yy,y); VecGetArray_Fast(zz,z); 80017ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 801cddf8d76SBarry Smith idx = a->j; 802cddf8d76SBarry Smith v = a->a; 803cddf8d76SBarry Smith ii = a->i; 80417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 805cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 80617ab2063SBarry Smith sum = y[i]; 807cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 808cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 80917ab2063SBarry Smith z[i] = sum; 81017ab2063SBarry Smith } 811416022c9SBarry Smith PLogFlops(2*a->nz); 81217ab2063SBarry Smith return 0; 81317ab2063SBarry Smith } 81417ab2063SBarry Smith 81517ab2063SBarry Smith /* 81617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 81717ab2063SBarry Smith */ 81817ab2063SBarry Smith 8195615d1e5SSatish Balay #undef __FUNC__ 8205615d1e5SSatish Balay #define __FUNC__ "MatMarkDiag_SeqAIJ" 821416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 82217ab2063SBarry Smith { 823416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 824416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 82517ab2063SBarry Smith 8260452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 827416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 828416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 829416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 830416022c9SBarry Smith if (a->j[j]+shift == i) { 83117ab2063SBarry Smith diag[i] = j - shift; 83217ab2063SBarry Smith break; 83317ab2063SBarry Smith } 83417ab2063SBarry Smith } 83517ab2063SBarry Smith } 836416022c9SBarry Smith a->diag = diag; 83717ab2063SBarry Smith return 0; 83817ab2063SBarry Smith } 83917ab2063SBarry Smith 8405615d1e5SSatish Balay #undef __FUNC__ 8415615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqAIJ" 84244cd7ae7SLois Curfman McInnes int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 84317ab2063SBarry Smith double fshift,int its,Vec xx) 84417ab2063SBarry Smith { 845416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 846416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 847d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 84817ab2063SBarry Smith 84990f02eecSBarry Smith VecGetArray_Fast(xx,x); VecGetArray_Fast(bb,b); 850416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 851416022c9SBarry Smith diag = a->diag; 852416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 85317ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 85417ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 85517ab2063SBarry Smith bs = b + shift; 85617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 857416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 858416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 859416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 860416022c9SBarry Smith v = a->a + diag[i] + (!shift); 86117ab2063SBarry Smith sum = b[i]*d/omega; 86217ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 86317ab2063SBarry Smith x[i] = sum; 86417ab2063SBarry Smith } 86517ab2063SBarry Smith return 0; 86617ab2063SBarry Smith } 86717ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 868e3372554SBarry Smith SETERRQ(1,0,"SOR_APPLY_LOWER is not done"); 86917ab2063SBarry Smith } 870416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 87117ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 87217ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 87317ab2063SBarry Smith 87417ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 87517ab2063SBarry Smith 87617ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 87717ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 87817ab2063SBarry Smith is the relaxation factor. 87917ab2063SBarry Smith */ 8800452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 88117ab2063SBarry Smith scale = (2.0/omega) - 1.0; 88217ab2063SBarry Smith 88317ab2063SBarry Smith /* x = (E + U)^{-1} b */ 88417ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 885416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 886416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 887416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 888416022c9SBarry Smith v = a->a + diag[i] + (!shift); 88917ab2063SBarry Smith sum = b[i]; 89017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 89117ab2063SBarry Smith x[i] = omega*(sum/d); 89217ab2063SBarry Smith } 89317ab2063SBarry Smith 89417ab2063SBarry Smith /* t = b - (2*E - D)x */ 895416022c9SBarry Smith v = a->a; 89617ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 89717ab2063SBarry Smith 89817ab2063SBarry Smith /* t = (E + L)^{-1}t */ 899416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 900416022c9SBarry Smith diag = a->diag; 90117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 902416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 903416022c9SBarry Smith n = diag[i] - a->i[i]; 904416022c9SBarry Smith idx = a->j + a->i[i] + shift; 905416022c9SBarry Smith v = a->a + a->i[i] + shift; 90617ab2063SBarry Smith sum = t[i]; 90717ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 90817ab2063SBarry Smith t[i] = omega*(sum/d); 90917ab2063SBarry Smith } 91017ab2063SBarry Smith 91117ab2063SBarry Smith /* x = x + t */ 91217ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 9130452661fSBarry Smith PetscFree(t); 91417ab2063SBarry Smith return 0; 91517ab2063SBarry Smith } 91617ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 91717ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 91817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 919416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 920416022c9SBarry Smith n = diag[i] - a->i[i]; 921416022c9SBarry Smith idx = a->j + a->i[i] + shift; 922416022c9SBarry Smith v = a->a + a->i[i] + shift; 92317ab2063SBarry Smith sum = b[i]; 92417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 92517ab2063SBarry Smith x[i] = omega*(sum/d); 92617ab2063SBarry Smith } 92717ab2063SBarry Smith xb = x; 92817ab2063SBarry Smith } 92917ab2063SBarry Smith else xb = b; 93017ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 93117ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 93217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 933416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 93417ab2063SBarry Smith } 93517ab2063SBarry Smith } 93617ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 93717ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 938416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 939416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 940416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 941416022c9SBarry Smith v = a->a + diag[i] + (!shift); 94217ab2063SBarry Smith sum = xb[i]; 94317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 94417ab2063SBarry Smith x[i] = omega*(sum/d); 94517ab2063SBarry Smith } 94617ab2063SBarry Smith } 94717ab2063SBarry Smith its--; 94817ab2063SBarry Smith } 94917ab2063SBarry Smith while (its--) { 95017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 95117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 952416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 953416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 954416022c9SBarry Smith idx = a->j + a->i[i] + shift; 955416022c9SBarry Smith v = a->a + a->i[i] + shift; 95617ab2063SBarry Smith sum = b[i]; 95717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 9587e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 95917ab2063SBarry Smith } 96017ab2063SBarry Smith } 96117ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 96217ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 963416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 964416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 965416022c9SBarry Smith idx = a->j + a->i[i] + shift; 966416022c9SBarry Smith v = a->a + a->i[i] + shift; 96717ab2063SBarry Smith sum = b[i]; 96817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 9697e33a6baSBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum + a->a[diag[i]+shift]*x[i])/d; 97017ab2063SBarry Smith } 97117ab2063SBarry Smith } 97217ab2063SBarry Smith } 97317ab2063SBarry Smith return 0; 97417ab2063SBarry Smith } 97517ab2063SBarry Smith 9765615d1e5SSatish Balay #undef __FUNC__ 9775615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqAIJ" 9784e220ebcSLois Curfman McInnes static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) 97917ab2063SBarry Smith { 980416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 9814e220ebcSLois Curfman McInnes 9824e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 9834e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 9844e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 9854e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 9864e220ebcSLois Curfman McInnes info->block_size = 1.0; 9874e220ebcSLois Curfman McInnes info->nz_allocated = (double)a->maxnz; 9884e220ebcSLois Curfman McInnes info->nz_used = (double)a->nz; 9894e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(a->maxnz - a->nz); 9904e220ebcSLois Curfman McInnes /* if (info->nz_unneeded != A->info.nz_unneeded) 9914e220ebcSLois Curfman McInnes printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 9924e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 9934e220ebcSLois Curfman McInnes info->mallocs = (double)a->reallocs; 9944e220ebcSLois Curfman McInnes info->memory = A->mem; 9954e220ebcSLois Curfman McInnes if (A->factor) { 9964e220ebcSLois Curfman McInnes info->fill_ratio_given = A->info.fill_ratio_given; 9974e220ebcSLois Curfman McInnes info->fill_ratio_needed = A->info.fill_ratio_needed; 9984e220ebcSLois Curfman McInnes info->factor_mallocs = A->info.factor_mallocs; 9994e220ebcSLois Curfman McInnes } else { 10004e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 10014e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 10024e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 10034e220ebcSLois Curfman McInnes } 100417ab2063SBarry Smith return 0; 100517ab2063SBarry Smith } 100617ab2063SBarry Smith 100717ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 100817ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 100917ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 101017ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 101117ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 101217ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 101317ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 101417ab2063SBarry Smith 10155615d1e5SSatish Balay #undef __FUNC__ 10165615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqAIJ" 101717ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 101817ab2063SBarry Smith { 1019416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1020416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 102117ab2063SBarry Smith 102277c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 102317ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 102417ab2063SBarry Smith if (diag) { 102517ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1026e3372554SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range"); 1027416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 1028416022c9SBarry Smith a->ilen[rows[i]] = 1; 1029416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 1030416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 103117ab2063SBarry Smith } 103217ab2063SBarry Smith else { 103317ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 103417ab2063SBarry Smith CHKERRQ(ierr); 103517ab2063SBarry Smith } 103617ab2063SBarry Smith } 103717ab2063SBarry Smith } 103817ab2063SBarry Smith else { 103917ab2063SBarry Smith for ( i=0; i<N; i++ ) { 1040e3372554SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,0,"row out of range"); 1041416022c9SBarry Smith a->ilen[rows[i]] = 0; 104217ab2063SBarry Smith } 104317ab2063SBarry Smith } 104417ab2063SBarry Smith ISRestoreIndices(is,&rows); 104543a90d84SBarry Smith ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 104617ab2063SBarry Smith return 0; 104717ab2063SBarry Smith } 104817ab2063SBarry Smith 10495615d1e5SSatish Balay #undef __FUNC__ 10505615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqAIJ" 1051416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 105217ab2063SBarry Smith { 1053416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1054416022c9SBarry Smith *m = a->m; *n = a->n; 105517ab2063SBarry Smith return 0; 105617ab2063SBarry Smith } 105717ab2063SBarry Smith 10585615d1e5SSatish Balay #undef __FUNC__ 10595615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqAIJ" 1060416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 106117ab2063SBarry Smith { 1062416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1063416022c9SBarry Smith *m = 0; *n = a->m; 106417ab2063SBarry Smith return 0; 106517ab2063SBarry Smith } 1066026e39d0SSatish Balay 10675615d1e5SSatish Balay #undef __FUNC__ 10685615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqAIJ" 10694e093b46SBarry Smith int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 107017ab2063SBarry Smith { 1071416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1072c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 107317ab2063SBarry Smith 1074e3372554SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,0,"Row out of range"); 107517ab2063SBarry Smith 1076416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 1077416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 107817ab2063SBarry Smith if (idx) { 1079416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 10804e093b46SBarry Smith if (*nz && shift) { 10810452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 108217ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 10834e093b46SBarry Smith } else if (*nz) { 10844e093b46SBarry Smith *idx = itmp; 108517ab2063SBarry Smith } 108617ab2063SBarry Smith else *idx = 0; 108717ab2063SBarry Smith } 108817ab2063SBarry Smith return 0; 108917ab2063SBarry Smith } 109017ab2063SBarry Smith 10915615d1e5SSatish Balay #undef __FUNC__ 10925615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqAIJ" 10934e093b46SBarry Smith int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 109417ab2063SBarry Smith { 10954e093b46SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 10964e093b46SBarry Smith if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 109717ab2063SBarry Smith return 0; 109817ab2063SBarry Smith } 109917ab2063SBarry Smith 11005615d1e5SSatish Balay #undef __FUNC__ 11015615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqAIJ" 1102cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 110317ab2063SBarry Smith { 1104416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1105416022c9SBarry Smith Scalar *v = a->a; 110617ab2063SBarry Smith double sum = 0.0; 1107416022c9SBarry Smith int i, j,shift = a->indexshift; 110817ab2063SBarry Smith 110917ab2063SBarry Smith if (type == NORM_FROBENIUS) { 1110416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 111117ab2063SBarry Smith #if defined(PETSC_COMPLEX) 111217ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 111317ab2063SBarry Smith #else 111417ab2063SBarry Smith sum += (*v)*(*v); v++; 111517ab2063SBarry Smith #endif 111617ab2063SBarry Smith } 111717ab2063SBarry Smith *norm = sqrt(sum); 111817ab2063SBarry Smith } 111917ab2063SBarry Smith else if (type == NORM_1) { 112017ab2063SBarry Smith double *tmp; 1121416022c9SBarry Smith int *jj = a->j; 11220452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 1123cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 112417ab2063SBarry Smith *norm = 0.0; 1125416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 1126a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 112717ab2063SBarry Smith } 1128416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 112917ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 113017ab2063SBarry Smith } 11310452661fSBarry Smith PetscFree(tmp); 113217ab2063SBarry Smith } 113317ab2063SBarry Smith else if (type == NORM_INFINITY) { 113417ab2063SBarry Smith *norm = 0.0; 1135416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 1136416022c9SBarry Smith v = a->a + a->i[j] + shift; 113717ab2063SBarry Smith sum = 0.0; 1138416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1139cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 114017ab2063SBarry Smith } 114117ab2063SBarry Smith if (sum > *norm) *norm = sum; 114217ab2063SBarry Smith } 114317ab2063SBarry Smith } 114417ab2063SBarry Smith else { 1145e3372554SBarry Smith SETERRQ(1,0,"No support for two norm yet"); 114617ab2063SBarry Smith } 114717ab2063SBarry Smith return 0; 114817ab2063SBarry Smith } 114917ab2063SBarry Smith 11505615d1e5SSatish Balay #undef __FUNC__ 11515615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqAIJ" 1152416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 115317ab2063SBarry Smith { 1154416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1155416022c9SBarry Smith Mat C; 1156416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 1157416022c9SBarry Smith int shift = a->indexshift; 1158d5d45c9bSBarry Smith Scalar *array = a->a; 115917ab2063SBarry Smith 11603638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 1161e3372554SBarry Smith SETERRQ(1,0,"Square matrix only for in-place"); 11620452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 1163cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 116417ab2063SBarry Smith if (shift) { 116517ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 116617ab2063SBarry Smith } 116717ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 1168416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 11690452661fSBarry Smith PetscFree(col); 117017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 117117ab2063SBarry Smith len = ai[i+1]-ai[i]; 1172416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 117317ab2063SBarry Smith array += len; aj += len; 117417ab2063SBarry Smith } 117517ab2063SBarry Smith if (shift) { 117617ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 117717ab2063SBarry Smith } 117817ab2063SBarry Smith 11796d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 11806d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 118117ab2063SBarry Smith 11823638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 1183416022c9SBarry Smith *B = C; 118417ab2063SBarry Smith } else { 1185416022c9SBarry Smith /* This isn't really an in-place transpose */ 11860452661fSBarry Smith PetscFree(a->a); 11870452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 11880452661fSBarry Smith if (a->diag) PetscFree(a->diag); 11890452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 11900452661fSBarry Smith if (a->imax) PetscFree(a->imax); 11910452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 11921073c447SSatish Balay if (a->inode.size) PetscFree(a->inode.size); 11930452661fSBarry Smith PetscFree(a); 1194416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 11950452661fSBarry Smith PetscHeaderDestroy(C); 119617ab2063SBarry Smith } 119717ab2063SBarry Smith return 0; 119817ab2063SBarry Smith } 119917ab2063SBarry Smith 12005615d1e5SSatish Balay #undef __FUNC__ 12015615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqAIJ" 1202f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 120317ab2063SBarry Smith { 1204416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 120517ab2063SBarry Smith Scalar *l,*r,x,*v; 1206d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 120717ab2063SBarry Smith 120817ab2063SBarry Smith if (ll) { 12093ea7c6a1SSatish Balay /* The local size is used so that VecMPI can be passed to this routine 12103ea7c6a1SSatish Balay by MatDiagonalScale_MPIAIJ */ 12119b1297e1SSatish Balay VecGetLocalSize_Fast(ll,m); 1212e3372554SBarry Smith if (m != a->m) SETERRQ(1,0,"Left scaling vector wrong length"); 121390f02eecSBarry Smith VecGetArray_Fast(ll,l); 1214416022c9SBarry Smith v = a->a; 121517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 121617ab2063SBarry Smith x = l[i]; 1217416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 121817ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 121917ab2063SBarry Smith } 122044cd7ae7SLois Curfman McInnes PLogFlops(nz); 122117ab2063SBarry Smith } 122217ab2063SBarry Smith if (rr) { 12239b1297e1SSatish Balay VecGetLocalSize_Fast(rr,n); 1224e3372554SBarry Smith if (n != a->n) SETERRQ(1,0,"Right scaling vector wrong length"); 122590f02eecSBarry Smith VecGetArray_Fast(rr,r); 1226416022c9SBarry Smith v = a->a; jj = a->j; 122717ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 122817ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 122917ab2063SBarry Smith } 123044cd7ae7SLois Curfman McInnes PLogFlops(nz); 123117ab2063SBarry Smith } 123217ab2063SBarry Smith return 0; 123317ab2063SBarry Smith } 123417ab2063SBarry Smith 12355615d1e5SSatish Balay #undef __FUNC__ 12365615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqAIJ" 1237cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 123817ab2063SBarry Smith { 1239db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 124002834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 124199141d43SSatish Balay int row,mat_i,*mat_j,tcol,first,step,*mat_ilen; 1242a2744918SBarry Smith register int sum,lensi; 124399141d43SSatish Balay int *irow, *icol, nrows, ncols, shift = a->indexshift,*ssmap; 124499141d43SSatish Balay int *starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 124599141d43SSatish Balay Scalar *a_new,*mat_a; 1246416022c9SBarry Smith Mat C; 124717ab2063SBarry Smith 1248b48a1e75SSatish Balay ierr = ISSorted(isrow,(PetscTruth*)&i); 1249e3372554SBarry Smith if (!i) SETERRQ(1,0,"ISrow is not sorted"); 125099141d43SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i); 1251e3372554SBarry Smith if (!i) SETERRQ(1,0,"IScol is not sorted"); 125299141d43SSatish Balay 125317ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 125417ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 125517ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 125617ab2063SBarry Smith 12577264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 125802834360SBarry Smith /* special case of contiguous rows */ 125957faeb66SBarry Smith lens = (int *) PetscMalloc((ncols+nrows+1)*sizeof(int)); CHKPTRQ(lens); 126002834360SBarry Smith starts = lens + ncols; 126102834360SBarry Smith /* loop over new rows determining lens and starting points */ 126202834360SBarry Smith for (i=0; i<nrows; i++) { 1263a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1264a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 126502834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1266d8ced48eSBarry Smith if (aj[k]+shift >= first) { 126702834360SBarry Smith starts[i] = k; 126802834360SBarry Smith break; 126902834360SBarry Smith } 127002834360SBarry Smith } 1271a2744918SBarry Smith sum = 0; 127202834360SBarry Smith while (k < kend) { 1273d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1274a2744918SBarry Smith sum++; 127502834360SBarry Smith } 1276a2744918SBarry Smith lens[i] = sum; 127702834360SBarry Smith } 127802834360SBarry Smith /* create submatrix */ 1279cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 128008480c60SBarry Smith int n_cols,n_rows; 128108480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1282e3372554SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,0,""); 1283d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 128408480c60SBarry Smith C = *B; 128508480c60SBarry Smith } 128608480c60SBarry Smith else { 128702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 128808480c60SBarry Smith } 1289db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1290db02288aSLois Curfman McInnes 129102834360SBarry Smith /* loop over rows inserting into submatrix */ 1292db02288aSLois Curfman McInnes a_new = c->a; 1293db02288aSLois Curfman McInnes j_new = c->j; 1294db02288aSLois Curfman McInnes i_new = c->i; 1295db02288aSLois Curfman McInnes i_new[0] = -shift; 129602834360SBarry Smith for (i=0; i<nrows; i++) { 1297a2744918SBarry Smith ii = starts[i]; 1298a2744918SBarry Smith lensi = lens[i]; 1299a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1300a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 130102834360SBarry Smith } 1302a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1303a2744918SBarry Smith a_new += lensi; 1304a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1305a2744918SBarry Smith c->ilen[i] = lensi; 130602834360SBarry Smith } 13070452661fSBarry Smith PetscFree(lens); 130802834360SBarry Smith } 130902834360SBarry Smith else { 131002834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 13110452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 131202834360SBarry Smith ssmap = smap + shift; 131399141d43SSatish Balay lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 1314cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 131517ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 131602834360SBarry Smith /* determine lens of each row */ 131702834360SBarry Smith for (i=0; i<nrows; i++) { 1318d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 131902834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 132002834360SBarry Smith lens[i] = 0; 132102834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1322d8ced48eSBarry Smith if (ssmap[aj[k]]) { 132302834360SBarry Smith lens[i]++; 132402834360SBarry Smith } 132502834360SBarry Smith } 132602834360SBarry Smith } 132717ab2063SBarry Smith /* Create and fill new matrix */ 1328a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 132999141d43SSatish Balay c = (Mat_SeqAIJ *)((*B)->data); 133099141d43SSatish Balay 1331e3372554SBarry Smith if (c->m != nrows || c->n != ncols) SETERRQ(1,0,""); 133299141d43SSatish Balay if (PetscMemcmp(c->ilen,lens, c->m *sizeof(int))) { 1333e3372554SBarry Smith SETERRQ(1,0,"Cannot reuse matrix. wrong no of nonzeros"); 133499141d43SSatish Balay } 133599141d43SSatish Balay PetscMemzero(c->ilen,c->m*sizeof(int)); 133608480c60SBarry Smith C = *B; 133708480c60SBarry Smith } 133808480c60SBarry Smith else { 133902834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 134008480c60SBarry Smith } 134199141d43SSatish Balay c = (Mat_SeqAIJ *)(C->data); 134217ab2063SBarry Smith for (i=0; i<nrows; i++) { 134399141d43SSatish Balay row = irow[i]; 134417ab2063SBarry Smith nznew = 0; 134599141d43SSatish Balay kstart = ai[row]+shift; 134699141d43SSatish Balay kend = kstart + a->ilen[row]; 134799141d43SSatish Balay mat_i = c->i[i]+shift; 134899141d43SSatish Balay mat_j = c->j + mat_i; 134999141d43SSatish Balay mat_a = c->a + mat_i; 135099141d43SSatish Balay mat_ilen = c->ilen + i; 135117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 135299141d43SSatish Balay if ((tcol=ssmap[a->j[k]])) { 135399141d43SSatish Balay *mat_j++ = tcol - (!shift); 135499141d43SSatish Balay *mat_a++ = a->a[k]; 135599141d43SSatish Balay (*mat_ilen)++; 135699141d43SSatish Balay 135717ab2063SBarry Smith } 135817ab2063SBarry Smith } 135917ab2063SBarry Smith } 136002834360SBarry Smith /* Free work space */ 136102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 136299141d43SSatish Balay PetscFree(smap); PetscFree(lens); 136302834360SBarry Smith } 13646d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 13656d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 136617ab2063SBarry Smith 136717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1368416022c9SBarry Smith *B = C; 136917ab2063SBarry Smith return 0; 137017ab2063SBarry Smith } 137117ab2063SBarry Smith 1372a871dcd8SBarry Smith /* 137363b91edcSBarry Smith note: This can only work for identity for row and col. It would 137463b91edcSBarry Smith be good to check this and otherwise generate an error. 1375a871dcd8SBarry Smith */ 13765615d1e5SSatish Balay #undef __FUNC__ 13775615d1e5SSatish Balay #define __FUNC__ "MatILUFactor_SeqAIJ" 137863b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1379a871dcd8SBarry Smith { 138063b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 138108480c60SBarry Smith int ierr; 138263b91edcSBarry Smith Mat outA; 138363b91edcSBarry Smith 1384e3372554SBarry Smith if (fill != 0) SETERRQ(1,0,"Only fill=0 supported"); 1385a871dcd8SBarry Smith 138663b91edcSBarry Smith outA = inA; 138763b91edcSBarry Smith inA->factor = FACTOR_LU; 138863b91edcSBarry Smith a->row = row; 138963b91edcSBarry Smith a->col = col; 139063b91edcSBarry Smith 13910452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 139263b91edcSBarry Smith 139308480c60SBarry Smith if (!a->diag) { 139408480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 139563b91edcSBarry Smith } 139663b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1397a871dcd8SBarry Smith return 0; 1398a871dcd8SBarry Smith } 1399a871dcd8SBarry Smith 1400f0b747eeSBarry Smith #include "pinclude/plapack.h" 14015615d1e5SSatish Balay #undef __FUNC__ 14025615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqAIJ" 1403f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1404f0b747eeSBarry Smith { 1405f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1406f0b747eeSBarry Smith int one = 1; 1407f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1408f0b747eeSBarry Smith PLogFlops(a->nz); 1409f0b747eeSBarry Smith return 0; 1410f0b747eeSBarry Smith } 1411f0b747eeSBarry Smith 14125615d1e5SSatish Balay #undef __FUNC__ 14135615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqAIJ" 1414cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1415cddf8d76SBarry Smith Mat **B) 1416cddf8d76SBarry Smith { 1417cddf8d76SBarry Smith int ierr,i; 1418cddf8d76SBarry Smith 1419cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 14200452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1421cddf8d76SBarry Smith } 1422cddf8d76SBarry Smith 1423cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1424905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1425cddf8d76SBarry Smith } 1426cddf8d76SBarry Smith return 0; 1427cddf8d76SBarry Smith } 1428cddf8d76SBarry Smith 14295615d1e5SSatish Balay #undef __FUNC__ 14305615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqAIJ" 14315a838052SSatish Balay static int MatGetBlockSize_SeqAIJ(Mat A, int *bs) 14325a838052SSatish Balay { 14335a838052SSatish Balay *bs = 1; 14345a838052SSatish Balay return 0; 14355a838052SSatish Balay } 14365a838052SSatish Balay 14375615d1e5SSatish Balay #undef __FUNC__ 14385615d1e5SSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqAIJ" 1439e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 14404dcbc457SBarry Smith { 1441e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 144206763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 14438a047759SSatish Balay int start, end, *ai, *aj; 144406763907SSatish Balay char *table; 14458a047759SSatish Balay shift = a->indexshift; 1446e4d965acSSatish Balay m = a->m; 1447e4d965acSSatish Balay ai = a->i; 14488a047759SSatish Balay aj = a->j+shift; 14498a047759SSatish Balay 1450e3372554SBarry Smith if (ov < 0) SETERRQ(1,0,"illegal overlap value used"); 145106763907SSatish Balay 145206763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 145306763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 145406763907SSatish Balay 1455e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1456e4d965acSSatish Balay /* Initialise the two local arrays */ 1457e4d965acSSatish Balay isz = 0; 145806763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1459e4d965acSSatish Balay 1460e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 14614dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 146277c4ece6SBarry Smith ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 1463e4d965acSSatish Balay 1464e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1465e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 146606763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 14674dcbc457SBarry Smith } 146806763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 146906763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1470e4d965acSSatish Balay 147104a348a9SBarry Smith k = 0; 147204a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 147304a348a9SBarry Smith n = isz; 147406763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1475e4d965acSSatish Balay row = nidx[k]; 1476e4d965acSSatish Balay start = ai[row]; 1477e4d965acSSatish Balay end = ai[row+1]; 147804a348a9SBarry Smith for ( l = start; l<end ; l++){ 14798a047759SSatish Balay val = aj[l] + shift; 148006763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1481e4d965acSSatish Balay } 1482e4d965acSSatish Balay } 1483e4d965acSSatish Balay } 1484537820f0SBarry Smith ierr = ISCreateGeneral(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1485e4d965acSSatish Balay } 148604a348a9SBarry Smith PetscFree(table); 148706763907SSatish Balay PetscFree(nidx); 1488e4d965acSSatish Balay return 0; 14894dcbc457SBarry Smith } 149017ab2063SBarry Smith 1491*0513a670SBarry Smith /* -------------------------------------------------------------- */ 1492*0513a670SBarry Smith #undef __FUNC__ 1493*0513a670SBarry Smith #define __FUNC__ "MatPermute_SeqAIJ" 1494*0513a670SBarry Smith int MatPermute_SeqAIJ(Mat A, IS rowp, IS colp, Mat *B) 1495*0513a670SBarry Smith { 1496*0513a670SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1497*0513a670SBarry Smith Scalar *vwork; 1498*0513a670SBarry Smith int i, ierr, nz, m = a->m, n = a->n, *cwork; 1499*0513a670SBarry Smith int *row,*col,*cnew,j,*lens; 1500*0513a670SBarry Smith 1501*0513a670SBarry Smith ierr = ISGetIndices(rowp,&row); CHKERRQ(ierr); 1502*0513a670SBarry Smith ierr = ISGetIndices(colp,&col); CHKERRQ(ierr); 1503*0513a670SBarry Smith 1504*0513a670SBarry Smith /* determine lengths of permuted rows */ 1505*0513a670SBarry Smith lens = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(lens); 1506*0513a670SBarry Smith for (i=0; i<m; i++ ) { 1507*0513a670SBarry Smith lens[row[i]] = a->i[i+1] - a->i[i]; 1508*0513a670SBarry Smith } 1509*0513a670SBarry Smith ierr = MatCreateSeqAIJ(A->comm,m,n,0,lens,B);CHKERRQ(ierr); 1510*0513a670SBarry Smith PetscFree(lens); 1511*0513a670SBarry Smith 1512*0513a670SBarry Smith cnew = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(cnew); 1513*0513a670SBarry Smith for (i=0; i<m; i++) { 1514*0513a670SBarry Smith ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 1515*0513a670SBarry Smith for (j=0; j<nz; j++ ) { cnew[j] = col[cwork[j]];} 1516*0513a670SBarry Smith ierr = MatSetValues(*B,1,&row[i],nz,cnew,vwork,INSERT_VALUES); CHKERRQ(ierr); 1517*0513a670SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 1518*0513a670SBarry Smith } 1519*0513a670SBarry Smith PetscFree(cnew); 1520*0513a670SBarry Smith ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1521*0513a670SBarry Smith ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 1522*0513a670SBarry Smith ierr = ISRestoreIndices(rowp,&row); CHKERRQ(ierr); 1523*0513a670SBarry Smith ierr = ISRestoreIndices(colp,&col); CHKERRQ(ierr); 1524*0513a670SBarry Smith return 0; 1525*0513a670SBarry Smith } 1526*0513a670SBarry Smith 15275615d1e5SSatish Balay #undef __FUNC__ 15285615d1e5SSatish Balay #define __FUNC__ "MatPrintHelp_SeqAIJ" 1529682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1530682d7d0cSBarry Smith { 1531682d7d0cSBarry Smith static int called = 0; 1532682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1533682d7d0cSBarry Smith 1534682d7d0cSBarry Smith if (called) return 0; else called = 1; 153577c4ece6SBarry Smith PetscPrintf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 153677c4ece6SBarry Smith PetscPrintf(comm," -mat_lu_pivotthreshold <threshold>\n"); 153777c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 153877c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_no_inode - Do not use inodes\n"); 153977c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1540682d7d0cSBarry Smith #if defined(HAVE_ESSL) 154177c4ece6SBarry Smith PetscPrintf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1542682d7d0cSBarry Smith #endif 1543682d7d0cSBarry Smith return 0; 1544682d7d0cSBarry Smith } 1545205e38a3SBarry Smith static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg); 1546a93ec695SBarry Smith extern int MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring); 1547a93ec695SBarry Smith extern int MatColoringPatch_SeqAIJ(Mat,int,int *,ISColoring *); 1548a93ec695SBarry Smith 1549682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 155017ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 155117ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1552416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1553416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 155417ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 155517ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 155617ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 155717ab2063SBarry Smith MatRelax_SeqAIJ, 155817ab2063SBarry Smith MatTranspose_SeqAIJ, 15597264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1560f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 156117ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 156217ab2063SBarry Smith MatCompress_SeqAIJ, 156317ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 156417ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 156517ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 156617ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 156717ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 15683d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1569cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 15707eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1571682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1572f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 15735a838052SSatish Balay MatScale_SeqAIJ,0,0, 15746945ee14SBarry Smith MatILUDTFactor_SeqAIJ, 15756945ee14SBarry Smith MatGetBlockSize_SeqAIJ, 15763b2fbd54SBarry Smith MatGetRowIJ_SeqAIJ, 15773b2fbd54SBarry Smith MatRestoreRowIJ_SeqAIJ, 15783b2fbd54SBarry Smith MatGetColumnIJ_SeqAIJ, 1579a93ec695SBarry Smith MatRestoreColumnIJ_SeqAIJ, 1580a93ec695SBarry Smith MatFDColoringCreate_SeqAIJ, 1581*0513a670SBarry Smith MatColoringPatch_SeqAIJ, 1582*0513a670SBarry Smith 0, 1583*0513a670SBarry Smith MatPermute_SeqAIJ}; 158417ab2063SBarry Smith 158517ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 158617ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 158717ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 158817ab2063SBarry Smith 15895615d1e5SSatish Balay #undef __FUNC__ 15905615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqAIJ" 159117ab2063SBarry Smith /*@C 1592682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 15930d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 15946e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 15952bd5e0b2SLois Curfman McInnes (or the array nzz). By setting these parameters accurately, performance 15962bd5e0b2SLois Curfman McInnes during matrix assembly can be increased by more than a factor of 50. 159717ab2063SBarry Smith 159817ab2063SBarry Smith Input Parameters: 159917ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 160017ab2063SBarry Smith . m - number of rows 160117ab2063SBarry Smith . n - number of columns 160217ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 16032bd5e0b2SLois Curfman McInnes . nzz - array containing the number of nonzeros in the various rows 16042bd5e0b2SLois Curfman McInnes (possibly different for each row) or PETSC_NULL 160517ab2063SBarry Smith 160617ab2063SBarry Smith Output Parameter: 1607416022c9SBarry Smith . A - the matrix 160817ab2063SBarry Smith 160917ab2063SBarry Smith Notes: 161017ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 161117ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 16120002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 161344cd7ae7SLois Curfman McInnes either one (as in Fortran) or zero. See the users' manual for details. 161417ab2063SBarry Smith 161517ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1616a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 16173d323bbdSBarry Smith allocation. For large problems you MUST preallocate memory or you 16186da5968aSLois Curfman McInnes will get TERRIBLE performance, see the users' manual chapter on matrices. 161917ab2063SBarry Smith 1620682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1621682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1622682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 16236c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 16246c7ebb05SLois Curfman McInnes 16256c7ebb05SLois Curfman McInnes Options Database Keys: 16266c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 16276e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 16286e62573dSLois Curfman McInnes $ (max limit=5) 16296e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 16306e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 16316e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 163217ab2063SBarry Smith 163317ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 163417ab2063SBarry Smith @*/ 1635416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 163617ab2063SBarry Smith { 1637416022c9SBarry Smith Mat B; 1638416022c9SBarry Smith Mat_SeqAIJ *b; 16396945ee14SBarry Smith int i, len, ierr, flg,size; 16406945ee14SBarry Smith 16416945ee14SBarry Smith MPI_Comm_size(comm,&size); 1642e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"Comm must be of size 1"); 1643d5d45c9bSBarry Smith 1644416022c9SBarry Smith *A = 0; 16450452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1646416022c9SBarry Smith PLogObjectCreate(B); 16470452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 164844cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqAIJ)); 1649416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1650416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1651416022c9SBarry Smith B->view = MatView_SeqAIJ; 1652416022c9SBarry Smith B->factor = 0; 1653416022c9SBarry Smith B->lupivotthreshold = 1.0; 165490f02eecSBarry Smith B->mapping = 0; 16557a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 165669957df2SSatish Balay &flg); CHKERRQ(ierr); 16577a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 16587a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 16597a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1660416022c9SBarry Smith b->row = 0; 1661416022c9SBarry Smith b->col = 0; 1662416022c9SBarry Smith b->indexshift = 0; 1663b810aeb4SBarry Smith b->reallocs = 0; 166469957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 166569957df2SSatish Balay if (flg) b->indexshift = -1; 166617ab2063SBarry Smith 166744cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 166844cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 16690452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1670b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 16717b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 16727b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1673416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 167417ab2063SBarry Smith nz = nz*m; 167517ab2063SBarry Smith } 167617ab2063SBarry Smith else { 167717ab2063SBarry Smith nz = 0; 1678416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 167917ab2063SBarry Smith } 168017ab2063SBarry Smith 168117ab2063SBarry Smith /* allocate the matrix space */ 1682416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 16830452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1684416022c9SBarry Smith b->j = (int *) (b->a + nz); 1685cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1686416022c9SBarry Smith b->i = b->j + nz; 1687416022c9SBarry Smith b->singlemalloc = 1; 168817ab2063SBarry Smith 1689416022c9SBarry Smith b->i[0] = -b->indexshift; 169017ab2063SBarry Smith for (i=1; i<m+1; i++) { 1691416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 169217ab2063SBarry Smith } 169317ab2063SBarry Smith 1694416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 16950452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1696416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1697416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 169817ab2063SBarry Smith 1699416022c9SBarry Smith b->nz = 0; 1700416022c9SBarry Smith b->maxnz = nz; 1701416022c9SBarry Smith b->sorted = 0; 1702416022c9SBarry Smith b->roworiented = 1; 1703416022c9SBarry Smith b->nonew = 0; 1704416022c9SBarry Smith b->diag = 0; 1705416022c9SBarry Smith b->solve_work = 0; 170671bd300dSLois Curfman McInnes b->spptr = 0; 1707754ec7b1SSatish Balay b->inode.node_count = 0; 1708754ec7b1SSatish Balay b->inode.size = 0; 17096c7ebb05SLois Curfman McInnes b->inode.limit = 5; 17106c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 17114e220ebcSLois Curfman McInnes B->info.nz_unneeded = (double)b->maxnz; 171217ab2063SBarry Smith 1713416022c9SBarry Smith *A = B; 17144e220ebcSLois Curfman McInnes 17154b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 17164b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 171769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 171869957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 17194b14c69eSBarry Smith #endif 172069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 172169957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 172269957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 172369957df2SSatish Balay if (flg) { 1724e3372554SBarry Smith if (!b->indexshift) SETERRQ(1,0,"need -mat_aij_oneindex with -mat_aij_dxml"); 1725416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 172617ab2063SBarry Smith } 172769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 172869957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 172917ab2063SBarry Smith return 0; 173017ab2063SBarry Smith } 173117ab2063SBarry Smith 17325615d1e5SSatish Balay #undef __FUNC__ 17335615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqAIJ" 17343d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 173517ab2063SBarry Smith { 1736416022c9SBarry Smith Mat C; 1737416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 173808480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 173917ab2063SBarry Smith 17404043dd9cSLois Curfman McInnes *B = 0; 17410452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1742416022c9SBarry Smith PLogObjectCreate(C); 17430452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 174441c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1745416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1746416022c9SBarry Smith C->view = MatView_SeqAIJ; 1747416022c9SBarry Smith C->factor = A->factor; 1748416022c9SBarry Smith c->row = 0; 1749416022c9SBarry Smith c->col = 0; 1750416022c9SBarry Smith c->indexshift = shift; 1751c456f294SBarry Smith C->assembled = PETSC_TRUE; 175217ab2063SBarry Smith 175344cd7ae7SLois Curfman McInnes c->m = C->m = a->m; 175444cd7ae7SLois Curfman McInnes c->n = C->n = a->n; 175544cd7ae7SLois Curfman McInnes C->M = a->m; 175644cd7ae7SLois Curfman McInnes C->N = a->n; 175717ab2063SBarry Smith 17580452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 17590452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 176017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1761416022c9SBarry Smith c->imax[i] = a->imax[i]; 1762416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 176317ab2063SBarry Smith } 176417ab2063SBarry Smith 176517ab2063SBarry Smith /* allocate the matrix space */ 1766416022c9SBarry Smith c->singlemalloc = 1; 1767416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 17680452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1769416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1770416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1771416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 177217ab2063SBarry Smith if (m > 0) { 1773416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 177408480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1775416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 177617ab2063SBarry Smith } 177708480c60SBarry Smith } 177817ab2063SBarry Smith 1779416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1780416022c9SBarry Smith c->sorted = a->sorted; 1781416022c9SBarry Smith c->roworiented = a->roworiented; 1782416022c9SBarry Smith c->nonew = a->nonew; 17837a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 178417ab2063SBarry Smith 1785416022c9SBarry Smith if (a->diag) { 17860452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1787416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 178817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1789416022c9SBarry Smith c->diag[i] = a->diag[i]; 179017ab2063SBarry Smith } 179117ab2063SBarry Smith } 1792416022c9SBarry Smith else c->diag = 0; 17936c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 17946c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1795754ec7b1SSatish Balay if (a->inode.size){ 1796daed632aSSatish Balay c->inode.size = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->inode.size); 1797754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1798daed632aSSatish Balay PetscMemcpy( c->inode.size, a->inode.size, (m+1)*sizeof(int)); 1799754ec7b1SSatish Balay } else { 1800754ec7b1SSatish Balay c->inode.size = 0; 1801754ec7b1SSatish Balay c->inode.node_count = 0; 1802754ec7b1SSatish Balay } 1803416022c9SBarry Smith c->nz = a->nz; 1804416022c9SBarry Smith c->maxnz = a->maxnz; 1805416022c9SBarry Smith c->solve_work = 0; 180676dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1807754ec7b1SSatish Balay 1808416022c9SBarry Smith *B = C; 180917ab2063SBarry Smith return 0; 181017ab2063SBarry Smith } 181117ab2063SBarry Smith 18125615d1e5SSatish Balay #undef __FUNC__ 18135615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqAIJ" 181419bcc07fSBarry Smith int MatLoad_SeqAIJ(Viewer viewer,MatType type,Mat *A) 181517ab2063SBarry Smith { 1816416022c9SBarry Smith Mat_SeqAIJ *a; 1817416022c9SBarry Smith Mat B; 181817699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1819bcd2baecSBarry Smith MPI_Comm comm; 182017ab2063SBarry Smith 182119bcc07fSBarry Smith PetscObjectGetComm((PetscObject) viewer,&comm); 182217699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 1823e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"view must have one processor"); 182490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 182577c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 1826e3372554SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,0,"not matrix object in file"); 182717ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 182817ab2063SBarry Smith 182917ab2063SBarry Smith /* read in row lengths */ 18300452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 183177c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 183217ab2063SBarry Smith 183317ab2063SBarry Smith /* create our matrix */ 1834416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1835416022c9SBarry Smith B = *A; 1836416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1837416022c9SBarry Smith shift = a->indexshift; 183817ab2063SBarry Smith 183917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 184077c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->j,nz,BINARY_INT); CHKERRQ(ierr); 184117ab2063SBarry Smith if (shift) { 184217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1843416022c9SBarry Smith a->j[i] += 1; 184417ab2063SBarry Smith } 184517ab2063SBarry Smith } 184617ab2063SBarry Smith 184717ab2063SBarry Smith /* read in nonzero values */ 184877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->a,nz,BINARY_SCALAR); CHKERRQ(ierr); 184917ab2063SBarry Smith 185017ab2063SBarry Smith /* set matrix "i" values */ 1851416022c9SBarry Smith a->i[0] = -shift; 185217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1853416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1854416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 185517ab2063SBarry Smith } 18560452661fSBarry Smith PetscFree(rowlengths); 185717ab2063SBarry Smith 18586d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 18596d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 186017ab2063SBarry Smith return 0; 186117ab2063SBarry Smith } 186217ab2063SBarry Smith 18635615d1e5SSatish Balay #undef __FUNC__ 18645615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqAIJ" 1865686e14cbSSatish Balay static int MatEqual_SeqAIJ(Mat A,Mat B, PetscTruth* flg) 18667264ac53SSatish Balay { 18677264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 18687264ac53SSatish Balay 1869e3372554SBarry Smith if (B->type !=MATSEQAIJ)SETERRQ(1,0,"Matrices must be same type"); 18707264ac53SSatish Balay 18717264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 18727264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1873bcd2baecSBarry Smith (a->indexshift != b->indexshift)) { 187477c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1875bcd2baecSBarry Smith } 18767264ac53SSatish Balay 18777264ac53SSatish Balay /* if the a->i are the same */ 1878bcd2baecSBarry Smith if (PetscMemcmp(a->i,b->i, (a->n+1)*sizeof(int))) { 187977c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 18807264ac53SSatish Balay } 18817264ac53SSatish Balay 18827264ac53SSatish Balay /* if a->j are the same */ 1883bcd2baecSBarry Smith if (PetscMemcmp(a->j, b->j, (a->nz)*sizeof(int))) { 188477c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 1885bcd2baecSBarry Smith } 1886bcd2baecSBarry Smith 1887bcd2baecSBarry Smith /* if a->a are the same */ 188819bcc07fSBarry Smith if (PetscMemcmp(a->a, b->a, (a->nz)*sizeof(Scalar))) { 188977c4ece6SBarry Smith *flg = PETSC_FALSE; return 0; 18907264ac53SSatish Balay } 189177c4ece6SBarry Smith *flg = PETSC_TRUE; 18927264ac53SSatish Balay return 0; 18937264ac53SSatish Balay 18947264ac53SSatish Balay } 1895