117ab2063SBarry Smith #ifndef lint 2*7264ac53SSatish Balay static char vcid[] = "$Id: aij.c,v 1.147 1996/02/13 23:27:29 bsmith Exp balay $"; 317ab2063SBarry Smith #endif 417ab2063SBarry Smith 5d5d45c9bSBarry Smith /* 6d5d45c9bSBarry Smith Defines the basic matrix operations for the AIJ (compressed row) 7d5d45c9bSBarry Smith matrix storage format. 8d5d45c9bSBarry Smith */ 917ab2063SBarry Smith #include "aij.h" 1017ab2063SBarry Smith #include "vec/vecimpl.h" 1117ab2063SBarry Smith #include "inline/spops.h" 12e4d965acSSatish Balay #include "petsc.h" 1306763907SSatish Balay #include "inline/bitarray.h" 1417ab2063SBarry Smith 157a743949SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int**,int**); 1617ab2063SBarry Smith 17416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1817ab2063SBarry Smith { 19416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 20a2744918SBarry Smith int ierr, *ia, *ja,n,*idx,i; 213d54f168SSatish Balay /*Viewer V1, V2;*/ 2217ab2063SBarry Smith 23a2744918SBarry Smith /* 24a2744918SBarry Smith this is tacky: In the future when we have written special factorization 25a2744918SBarry Smith and solve routines for the identity permutation we should use a 26a2744918SBarry Smith stride index set instead of the general one. 27a2744918SBarry Smith */ 28a2744918SBarry Smith if (type == ORDER_NATURAL) { 29a2744918SBarry Smith n = a->n; 30a2744918SBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 31a2744918SBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 32a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 33a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 34a2744918SBarry Smith PetscFree(idx); 35a2744918SBarry Smith ISSetPermutation(*rperm); 36a2744918SBarry Smith ISSetPermutation(*cperm); 37a2744918SBarry Smith ISSetIdentity(*rperm); 38a2744918SBarry Smith ISSetIdentity(*cperm); 39a2744918SBarry Smith return 0; 40a2744918SBarry Smith } 41a2744918SBarry Smith 427a743949SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,a->indexshift, &ia, &ja);CHKERRQ(ierr); 43416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 44f31bba2fSSatish Balay 450452661fSBarry Smith PetscFree(ia); PetscFree(ja); 4617ab2063SBarry Smith return 0; 4717ab2063SBarry Smith } 4817ab2063SBarry Smith 49227d817aSBarry Smith #define CHUNKSIZE 15 5017ab2063SBarry Smith 5117ab2063SBarry Smith /* This version has row oriented v */ 52416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 5317ab2063SBarry Smith { 54416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 55416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 564b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 57d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 58416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 5917ab2063SBarry Smith 6017ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 61416022c9SBarry Smith row = im[k]; 6217ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 63416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 6417ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 6517ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 66416022c9SBarry Smith low = 0; 6717ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 68416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 69416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 704b0e389bSBarry Smith col = in[l] - shift; 714b0e389bSBarry Smith if (roworiented) { 724b0e389bSBarry Smith value = *v++; 734b0e389bSBarry Smith } 744b0e389bSBarry Smith else { 754b0e389bSBarry Smith value = v[k + l*m]; 764b0e389bSBarry Smith } 77416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 78416022c9SBarry Smith while (high-low > 5) { 79416022c9SBarry Smith t = (low+high)/2; 80416022c9SBarry Smith if (rp[t] > col) high = t; 81416022c9SBarry Smith else low = t; 8217ab2063SBarry Smith } 83416022c9SBarry Smith for ( i=low; i<high; i++ ) { 8417ab2063SBarry Smith if (rp[i] > col) break; 8517ab2063SBarry Smith if (rp[i] == col) { 86416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 8717ab2063SBarry Smith else ap[i] = value; 8817ab2063SBarry Smith goto noinsert; 8917ab2063SBarry Smith } 9017ab2063SBarry Smith } 9117ab2063SBarry Smith if (nonew) goto noinsert; 9217ab2063SBarry Smith if (nrow >= rmax) { 9317ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 94416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 9517ab2063SBarry Smith Scalar *new_a; 9617ab2063SBarry Smith 9717ab2063SBarry Smith /* malloc new storage space */ 98416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 990452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 10017ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 10117ab2063SBarry Smith new_i = new_j + new_nz; 10217ab2063SBarry Smith 10317ab2063SBarry Smith /* copy over old data into new slots */ 10417ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 105416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 106416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 107416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 108416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 10917ab2063SBarry Smith len*sizeof(int)); 110416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 111416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 11217ab2063SBarry Smith len*sizeof(Scalar)); 11317ab2063SBarry Smith /* free up old matrix storage */ 1140452661fSBarry Smith PetscFree(a->a); 1150452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 116416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 117416022c9SBarry Smith a->singlemalloc = 1; 11817ab2063SBarry Smith 11917ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 120416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 121416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 122416022c9SBarry Smith a->maxnz += CHUNKSIZE; 123b810aeb4SBarry Smith a->reallocs++; 12417ab2063SBarry Smith } 125416022c9SBarry Smith N = nrow++ - 1; a->nz++; 126416022c9SBarry Smith /* shift up all the later entries in this row */ 127416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 12817ab2063SBarry Smith rp[ii+1] = rp[ii]; 12917ab2063SBarry Smith ap[ii+1] = ap[ii]; 13017ab2063SBarry Smith } 13117ab2063SBarry Smith rp[i] = col; 13217ab2063SBarry Smith ap[i] = value; 13317ab2063SBarry Smith noinsert:; 134416022c9SBarry Smith low = i + 1; 13517ab2063SBarry Smith } 13617ab2063SBarry Smith ailen[row] = nrow; 13717ab2063SBarry Smith } 13817ab2063SBarry Smith return 0; 13917ab2063SBarry Smith } 14017ab2063SBarry Smith 1417eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 1427eb43aa7SLois Curfman McInnes { 1437eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 144b49de8d1SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 1457eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 1467eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 1477eb43aa7SLois Curfman McInnes 1487eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 1497eb43aa7SLois Curfman McInnes row = im[k]; 1507eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 1517eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 1527eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 1537eb43aa7SLois Curfman McInnes nrow = ailen[row]; 1547eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 1557eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 1567eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 1577eb43aa7SLois Curfman McInnes col = in[l] - shift; 1587eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 1597eb43aa7SLois Curfman McInnes while (high-low > 5) { 1607eb43aa7SLois Curfman McInnes t = (low+high)/2; 1617eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 1627eb43aa7SLois Curfman McInnes else low = t; 1637eb43aa7SLois Curfman McInnes } 1647eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 1657eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 1667eb43aa7SLois Curfman McInnes if (rp[i] == col) { 167b49de8d1SLois Curfman McInnes *v++ = ap[i]; 1687eb43aa7SLois Curfman McInnes goto finished; 1697eb43aa7SLois Curfman McInnes } 1707eb43aa7SLois Curfman McInnes } 171b49de8d1SLois Curfman McInnes *v++ = zero; 1727eb43aa7SLois Curfman McInnes finished:; 1737eb43aa7SLois Curfman McInnes } 1747eb43aa7SLois Curfman McInnes } 1757eb43aa7SLois Curfman McInnes return 0; 1767eb43aa7SLois Curfman McInnes } 1777eb43aa7SLois Curfman McInnes 17817ab2063SBarry Smith #include "draw.h" 17917ab2063SBarry Smith #include "pinclude/pviewer.h" 180416022c9SBarry Smith #include "sysio.h" 18117ab2063SBarry Smith 182416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 18317ab2063SBarry Smith { 184416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 185416022c9SBarry Smith int i, fd, *col_lens, ierr; 18617ab2063SBarry Smith 187416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1880452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 189416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 190416022c9SBarry Smith col_lens[1] = a->m; 191416022c9SBarry Smith col_lens[2] = a->n; 192416022c9SBarry Smith col_lens[3] = a->nz; 193416022c9SBarry Smith 194416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 195416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 196416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 19717ab2063SBarry Smith } 198416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 1990452661fSBarry Smith PetscFree(col_lens); 200416022c9SBarry Smith 201416022c9SBarry Smith /* store column indices (zero start index) */ 202416022c9SBarry Smith if (a->indexshift) { 203416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 20417ab2063SBarry Smith } 205416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 206416022c9SBarry Smith if (a->indexshift) { 207416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 20817ab2063SBarry Smith } 209416022c9SBarry Smith 210416022c9SBarry Smith /* store nonzero values */ 211416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 21217ab2063SBarry Smith return 0; 21317ab2063SBarry Smith } 214416022c9SBarry Smith 215416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 216416022c9SBarry Smith { 217416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 218416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 21917ab2063SBarry Smith FILE *fd; 22017ab2063SBarry Smith char *outputname; 22117ab2063SBarry Smith 222227d817aSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 223416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 224416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 22517ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 22608480c60SBarry Smith /* no need to print additional information */ ; 22717ab2063SBarry Smith } 22817ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 22917ab2063SBarry Smith int nz, nzalloc, mem; 230416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 231416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 23217ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 23317ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 23417ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 23517ab2063SBarry Smith 23617ab2063SBarry Smith for (i=0; i<m; i++) { 237416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 23817ab2063SBarry Smith #if defined(PETSC_COMPLEX) 2397a743949SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]), 240416022c9SBarry Smith imag(a->a[j])); 24117ab2063SBarry Smith #else 2427a743949SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 24317ab2063SBarry Smith #endif 24417ab2063SBarry Smith } 24517ab2063SBarry Smith } 24617ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 24717ab2063SBarry Smith } 24817ab2063SBarry Smith else { 24917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 25017ab2063SBarry Smith fprintf(fd,"row %d:",i); 251416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 25217ab2063SBarry Smith #if defined(PETSC_COMPLEX) 253416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 254416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 25517ab2063SBarry Smith } 25617ab2063SBarry Smith else { 257416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 25817ab2063SBarry Smith } 25917ab2063SBarry Smith #else 260416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 26117ab2063SBarry Smith #endif 26217ab2063SBarry Smith } 26317ab2063SBarry Smith fprintf(fd,"\n"); 26417ab2063SBarry Smith } 26517ab2063SBarry Smith } 26617ab2063SBarry Smith fflush(fd); 267416022c9SBarry Smith return 0; 268416022c9SBarry Smith } 269416022c9SBarry Smith 270416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 271416022c9SBarry Smith { 272416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 273cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 274cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 275d7e8b826SBarry Smith Draw draw = (Draw) viewer; 276cddf8d76SBarry Smith DrawButton button; 277cddf8d76SBarry Smith 278416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 279416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 280416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 281416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 282cddf8d76SBarry Smith color = DRAW_BLUE; 283416022c9SBarry Smith for ( i=0; i<m; i++ ) { 284cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 285416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 286cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 287cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 288cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 289cddf8d76SBarry Smith #else 290cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 291cddf8d76SBarry Smith #endif 292cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 293cddf8d76SBarry Smith } 294cddf8d76SBarry Smith } 295cddf8d76SBarry Smith color = DRAW_CYAN; 296cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 297cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 298cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 299cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 300cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 301cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 302cddf8d76SBarry Smith } 303cddf8d76SBarry Smith } 304cddf8d76SBarry Smith color = DRAW_RED; 305cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 306cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 307cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 308cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 309cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 310cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 311cddf8d76SBarry Smith #else 312cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 313cddf8d76SBarry Smith #endif 314cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 315416022c9SBarry Smith } 316416022c9SBarry Smith } 317416022c9SBarry Smith DrawFlush(draw); 318cddf8d76SBarry Smith DrawGetPause(draw,&pause); 319cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 320cddf8d76SBarry Smith 321cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 322cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 323cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 324cddf8d76SBarry Smith DrawClear(draw); 325cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 326cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 327cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 328cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 329cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 330cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 331cddf8d76SBarry Smith w *= scale; h *= scale; 332cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 333cddf8d76SBarry Smith color = DRAW_BLUE; 334cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 335cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 336cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 337cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 338cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 339cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 340cddf8d76SBarry Smith #else 341cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 342cddf8d76SBarry Smith #endif 343cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 344cddf8d76SBarry Smith } 345cddf8d76SBarry Smith } 346cddf8d76SBarry Smith color = DRAW_CYAN; 347cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 348cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 349cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 350cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 351cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 352cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 353cddf8d76SBarry Smith } 354cddf8d76SBarry Smith } 355cddf8d76SBarry Smith color = DRAW_RED; 356cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 357cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 358cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 359cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 360cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 361cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 362cddf8d76SBarry Smith #else 363cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 364cddf8d76SBarry Smith #endif 365cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 366cddf8d76SBarry Smith } 367cddf8d76SBarry Smith } 368cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 369cddf8d76SBarry Smith } 370416022c9SBarry Smith return 0; 371416022c9SBarry Smith } 372416022c9SBarry Smith 373416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 374416022c9SBarry Smith { 375416022c9SBarry Smith Mat A = (Mat) obj; 376416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 377416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 378416022c9SBarry Smith 379416022c9SBarry Smith if (!viewer) { 380416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 381416022c9SBarry Smith } 382416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 383416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 384416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 385416022c9SBarry Smith } 386416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 387416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 388416022c9SBarry Smith } 389416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 390416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 391416022c9SBarry Smith } 392416022c9SBarry Smith } 393416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 394416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 395416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 39617ab2063SBarry Smith } 39717ab2063SBarry Smith return 0; 39817ab2063SBarry Smith } 399c456f294SBarry Smith extern int Mat_AIJ_CheckInode(Mat); 400416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 40117ab2063SBarry Smith { 402416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 40341c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 404416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 405416022c9SBarry Smith Scalar *aa = a->a, *ap; 40617ab2063SBarry Smith 40717ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 40817ab2063SBarry Smith 40917ab2063SBarry Smith for ( i=1; i<m; i++ ) { 410416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 41117ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 41217ab2063SBarry Smith if (fshift) { 413416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 41417ab2063SBarry Smith N = ailen[i]; 41517ab2063SBarry Smith for ( j=0; j<N; j++ ) { 41617ab2063SBarry Smith ip[j-fshift] = ip[j]; 41717ab2063SBarry Smith ap[j-fshift] = ap[j]; 41817ab2063SBarry Smith } 41917ab2063SBarry Smith } 42017ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 42117ab2063SBarry Smith } 42217ab2063SBarry Smith if (m) { 42317ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 42417ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 42517ab2063SBarry Smith } 42617ab2063SBarry Smith /* reset ilen and imax for each row */ 42717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 42817ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 42917ab2063SBarry Smith } 430416022c9SBarry Smith a->nz = ai[m] + shift; 43117ab2063SBarry Smith 43217ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 433416022c9SBarry Smith if (fshift && a->diag) { 4340452661fSBarry Smith PetscFree(a->diag); 435416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 436416022c9SBarry Smith a->diag = 0; 43717ab2063SBarry Smith } 438b810aeb4SBarry Smith PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n", 439b810aeb4SBarry Smith fshift,a->nz,m); 440b810aeb4SBarry Smith PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n", 441b810aeb4SBarry Smith a->reallocs); 44276dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 44341c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 44417ab2063SBarry Smith return 0; 44517ab2063SBarry Smith } 44617ab2063SBarry Smith 447416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 44817ab2063SBarry Smith { 449416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 450cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 45117ab2063SBarry Smith return 0; 45217ab2063SBarry Smith } 453416022c9SBarry Smith 45417ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 45517ab2063SBarry Smith { 456416022c9SBarry Smith Mat A = (Mat) obj; 457416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 458d5d45c9bSBarry Smith 45917ab2063SBarry Smith #if defined(PETSC_LOG) 460416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 46117ab2063SBarry Smith #endif 4620452661fSBarry Smith PetscFree(a->a); 4630452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4640452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4650452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4660452661fSBarry Smith if (a->imax) PetscFree(a->imax); 4670452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 46876dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 4690452661fSBarry Smith PetscFree(a); 470416022c9SBarry Smith PLogObjectDestroy(A); 4710452661fSBarry Smith PetscHeaderDestroy(A); 47217ab2063SBarry Smith return 0; 47317ab2063SBarry Smith } 47417ab2063SBarry Smith 475416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 47617ab2063SBarry Smith { 47717ab2063SBarry Smith return 0; 47817ab2063SBarry Smith } 47917ab2063SBarry Smith 480416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 48117ab2063SBarry Smith { 482416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 483416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 4844b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 485416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 486416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 487416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 488e2f28af5SBarry Smith else if (op == ROWS_SORTED || 489e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 490e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 491e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 492df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 493df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4944d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 4956c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_1) a->inode.limit = 1; 4966c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_2) a->inode.limit = 2; 4976c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_3) a->inode.limit = 3; 4986c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_4) a->inode.limit = 4; 4996c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_5) a->inode.limit = 5; 500e2f28af5SBarry Smith else 5014d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 50217ab2063SBarry Smith return 0; 50317ab2063SBarry Smith } 50417ab2063SBarry Smith 505416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 50617ab2063SBarry Smith { 507416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 508416022c9SBarry Smith int i,j, n,shift = a->indexshift; 50917ab2063SBarry Smith Scalar *x, zero = 0.0; 51017ab2063SBarry Smith 51117ab2063SBarry Smith VecSet(&zero,v); 51217ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 513416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 514416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 515416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 516416022c9SBarry Smith if (a->j[j]+shift == i) { 517416022c9SBarry Smith x[i] = a->a[j]; 51817ab2063SBarry Smith break; 51917ab2063SBarry Smith } 52017ab2063SBarry Smith } 52117ab2063SBarry Smith } 52217ab2063SBarry Smith return 0; 52317ab2063SBarry Smith } 52417ab2063SBarry Smith 52517ab2063SBarry Smith /* -------------------------------------------------------*/ 52617ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 52717ab2063SBarry Smith /* -------------------------------------------------------*/ 528416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 52917ab2063SBarry Smith { 530416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 53117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 532416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 53317ab2063SBarry Smith 53417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 535cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 53617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 53717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 538416022c9SBarry Smith idx = a->j + a->i[i] + shift; 539416022c9SBarry Smith v = a->a + a->i[i] + shift; 540416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 54117ab2063SBarry Smith alpha = x[i]; 54217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 54317ab2063SBarry Smith } 544416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 54517ab2063SBarry Smith return 0; 54617ab2063SBarry Smith } 547d5d45c9bSBarry Smith 548416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 54917ab2063SBarry Smith { 550416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 55117ab2063SBarry Smith Scalar *x, *y, *v, alpha; 552416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 55317ab2063SBarry Smith 55417ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 55517ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 55617ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 55717ab2063SBarry Smith for ( i=0; i<m; i++ ) { 558416022c9SBarry Smith idx = a->j + a->i[i] + shift; 559416022c9SBarry Smith v = a->a + a->i[i] + shift; 560416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 56117ab2063SBarry Smith alpha = x[i]; 56217ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 56317ab2063SBarry Smith } 56417ab2063SBarry Smith return 0; 56517ab2063SBarry Smith } 56617ab2063SBarry Smith 567416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 56817ab2063SBarry Smith { 569416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 57017ab2063SBarry Smith Scalar *x, *y, *v, sum; 571416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 57217ab2063SBarry Smith 57317ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 57417ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 575416022c9SBarry Smith idx = a->j; 576416022c9SBarry Smith v = a->a; 577416022c9SBarry Smith ii = a->i; 57817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 579416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 58017ab2063SBarry Smith sum = 0.0; 58117ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 58217ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 583416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 58417ab2063SBarry Smith y[i] = sum; 58517ab2063SBarry Smith } 586416022c9SBarry Smith PLogFlops(2*a->nz - m); 58717ab2063SBarry Smith return 0; 58817ab2063SBarry Smith } 58917ab2063SBarry Smith 590416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 59117ab2063SBarry Smith { 592416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 59317ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 594cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 59517ab2063SBarry Smith 59617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 59717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 598cddf8d76SBarry Smith idx = a->j; 599cddf8d76SBarry Smith v = a->a; 600cddf8d76SBarry Smith ii = a->i; 60117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 602cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 60317ab2063SBarry Smith sum = y[i]; 604cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 605cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 60617ab2063SBarry Smith z[i] = sum; 60717ab2063SBarry Smith } 608416022c9SBarry Smith PLogFlops(2*a->nz); 60917ab2063SBarry Smith return 0; 61017ab2063SBarry Smith } 61117ab2063SBarry Smith 61217ab2063SBarry Smith /* 61317ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 61417ab2063SBarry Smith */ 61517ab2063SBarry Smith 616416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 61717ab2063SBarry Smith { 618416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 619416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 62017ab2063SBarry Smith 6210452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 622416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 623416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 624416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 625416022c9SBarry Smith if (a->j[j]+shift == i) { 62617ab2063SBarry Smith diag[i] = j - shift; 62717ab2063SBarry Smith break; 62817ab2063SBarry Smith } 62917ab2063SBarry Smith } 63017ab2063SBarry Smith } 631416022c9SBarry Smith a->diag = diag; 63217ab2063SBarry Smith return 0; 63317ab2063SBarry Smith } 63417ab2063SBarry Smith 635416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 63617ab2063SBarry Smith double fshift,int its,Vec xx) 63717ab2063SBarry Smith { 638416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 639416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 640d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 64117ab2063SBarry Smith 64217ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 643416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 644416022c9SBarry Smith diag = a->diag; 645416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 64617ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 64717ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 64817ab2063SBarry Smith bs = b + shift; 64917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 650416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 651416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 652416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 653416022c9SBarry Smith v = a->a + diag[i] + (!shift); 65417ab2063SBarry Smith sum = b[i]*d/omega; 65517ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 65617ab2063SBarry Smith x[i] = sum; 65717ab2063SBarry Smith } 65817ab2063SBarry Smith return 0; 65917ab2063SBarry Smith } 66017ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 66117ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 66217ab2063SBarry Smith } 663416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 66417ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 66517ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 66617ab2063SBarry Smith 66717ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 66817ab2063SBarry Smith 66917ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 67017ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 67117ab2063SBarry Smith is the relaxation factor. 67217ab2063SBarry Smith */ 6730452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 67417ab2063SBarry Smith scale = (2.0/omega) - 1.0; 67517ab2063SBarry Smith 67617ab2063SBarry Smith /* x = (E + U)^{-1} b */ 67717ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 678416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 679416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 680416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 681416022c9SBarry Smith v = a->a + diag[i] + (!shift); 68217ab2063SBarry Smith sum = b[i]; 68317ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 68417ab2063SBarry Smith x[i] = omega*(sum/d); 68517ab2063SBarry Smith } 68617ab2063SBarry Smith 68717ab2063SBarry Smith /* t = b - (2*E - D)x */ 688416022c9SBarry Smith v = a->a; 68917ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 69017ab2063SBarry Smith 69117ab2063SBarry Smith /* t = (E + L)^{-1}t */ 692416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 693416022c9SBarry Smith diag = a->diag; 69417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 695416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 696416022c9SBarry Smith n = diag[i] - a->i[i]; 697416022c9SBarry Smith idx = a->j + a->i[i] + shift; 698416022c9SBarry Smith v = a->a + a->i[i] + shift; 69917ab2063SBarry Smith sum = t[i]; 70017ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 70117ab2063SBarry Smith t[i] = omega*(sum/d); 70217ab2063SBarry Smith } 70317ab2063SBarry Smith 70417ab2063SBarry Smith /* x = x + t */ 70517ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7060452661fSBarry Smith PetscFree(t); 70717ab2063SBarry Smith return 0; 70817ab2063SBarry Smith } 70917ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 71017ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 71117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 712416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 713416022c9SBarry Smith n = diag[i] - a->i[i]; 714416022c9SBarry Smith idx = a->j + a->i[i] + shift; 715416022c9SBarry Smith v = a->a + a->i[i] + shift; 71617ab2063SBarry Smith sum = b[i]; 71717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 71817ab2063SBarry Smith x[i] = omega*(sum/d); 71917ab2063SBarry Smith } 72017ab2063SBarry Smith xb = x; 72117ab2063SBarry Smith } 72217ab2063SBarry Smith else xb = b; 72317ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 72417ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 72517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 726416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 72717ab2063SBarry Smith } 72817ab2063SBarry Smith } 72917ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 73017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 731416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 732416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 733416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 734416022c9SBarry Smith v = a->a + diag[i] + (!shift); 73517ab2063SBarry Smith sum = xb[i]; 73617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 73717ab2063SBarry Smith x[i] = omega*(sum/d); 73817ab2063SBarry Smith } 73917ab2063SBarry Smith } 74017ab2063SBarry Smith its--; 74117ab2063SBarry Smith } 74217ab2063SBarry Smith while (its--) { 74317ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 74417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 745416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 746416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 747416022c9SBarry Smith idx = a->j + a->i[i] + shift; 748416022c9SBarry Smith v = a->a + a->i[i] + shift; 74917ab2063SBarry Smith sum = b[i]; 75017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 75117ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 75217ab2063SBarry Smith } 75317ab2063SBarry Smith } 75417ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 75517ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 756416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 757416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 758416022c9SBarry Smith idx = a->j + a->i[i] + shift; 759416022c9SBarry Smith v = a->a + a->i[i] + shift; 76017ab2063SBarry Smith sum = b[i]; 76117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 76217ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 76317ab2063SBarry Smith } 76417ab2063SBarry Smith } 76517ab2063SBarry Smith } 76617ab2063SBarry Smith return 0; 76717ab2063SBarry Smith } 76817ab2063SBarry Smith 769d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 77017ab2063SBarry Smith { 771416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 772416022c9SBarry Smith *nz = a->nz; 773416022c9SBarry Smith *nzalloc = a->maxnz; 774416022c9SBarry Smith *mem = (int)A->mem; 77517ab2063SBarry Smith return 0; 77617ab2063SBarry Smith } 77717ab2063SBarry Smith 77817ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 77917ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 78017ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 78117ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 78217ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78317ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 78417ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78517ab2063SBarry Smith 78617ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 78717ab2063SBarry Smith { 788416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 789416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 79017ab2063SBarry Smith 79117ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 79217ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 79317ab2063SBarry Smith if (diag) { 79417ab2063SBarry Smith for ( i=0; i<N; i++ ) { 795416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 796416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 797416022c9SBarry Smith a->ilen[rows[i]] = 1; 798416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 799416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 80017ab2063SBarry Smith } 80117ab2063SBarry Smith else { 80217ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 80317ab2063SBarry Smith CHKERRQ(ierr); 80417ab2063SBarry Smith } 80517ab2063SBarry Smith } 80617ab2063SBarry Smith } 80717ab2063SBarry Smith else { 80817ab2063SBarry Smith for ( i=0; i<N; i++ ) { 809416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 810416022c9SBarry Smith a->ilen[rows[i]] = 0; 81117ab2063SBarry Smith } 81217ab2063SBarry Smith } 81317ab2063SBarry Smith ISRestoreIndices(is,&rows); 81417ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 81517ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 81617ab2063SBarry Smith return 0; 81717ab2063SBarry Smith } 81817ab2063SBarry Smith 819416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 82017ab2063SBarry Smith { 821416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 822416022c9SBarry Smith *m = a->m; *n = a->n; 82317ab2063SBarry Smith return 0; 82417ab2063SBarry Smith } 82517ab2063SBarry Smith 826416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 82717ab2063SBarry Smith { 828416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 829416022c9SBarry Smith *m = 0; *n = a->m; 83017ab2063SBarry Smith return 0; 83117ab2063SBarry Smith } 832416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 83317ab2063SBarry Smith { 834416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 835c456f294SBarry Smith int *itmp,i,shift = a->indexshift; 83617ab2063SBarry Smith 837416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 83817ab2063SBarry Smith 839416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 840416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 84117ab2063SBarry Smith if (idx) { 84217ab2063SBarry Smith if (*nz) { 843416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8440452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 84517ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 84617ab2063SBarry Smith } 84717ab2063SBarry Smith else *idx = 0; 84817ab2063SBarry Smith } 84917ab2063SBarry Smith return 0; 85017ab2063SBarry Smith } 85117ab2063SBarry Smith 852416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 85317ab2063SBarry Smith { 8540452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 85517ab2063SBarry Smith return 0; 85617ab2063SBarry Smith } 85717ab2063SBarry Smith 858cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 85917ab2063SBarry Smith { 860416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 861416022c9SBarry Smith Scalar *v = a->a; 86217ab2063SBarry Smith double sum = 0.0; 863416022c9SBarry Smith int i, j,shift = a->indexshift; 86417ab2063SBarry Smith 86517ab2063SBarry Smith if (type == NORM_FROBENIUS) { 866416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 86717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 86817ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 86917ab2063SBarry Smith #else 87017ab2063SBarry Smith sum += (*v)*(*v); v++; 87117ab2063SBarry Smith #endif 87217ab2063SBarry Smith } 87317ab2063SBarry Smith *norm = sqrt(sum); 87417ab2063SBarry Smith } 87517ab2063SBarry Smith else if (type == NORM_1) { 87617ab2063SBarry Smith double *tmp; 877416022c9SBarry Smith int *jj = a->j; 8780452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 879cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 88017ab2063SBarry Smith *norm = 0.0; 881416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 882a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 88317ab2063SBarry Smith } 884416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 88517ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 88617ab2063SBarry Smith } 8870452661fSBarry Smith PetscFree(tmp); 88817ab2063SBarry Smith } 88917ab2063SBarry Smith else if (type == NORM_INFINITY) { 89017ab2063SBarry Smith *norm = 0.0; 891416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 892416022c9SBarry Smith v = a->a + a->i[j] + shift; 89317ab2063SBarry Smith sum = 0.0; 894416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 895cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 89617ab2063SBarry Smith } 89717ab2063SBarry Smith if (sum > *norm) *norm = sum; 89817ab2063SBarry Smith } 89917ab2063SBarry Smith } 90017ab2063SBarry Smith else { 90148d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 90217ab2063SBarry Smith } 90317ab2063SBarry Smith return 0; 90417ab2063SBarry Smith } 90517ab2063SBarry Smith 906416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 90717ab2063SBarry Smith { 908416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 909416022c9SBarry Smith Mat C; 910416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 911416022c9SBarry Smith int shift = a->indexshift; 912d5d45c9bSBarry Smith Scalar *array = a->a; 91317ab2063SBarry Smith 9143638b69dSLois Curfman McInnes if (B == PETSC_NULL && m != a->n) 9153638b69dSLois Curfman McInnes SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place"); 9160452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 917cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 91817ab2063SBarry Smith if (shift) { 91917ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 92017ab2063SBarry Smith } 92117ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 922416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 9230452661fSBarry Smith PetscFree(col); 92417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 92517ab2063SBarry Smith len = ai[i+1]-ai[i]; 926416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 92717ab2063SBarry Smith array += len; aj += len; 92817ab2063SBarry Smith } 92917ab2063SBarry Smith if (shift) { 93017ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 93117ab2063SBarry Smith } 93217ab2063SBarry Smith 933416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 934416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 93517ab2063SBarry Smith 9363638b69dSLois Curfman McInnes if (B != PETSC_NULL) { 937416022c9SBarry Smith *B = C; 93817ab2063SBarry Smith } else { 939416022c9SBarry Smith /* This isn't really an in-place transpose */ 9400452661fSBarry Smith PetscFree(a->a); 9410452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9420452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9430452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9440452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9450452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9460452661fSBarry Smith PetscFree(a); 947416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9480452661fSBarry Smith PetscHeaderDestroy(C); 94917ab2063SBarry Smith } 95017ab2063SBarry Smith return 0; 95117ab2063SBarry Smith } 95217ab2063SBarry Smith 953f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 95417ab2063SBarry Smith { 955416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 95617ab2063SBarry Smith Scalar *l,*r,x,*v; 957d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 95817ab2063SBarry Smith 95917ab2063SBarry Smith if (ll) { 96017ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 961f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 962416022c9SBarry Smith v = a->a; 96317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 96417ab2063SBarry Smith x = l[i]; 965416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 96617ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 96717ab2063SBarry Smith } 96817ab2063SBarry Smith } 96917ab2063SBarry Smith if (rr) { 97017ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 971f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 972416022c9SBarry Smith v = a->a; jj = a->j; 97317ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 97417ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 97517ab2063SBarry Smith } 97617ab2063SBarry Smith } 97717ab2063SBarry Smith return 0; 97817ab2063SBarry Smith } 97917ab2063SBarry Smith 980cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 98117ab2063SBarry Smith { 982db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 98302834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 984a2744918SBarry Smith register int sum,lensi; 98502834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 986a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 987db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 988416022c9SBarry Smith Mat C; 98917ab2063SBarry Smith 99017ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 99117ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 99217ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 993*7264ac53SSatish Balay ierr = SYIsort(nrows, irow); CHKERRQ(ierr); 99417ab2063SBarry Smith 995*7264ac53SSatish Balay if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 99602834360SBarry Smith /* special case of contiguous rows */ 9970452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 99802834360SBarry Smith starts = lens + ncols; 99902834360SBarry Smith /* loop over new rows determining lens and starting points */ 100002834360SBarry Smith for (i=0; i<nrows; i++) { 1001a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1002a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 100302834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1004d8ced48eSBarry Smith if (aj[k]+shift >= first) { 100502834360SBarry Smith starts[i] = k; 100602834360SBarry Smith break; 100702834360SBarry Smith } 100802834360SBarry Smith } 1009a2744918SBarry Smith sum = 0; 101002834360SBarry Smith while (k < kend) { 1011d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1012a2744918SBarry Smith sum++; 101302834360SBarry Smith } 1014a2744918SBarry Smith lens[i] = sum; 101502834360SBarry Smith } 101602834360SBarry Smith /* create submatrix */ 1017cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 101808480c60SBarry Smith int n_cols,n_rows; 101908480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 102008480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1021d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 102208480c60SBarry Smith C = *B; 102308480c60SBarry Smith } 102408480c60SBarry Smith else { 102502834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 102608480c60SBarry Smith } 1027db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1028db02288aSLois Curfman McInnes 102902834360SBarry Smith /* loop over rows inserting into submatrix */ 1030db02288aSLois Curfman McInnes a_new = c->a; 1031db02288aSLois Curfman McInnes j_new = c->j; 1032db02288aSLois Curfman McInnes i_new = c->i; 1033db02288aSLois Curfman McInnes i_new[0] = -shift; 103402834360SBarry Smith for (i=0; i<nrows; i++) { 1035a2744918SBarry Smith ii = starts[i]; 1036a2744918SBarry Smith lensi = lens[i]; 1037a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1038a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 103902834360SBarry Smith } 1040a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1041a2744918SBarry Smith a_new += lensi; 1042a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1043a2744918SBarry Smith c->ilen[i] = lensi; 104402834360SBarry Smith } 10450452661fSBarry Smith PetscFree(lens); 104602834360SBarry Smith } 104702834360SBarry Smith else { 104802834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 1049*7264ac53SSatish Balay ierr = SYIsort(ncols,icol); CHKERRQ(ierr); 10500452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 105102834360SBarry Smith ssmap = smap + shift; 10527eb43aa7SLois Curfman McInnes cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 105302834360SBarry Smith lens = cwork + ncols; 10540452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1055cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 105617ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 105702834360SBarry Smith /* determine lens of each row */ 105802834360SBarry Smith for (i=0; i<nrows; i++) { 1059d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 106002834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 106102834360SBarry Smith lens[i] = 0; 106202834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1063d8ced48eSBarry Smith if (ssmap[aj[k]]) { 106402834360SBarry Smith lens[i]++; 106502834360SBarry Smith } 106602834360SBarry Smith } 106702834360SBarry Smith } 106817ab2063SBarry Smith /* Create and fill new matrix */ 1069a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 107008480c60SBarry Smith int n_cols,n_rows; 107108480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 107208480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1073d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 107408480c60SBarry Smith C = *B; 107508480c60SBarry Smith } 107608480c60SBarry Smith else { 107702834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 107808480c60SBarry Smith } 107917ab2063SBarry Smith for (i=0; i<nrows; i++) { 108017ab2063SBarry Smith nznew = 0; 1081d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 1082416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 108317ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 108402834360SBarry Smith if (ssmap[a->j[k]]) { 108502834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1086416022c9SBarry Smith vwork[nznew++] = a->a[k]; 108717ab2063SBarry Smith } 108817ab2063SBarry Smith } 1089416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 109017ab2063SBarry Smith } 109102834360SBarry Smith /* Free work space */ 109202834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 10930452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 109402834360SBarry Smith } 1095416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1096416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 109717ab2063SBarry Smith 109817ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1099416022c9SBarry Smith *B = C; 110017ab2063SBarry Smith return 0; 110117ab2063SBarry Smith } 110217ab2063SBarry Smith 1103a871dcd8SBarry Smith /* 110463b91edcSBarry Smith note: This can only work for identity for row and col. It would 110563b91edcSBarry Smith be good to check this and otherwise generate an error. 1106a871dcd8SBarry Smith */ 110763b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1108a871dcd8SBarry Smith { 110963b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 111008480c60SBarry Smith int ierr; 111163b91edcSBarry Smith Mat outA; 111263b91edcSBarry Smith 1113a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1114a871dcd8SBarry Smith 111563b91edcSBarry Smith outA = inA; 111663b91edcSBarry Smith inA->factor = FACTOR_LU; 111763b91edcSBarry Smith a->row = row; 111863b91edcSBarry Smith a->col = col; 111963b91edcSBarry Smith 11200452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 112163b91edcSBarry Smith 112208480c60SBarry Smith if (!a->diag) { 112308480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 112463b91edcSBarry Smith } 112563b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1126a871dcd8SBarry Smith return 0; 1127a871dcd8SBarry Smith } 1128a871dcd8SBarry Smith 1129f0b747eeSBarry Smith #include "pinclude/plapack.h" 1130f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1131f0b747eeSBarry Smith { 1132f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1133f0b747eeSBarry Smith int one = 1; 1134f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1135f0b747eeSBarry Smith PLogFlops(a->nz); 1136f0b747eeSBarry Smith return 0; 1137f0b747eeSBarry Smith } 1138f0b747eeSBarry Smith 1139cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1140cddf8d76SBarry Smith Mat **B) 1141cddf8d76SBarry Smith { 1142cddf8d76SBarry Smith int ierr,i; 1143cddf8d76SBarry Smith 1144cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 11450452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1146cddf8d76SBarry Smith } 1147cddf8d76SBarry Smith 1148cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1149cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1150cddf8d76SBarry Smith } 1151cddf8d76SBarry Smith return 0; 1152cddf8d76SBarry Smith } 1153cddf8d76SBarry Smith 1154e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 11554dcbc457SBarry Smith { 1156e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 115706763907SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 11588a047759SSatish Balay int start, end, *ai, *aj; 115906763907SSatish Balay char *table; 11608a047759SSatish Balay shift = a->indexshift; 1161e4d965acSSatish Balay m = a->m; 1162e4d965acSSatish Balay ai = a->i; 11638a047759SSatish Balay aj = a->j+shift; 11648a047759SSatish Balay 1165e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 116606763907SSatish Balay 116706763907SSatish Balay table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 116806763907SSatish Balay nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 116906763907SSatish Balay 1170e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1171e4d965acSSatish Balay /* Initialise the two local arrays */ 1172e4d965acSSatish Balay isz = 0; 117306763907SSatish Balay PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1174e4d965acSSatish Balay 1175e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 11764dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1177e4d965acSSatish Balay ierr = ISGetLocalSize(is[i],&n); CHKERRQ(ierr); 1178e4d965acSSatish Balay 1179e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1180e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 118106763907SSatish Balay if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 11824dcbc457SBarry Smith } 118306763907SSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 118406763907SSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1185e4d965acSSatish Balay 118604a348a9SBarry Smith k = 0; 118704a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 118804a348a9SBarry Smith n = isz; 118906763907SSatish Balay for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1190e4d965acSSatish Balay row = nidx[k]; 1191e4d965acSSatish Balay start = ai[row]; 1192e4d965acSSatish Balay end = ai[row+1]; 119304a348a9SBarry Smith for ( l = start; l<end ; l++){ 11948a047759SSatish Balay val = aj[l] + shift; 119506763907SSatish Balay if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1196e4d965acSSatish Balay } 1197e4d965acSSatish Balay } 1198e4d965acSSatish Balay } 1199e4d965acSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1200e4d965acSSatish Balay } 120104a348a9SBarry Smith PetscFree(table); 120206763907SSatish Balay PetscFree(nidx); 1203e4d965acSSatish Balay return 0; 12044dcbc457SBarry Smith } 120517ab2063SBarry Smith 1206682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1207682d7d0cSBarry Smith { 1208682d7d0cSBarry Smith static int called = 0; 1209682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1210682d7d0cSBarry Smith 1211682d7d0cSBarry Smith if (called) return 0; else called = 1; 12122a7368beSLois Curfman McInnes MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 1213682d7d0cSBarry Smith MPIU_printf(comm," -mat_lu_pivotthreshold <threshold>\n"); 12142a7368beSLois Curfman McInnes MPIU_printf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 1215682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_no_inode - Do not use inodes\n"); 1216682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1217682d7d0cSBarry Smith #if defined(HAVE_ESSL) 1218682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1219682d7d0cSBarry Smith #endif 1220682d7d0cSBarry Smith return 0; 1221682d7d0cSBarry Smith } 1222*7264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg); 1223682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 122417ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 122517ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1226416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1227416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 122817ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 122917ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 123017ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 123117ab2063SBarry Smith MatRelax_SeqAIJ, 123217ab2063SBarry Smith MatTranspose_SeqAIJ, 1233*7264ac53SSatish Balay MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1234f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 123517ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 123617ab2063SBarry Smith MatCompress_SeqAIJ, 123717ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 123817ab2063SBarry Smith MatGetReordering_SeqAIJ, 123917ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 124017ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 124117ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 124217ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1243416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 12443d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1245cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 12467eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1247682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1248f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1249f0b747eeSBarry Smith MatScale_SeqAIJ}; 125017ab2063SBarry Smith 125117ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 125217ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 125317ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 125417ab2063SBarry Smith 125517ab2063SBarry Smith /*@C 1256682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 12570d15e28bSLois Curfman McInnes (the default parallel PETSc format). For good matrix assembly performance 12586e62573dSLois Curfman McInnes the user should preallocate the matrix storage by setting the parameter nz 12590d15e28bSLois Curfman McInnes (or nzz). By setting these parameters accurately, performance can be 12600d15e28bSLois Curfman McInnes increased by more than a factor of 50. 126117ab2063SBarry Smith 126217ab2063SBarry Smith Input Parameters: 126317ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 126417ab2063SBarry Smith . m - number of rows 126517ab2063SBarry Smith . n - number of columns 126617ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 126717ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 126817ab2063SBarry Smith 126917ab2063SBarry Smith Output Parameter: 1270416022c9SBarry Smith . A - the matrix 127117ab2063SBarry Smith 127217ab2063SBarry Smith Notes: 127317ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 127417ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 12750002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 12760002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 127717ab2063SBarry Smith 127817ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1279a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 12800d15e28bSLois Curfman McInnes allocation. For additional details, see the users manual chapter on 12810d15e28bSLois Curfman McInnes matrices and the file $(PETSC_DIR)/Performance. 128217ab2063SBarry Smith 1283682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1284682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1285682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 12866c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 12876c7ebb05SLois Curfman McInnes 12886c7ebb05SLois Curfman McInnes Options Database Keys: 12896c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 12906e62573dSLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit. 12916e62573dSLois Curfman McInnes $ (max limit=5) 12926e62573dSLois Curfman McInnes $ -mat_aij_oneindex - Internally use indexing starting at 1 12936e62573dSLois Curfman McInnes $ rather than 0. Note: When calling MatSetValues(), 12946e62573dSLois Curfman McInnes $ the user still MUST index entries starting at 0! 129517ab2063SBarry Smith 129617ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 129717ab2063SBarry Smith @*/ 1298416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 129917ab2063SBarry Smith { 1300416022c9SBarry Smith Mat B; 1301416022c9SBarry Smith Mat_SeqAIJ *b; 130269957df2SSatish Balay int i,len,ierr, flg; 1303d5d45c9bSBarry Smith 1304416022c9SBarry Smith *A = 0; 13050452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1306416022c9SBarry Smith PLogObjectCreate(B); 13070452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1308416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1309416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1310416022c9SBarry Smith B->view = MatView_SeqAIJ; 1311416022c9SBarry Smith B->factor = 0; 1312416022c9SBarry Smith B->lupivotthreshold = 1.0; 13137a743949SBarry Smith ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 131469957df2SSatish Balay &flg); CHKERRQ(ierr); 13157a743949SBarry Smith b->ilu_preserve_row_sums = PETSC_FALSE; 13167a743949SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 13177a743949SBarry Smith (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1318416022c9SBarry Smith b->row = 0; 1319416022c9SBarry Smith b->col = 0; 1320416022c9SBarry Smith b->indexshift = 0; 1321b810aeb4SBarry Smith b->reallocs = 0; 132269957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 132369957df2SSatish Balay if (flg) b->indexshift = -1; 132417ab2063SBarry Smith 1325416022c9SBarry Smith b->m = m; 1326416022c9SBarry Smith b->n = n; 13270452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1328b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 13297b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 13307b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1331416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 133217ab2063SBarry Smith nz = nz*m; 133317ab2063SBarry Smith } 133417ab2063SBarry Smith else { 133517ab2063SBarry Smith nz = 0; 1336416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 133717ab2063SBarry Smith } 133817ab2063SBarry Smith 133917ab2063SBarry Smith /* allocate the matrix space */ 1340416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 13410452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1342416022c9SBarry Smith b->j = (int *) (b->a + nz); 1343cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1344416022c9SBarry Smith b->i = b->j + nz; 1345416022c9SBarry Smith b->singlemalloc = 1; 134617ab2063SBarry Smith 1347416022c9SBarry Smith b->i[0] = -b->indexshift; 134817ab2063SBarry Smith for (i=1; i<m+1; i++) { 1349416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 135017ab2063SBarry Smith } 135117ab2063SBarry Smith 1352416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 13530452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1354416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1355416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 135617ab2063SBarry Smith 1357416022c9SBarry Smith b->nz = 0; 1358416022c9SBarry Smith b->maxnz = nz; 1359416022c9SBarry Smith b->sorted = 0; 1360416022c9SBarry Smith b->roworiented = 1; 1361416022c9SBarry Smith b->nonew = 0; 1362416022c9SBarry Smith b->diag = 0; 1363416022c9SBarry Smith b->solve_work = 0; 136471bd300dSLois Curfman McInnes b->spptr = 0; 1365754ec7b1SSatish Balay b->inode.node_count = 0; 1366754ec7b1SSatish Balay b->inode.size = 0; 13676c7ebb05SLois Curfman McInnes b->inode.limit = 5; 13686c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 136917ab2063SBarry Smith 1370416022c9SBarry Smith *A = B; 13714b14c69eSBarry Smith /* SuperLU is not currently supported through PETSc */ 13724b14c69eSBarry Smith #if defined(HAVE_SUPERLU) 137369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 137469957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 13754b14c69eSBarry Smith #endif 137669957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 137769957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 137869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 137969957df2SSatish Balay if (flg) { 1380416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1381416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 138217ab2063SBarry Smith } 138369957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 138469957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 138517ab2063SBarry Smith return 0; 138617ab2063SBarry Smith } 138717ab2063SBarry Smith 13883d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 138917ab2063SBarry Smith { 1390416022c9SBarry Smith Mat C; 1391416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 139208480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 139317ab2063SBarry Smith 13944043dd9cSLois Curfman McInnes *B = 0; 13950452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1396416022c9SBarry Smith PLogObjectCreate(C); 13970452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 139841c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1399416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1400416022c9SBarry Smith C->view = MatView_SeqAIJ; 1401416022c9SBarry Smith C->factor = A->factor; 1402416022c9SBarry Smith c->row = 0; 1403416022c9SBarry Smith c->col = 0; 1404416022c9SBarry Smith c->indexshift = shift; 1405c456f294SBarry Smith C->assembled = PETSC_TRUE; 140617ab2063SBarry Smith 1407416022c9SBarry Smith c->m = a->m; 1408416022c9SBarry Smith c->n = a->n; 140917ab2063SBarry Smith 14100452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 14110452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 141217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1413416022c9SBarry Smith c->imax[i] = a->imax[i]; 1414416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 141517ab2063SBarry Smith } 141617ab2063SBarry Smith 141717ab2063SBarry Smith /* allocate the matrix space */ 1418416022c9SBarry Smith c->singlemalloc = 1; 1419416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 14200452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1421416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1422416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1423416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 142417ab2063SBarry Smith if (m > 0) { 1425416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 142608480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1427416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 142817ab2063SBarry Smith } 142908480c60SBarry Smith } 143017ab2063SBarry Smith 1431416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1432416022c9SBarry Smith c->sorted = a->sorted; 1433416022c9SBarry Smith c->roworiented = a->roworiented; 1434416022c9SBarry Smith c->nonew = a->nonew; 14357a743949SBarry Smith c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 143617ab2063SBarry Smith 1437416022c9SBarry Smith if (a->diag) { 14380452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1439416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 144017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1441416022c9SBarry Smith c->diag[i] = a->diag[i]; 144217ab2063SBarry Smith } 144317ab2063SBarry Smith } 1444416022c9SBarry Smith else c->diag = 0; 14456c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 14466c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1447754ec7b1SSatish Balay if (a->inode.size){ 1448754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1449754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1450754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1451754ec7b1SSatish Balay } else { 1452754ec7b1SSatish Balay c->inode.size = 0; 1453754ec7b1SSatish Balay c->inode.node_count = 0; 1454754ec7b1SSatish Balay } 1455416022c9SBarry Smith c->nz = a->nz; 1456416022c9SBarry Smith c->maxnz = a->maxnz; 1457416022c9SBarry Smith c->solve_work = 0; 145876dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1459754ec7b1SSatish Balay 1460416022c9SBarry Smith *B = C; 146117ab2063SBarry Smith return 0; 146217ab2063SBarry Smith } 146317ab2063SBarry Smith 1464416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 146517ab2063SBarry Smith { 1466416022c9SBarry Smith Mat_SeqAIJ *a; 1467416022c9SBarry Smith Mat B; 146817699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 146917ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 147017ab2063SBarry Smith MPI_Comm comm = vobj->comm; 147117ab2063SBarry Smith 147217699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 147317699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 147417ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1475416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 147648d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 147717ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 147817ab2063SBarry Smith 147917ab2063SBarry Smith /* read in row lengths */ 14800452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1481416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 148217ab2063SBarry Smith 148317ab2063SBarry Smith /* create our matrix */ 1484416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1485416022c9SBarry Smith B = *A; 1486416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1487416022c9SBarry Smith shift = a->indexshift; 148817ab2063SBarry Smith 148917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1490416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 149117ab2063SBarry Smith if (shift) { 149217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1493416022c9SBarry Smith a->j[i] += 1; 149417ab2063SBarry Smith } 149517ab2063SBarry Smith } 149617ab2063SBarry Smith 149717ab2063SBarry Smith /* read in nonzero values */ 1498416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 149917ab2063SBarry Smith 150017ab2063SBarry Smith /* set matrix "i" values */ 1501416022c9SBarry Smith a->i[0] = -shift; 150217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1503416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1504416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 150517ab2063SBarry Smith } 15060452661fSBarry Smith PetscFree(rowlengths); 150717ab2063SBarry Smith 1508416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1509416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 151017ab2063SBarry Smith return 0; 151117ab2063SBarry Smith } 151217ab2063SBarry Smith 151317ab2063SBarry Smith 151417ab2063SBarry Smith 1515*7264ac53SSatish Balay /* 1516*7264ac53SSatish Balay 1517*7264ac53SSatish Balay flg =0 if error; 1518*7264ac53SSatish Balay flg =1 if both are equal; 1519*7264ac53SSatish Balay */ 1520*7264ac53SSatish Balay int MatEqual_SeqAIJ(Mat A,Mat B, int* flg) 1521*7264ac53SSatish Balay { 1522*7264ac53SSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 1523*7264ac53SSatish Balay 1524*7264ac53SSatish Balay if (B->type !=MATSEQAIJ) SETERRQ(1,"MatEqual_SeqAIJ:Both matrices should be of type MATSEQAIJ"); 1525*7264ac53SSatish Balay 1526*7264ac53SSatish Balay /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 1527*7264ac53SSatish Balay if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1528*7264ac53SSatish Balay (a->indexshift != b->indexshift)) { *flg =0 ; return 0; } 1529*7264ac53SSatish Balay 1530*7264ac53SSatish Balay /* if the a->i are the same */ 1531*7264ac53SSatish Balay if(PetscMemcmp((char *)a->i, (char*)b->i, (a->n+1)*sizeof(int))) { *flg =0 ; return 0;} 1532*7264ac53SSatish Balay 1533*7264ac53SSatish Balay /* if a->j are the same */ 1534*7264ac53SSatish Balay if(PetscMemcmp((char *)a->j, (char*)b->j, (a->nz)*sizeof(int))) { 1535*7264ac53SSatish Balay *flg =0 ; return 0; 1536*7264ac53SSatish Balay } 1537*7264ac53SSatish Balay 1538*7264ac53SSatish Balay /* if a->j are the same */ 1539*7264ac53SSatish Balay if(PetscMemcmp((char *)a->a, (char*)b->a, (a->nz)*sizeof(Scalar))) { 1540*7264ac53SSatish Balay *flg =0 ; return 0; 1541*7264ac53SSatish Balay } 1542*7264ac53SSatish Balay *flg =1 ; 1543*7264ac53SSatish Balay return 0; 1544*7264ac53SSatish Balay 1545*7264ac53SSatish Balay } 1546