117ab2063SBarry Smith #ifndef lint 2*7eb43aa7SLois Curfman McInnes static char vcid[] = "$Id: aij.c,v 1.120 1995/11/21 21:07:55 balay Exp curfman $"; 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" 1217ab2063SBarry Smith 1317ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**); 1417ab2063SBarry Smith 15416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1617ab2063SBarry Smith { 17416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 18a2744918SBarry Smith int ierr, *ia, *ja,n,*idx,i; 193d54f168SSatish Balay /*Viewer V1, V2;*/ 2017ab2063SBarry Smith 21416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix"); 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 42416022c9SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr); 43416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 4444f6d32bSSatish Balay /* ISView(*rperm, STDOUT_VIEWER_SELF);*/ 45f31bba2fSSatish Balay 463d54f168SSatish Balay /*ViewerFileOpenASCII(MPI_COMM_SELF,"row_is_orig", &V1); 47f31bba2fSSatish Balay ViewerFileOpenASCII(MPI_COMM_SELF,"col_is_orig", &V2); 48f31bba2fSSatish Balay ISView(*rperm,V1); 49f31bba2fSSatish Balay ISView(*cperm,V2); 50f31bba2fSSatish Balay ViewerDestroy(V1); 513d54f168SSatish Balay ViewerDestroy(V2);*/ 52f31bba2fSSatish Balay 530452661fSBarry Smith PetscFree(ia); PetscFree(ja); 5417ab2063SBarry Smith return 0; 5517ab2063SBarry Smith } 5617ab2063SBarry Smith 57416022c9SBarry Smith #define CHUNKSIZE 10 5817ab2063SBarry Smith 5917ab2063SBarry Smith /* This version has row oriented v */ 60416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 6117ab2063SBarry Smith { 62416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 63416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 64416022c9SBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen; 65d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 66416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 6717ab2063SBarry Smith 6817ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 69416022c9SBarry Smith row = im[k]; 7017ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 71416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 7217ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 7317ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 74416022c9SBarry Smith low = 0; 7517ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 76416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 77416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 78416022c9SBarry Smith col = in[l] - shift; value = *v++; 79416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 80416022c9SBarry Smith while (high-low > 5) { 81416022c9SBarry Smith t = (low+high)/2; 82416022c9SBarry Smith if (rp[t] > col) high = t; 83416022c9SBarry Smith else low = t; 8417ab2063SBarry Smith } 85416022c9SBarry Smith for ( i=low; i<high; i++ ) { 8617ab2063SBarry Smith if (rp[i] > col) break; 8717ab2063SBarry Smith if (rp[i] == col) { 88416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 8917ab2063SBarry Smith else ap[i] = value; 9017ab2063SBarry Smith goto noinsert; 9117ab2063SBarry Smith } 9217ab2063SBarry Smith } 9317ab2063SBarry Smith if (nonew) goto noinsert; 9417ab2063SBarry Smith if (nrow >= rmax) { 9517ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 96416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 9717ab2063SBarry Smith Scalar *new_a; 9817ab2063SBarry Smith 9917ab2063SBarry Smith /* malloc new storage space */ 100416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1010452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 10217ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 10317ab2063SBarry Smith new_i = new_j + new_nz; 10417ab2063SBarry Smith 10517ab2063SBarry Smith /* copy over old data into new slots */ 10617ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 107416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 108416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 109416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 110416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 11117ab2063SBarry Smith len*sizeof(int)); 112416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 113416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 11417ab2063SBarry Smith len*sizeof(Scalar)); 11517ab2063SBarry Smith /* free up old matrix storage */ 1160452661fSBarry Smith PetscFree(a->a); 1170452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 118416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 119416022c9SBarry Smith a->singlemalloc = 1; 12017ab2063SBarry Smith 12117ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 122416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 123416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 124416022c9SBarry Smith a->maxnz += CHUNKSIZE; 12517ab2063SBarry Smith } 126416022c9SBarry Smith N = nrow++ - 1; a->nz++; 127416022c9SBarry Smith /* shift up all the later entries in this row */ 128416022c9SBarry Smith for ( ii=N; ii>=i; ii-- ) { 12917ab2063SBarry Smith rp[ii+1] = rp[ii]; 13017ab2063SBarry Smith ap[ii+1] = ap[ii]; 13117ab2063SBarry Smith } 13217ab2063SBarry Smith rp[i] = col; 13317ab2063SBarry Smith ap[i] = value; 13417ab2063SBarry Smith noinsert:; 135416022c9SBarry Smith low = i + 1; 13617ab2063SBarry Smith } 13717ab2063SBarry Smith ailen[row] = nrow; 13817ab2063SBarry Smith } 13917ab2063SBarry Smith return 0; 14017ab2063SBarry Smith } 14117ab2063SBarry Smith 142*7eb43aa7SLois Curfman McInnes static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 143*7eb43aa7SLois Curfman McInnes { 144*7eb43aa7SLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 145*7eb43aa7SLois Curfman McInnes int *rp, k, low, high, t, row, nrow, i, col, l, ict, *aj = a->j; 146*7eb43aa7SLois Curfman McInnes int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 147*7eb43aa7SLois Curfman McInnes Scalar *ap, *aa = a->a, zero = 0.0; 148*7eb43aa7SLois Curfman McInnes 149*7eb43aa7SLois Curfman McInnes if (!a->assembled) SETERRQ(1,"MatGetValues_SeqAIJ:Not for unassembled matrix"); 150*7eb43aa7SLois Curfman McInnes ict = 0; 151*7eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 152*7eb43aa7SLois Curfman McInnes row = im[k]; 153*7eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 154*7eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 155*7eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 156*7eb43aa7SLois Curfman McInnes nrow = ailen[row]; 157*7eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 158*7eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 159*7eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 160*7eb43aa7SLois Curfman McInnes col = in[l] - shift; 161*7eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 162*7eb43aa7SLois Curfman McInnes while (high-low > 5) { 163*7eb43aa7SLois Curfman McInnes t = (low+high)/2; 164*7eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 165*7eb43aa7SLois Curfman McInnes else low = t; 166*7eb43aa7SLois Curfman McInnes } 167*7eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 168*7eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 169*7eb43aa7SLois Curfman McInnes if (rp[i] == col) { 170*7eb43aa7SLois Curfman McInnes v[ict] = ap[i]; 171*7eb43aa7SLois Curfman McInnes goto finished; 172*7eb43aa7SLois Curfman McInnes } 173*7eb43aa7SLois Curfman McInnes } 174*7eb43aa7SLois Curfman McInnes v[ict] = zero; 175*7eb43aa7SLois Curfman McInnes finished:; 176*7eb43aa7SLois Curfman McInnes ict++; 177*7eb43aa7SLois Curfman McInnes } 178*7eb43aa7SLois Curfman McInnes } 179*7eb43aa7SLois Curfman McInnes return 0; 180*7eb43aa7SLois Curfman McInnes } 181*7eb43aa7SLois Curfman McInnes 18217ab2063SBarry Smith #include "draw.h" 18317ab2063SBarry Smith #include "pinclude/pviewer.h" 184416022c9SBarry Smith #include "sysio.h" 18517ab2063SBarry Smith 186416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 18717ab2063SBarry Smith { 188416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 189416022c9SBarry Smith int i, fd, *col_lens, ierr; 19017ab2063SBarry Smith 191416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1920452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 193416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 194416022c9SBarry Smith col_lens[1] = a->m; 195416022c9SBarry Smith col_lens[2] = a->n; 196416022c9SBarry Smith col_lens[3] = a->nz; 197416022c9SBarry Smith 198416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 199416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 200416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 20117ab2063SBarry Smith } 202416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 2030452661fSBarry Smith PetscFree(col_lens); 204416022c9SBarry Smith 205416022c9SBarry Smith /* store column indices (zero start index) */ 206416022c9SBarry Smith if (a->indexshift) { 207416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 20817ab2063SBarry Smith } 209416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 210416022c9SBarry Smith if (a->indexshift) { 211416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 21217ab2063SBarry Smith } 213416022c9SBarry Smith 214416022c9SBarry Smith /* store nonzero values */ 215416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 21617ab2063SBarry Smith return 0; 21717ab2063SBarry Smith } 218416022c9SBarry Smith 219416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 220416022c9SBarry Smith { 221416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 222416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 22317ab2063SBarry Smith FILE *fd; 22417ab2063SBarry Smith char *outputname; 22517ab2063SBarry Smith 226416022c9SBarry Smith ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 227416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 228416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 22917ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 23008480c60SBarry Smith /* no need to print additional information */ ; 23117ab2063SBarry Smith } 23217ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 23317ab2063SBarry Smith int nz, nzalloc, mem; 234416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 235416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 23617ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 23717ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 23817ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 23917ab2063SBarry Smith 24017ab2063SBarry Smith for (i=0; i<m; i++) { 241416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 24217ab2063SBarry Smith #if defined(PETSC_COMPLEX) 243416022c9SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j],real(a->a[j]), 244416022c9SBarry Smith imag(a->a[j])); 24517ab2063SBarry Smith #else 246416022c9SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j], a->a[j]); 24717ab2063SBarry Smith #endif 24817ab2063SBarry Smith } 24917ab2063SBarry Smith } 25017ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 25117ab2063SBarry Smith } 25217ab2063SBarry Smith else { 25317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 25417ab2063SBarry Smith fprintf(fd,"row %d:",i); 255416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 25617ab2063SBarry Smith #if defined(PETSC_COMPLEX) 257416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 258416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 25917ab2063SBarry Smith } 26017ab2063SBarry Smith else { 261416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 26217ab2063SBarry Smith } 26317ab2063SBarry Smith #else 264416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 26517ab2063SBarry Smith #endif 26617ab2063SBarry Smith } 26717ab2063SBarry Smith fprintf(fd,"\n"); 26817ab2063SBarry Smith } 26917ab2063SBarry Smith } 27017ab2063SBarry Smith fflush(fd); 271416022c9SBarry Smith return 0; 272416022c9SBarry Smith } 273416022c9SBarry Smith 274416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 275416022c9SBarry Smith { 276416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 277cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 278cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 279d7e8b826SBarry Smith Draw draw = (Draw) viewer; 280cddf8d76SBarry Smith DrawButton button; 281cddf8d76SBarry Smith 282416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 283416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 284416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 285416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 286cddf8d76SBarry Smith color = DRAW_BLUE; 287416022c9SBarry Smith for ( i=0; i<m; i++ ) { 288cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 289416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 290cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 291cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 292cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 293cddf8d76SBarry Smith #else 294cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 295cddf8d76SBarry Smith #endif 296cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 297cddf8d76SBarry Smith } 298cddf8d76SBarry Smith } 299cddf8d76SBarry Smith color = DRAW_CYAN; 300cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 301cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 302cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 303cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 304cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 305cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 306cddf8d76SBarry Smith } 307cddf8d76SBarry Smith } 308cddf8d76SBarry Smith color = DRAW_RED; 309cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 310cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 311cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 312cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 313cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 314cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 315cddf8d76SBarry Smith #else 316cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 317cddf8d76SBarry Smith #endif 318cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 319416022c9SBarry Smith } 320416022c9SBarry Smith } 321416022c9SBarry Smith DrawFlush(draw); 322cddf8d76SBarry Smith DrawGetPause(draw,&pause); 323cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 324cddf8d76SBarry Smith 325cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 326cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 327cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 328cddf8d76SBarry Smith DrawClear(draw); 329cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 330cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 331cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 332cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 333cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 334cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 335cddf8d76SBarry Smith w *= scale; h *= scale; 336cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 337cddf8d76SBarry Smith color = DRAW_BLUE; 338cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 339cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 340cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 341cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 342cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 343cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 344cddf8d76SBarry Smith #else 345cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 346cddf8d76SBarry Smith #endif 347cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 348cddf8d76SBarry Smith } 349cddf8d76SBarry Smith } 350cddf8d76SBarry Smith color = DRAW_CYAN; 351cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 352cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 353cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 354cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 355cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 356cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 357cddf8d76SBarry Smith } 358cddf8d76SBarry Smith } 359cddf8d76SBarry Smith color = DRAW_RED; 360cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 361cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 362cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 363cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 364cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 365cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 366cddf8d76SBarry Smith #else 367cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 368cddf8d76SBarry Smith #endif 369cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 370cddf8d76SBarry Smith } 371cddf8d76SBarry Smith } 372cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 373cddf8d76SBarry Smith } 374416022c9SBarry Smith return 0; 375416022c9SBarry Smith } 376416022c9SBarry Smith 377416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 378416022c9SBarry Smith { 379416022c9SBarry Smith Mat A = (Mat) obj; 380416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 381416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 382416022c9SBarry Smith 383416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix"); 384416022c9SBarry Smith if (!viewer) { 385416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 386416022c9SBarry Smith } 387416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 388416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 389416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 390416022c9SBarry Smith } 391416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 392416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 393416022c9SBarry Smith } 394416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 395416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 396416022c9SBarry Smith } 397416022c9SBarry Smith } 398416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 399416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 400416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 40117ab2063SBarry Smith } 40217ab2063SBarry Smith return 0; 40317ab2063SBarry Smith } 40441c01911SSatish Balay int Mat_AIJ_CheckInode(Mat); 405416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 40617ab2063SBarry Smith { 407416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 40841c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 409416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 410416022c9SBarry Smith Scalar *aa = a->a, *ap; 41117ab2063SBarry Smith 41217ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 41317ab2063SBarry Smith 41417ab2063SBarry Smith for ( i=1; i<m; i++ ) { 415416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 41617ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 41717ab2063SBarry Smith if (fshift) { 418416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 41917ab2063SBarry Smith N = ailen[i]; 42017ab2063SBarry Smith for ( j=0; j<N; j++ ) { 42117ab2063SBarry Smith ip[j-fshift] = ip[j]; 42217ab2063SBarry Smith ap[j-fshift] = ap[j]; 42317ab2063SBarry Smith } 42417ab2063SBarry Smith } 42517ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 42617ab2063SBarry Smith } 42717ab2063SBarry Smith if (m) { 42817ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 42917ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 43017ab2063SBarry Smith } 43117ab2063SBarry Smith /* reset ilen and imax for each row */ 43217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 43317ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 43417ab2063SBarry Smith } 435416022c9SBarry Smith a->nz = ai[m] + shift; 43617ab2063SBarry Smith 43717ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 438416022c9SBarry Smith if (fshift && a->diag) { 4390452661fSBarry Smith PetscFree(a->diag); 440416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 441416022c9SBarry Smith a->diag = 0; 44217ab2063SBarry Smith } 44376dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 44441c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 445416022c9SBarry Smith a->assembled = 1; 44617ab2063SBarry Smith return 0; 44717ab2063SBarry Smith } 44817ab2063SBarry Smith 449416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 45017ab2063SBarry Smith { 451416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 452cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 45317ab2063SBarry Smith return 0; 45417ab2063SBarry Smith } 455416022c9SBarry Smith 45617ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 45717ab2063SBarry Smith { 458416022c9SBarry Smith Mat A = (Mat) obj; 459416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 460d5d45c9bSBarry Smith 46117ab2063SBarry Smith #if defined(PETSC_LOG) 462416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 46317ab2063SBarry Smith #endif 4640452661fSBarry Smith PetscFree(a->a); 4650452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4660452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4670452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4680452661fSBarry Smith if (a->imax) PetscFree(a->imax); 4690452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 47076dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 4710452661fSBarry Smith PetscFree(a); 472416022c9SBarry Smith PLogObjectDestroy(A); 4730452661fSBarry Smith PetscHeaderDestroy(A); 47417ab2063SBarry Smith return 0; 47517ab2063SBarry Smith } 47617ab2063SBarry Smith 477416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 47817ab2063SBarry Smith { 47917ab2063SBarry Smith return 0; 48017ab2063SBarry Smith } 48117ab2063SBarry Smith 482416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 48317ab2063SBarry Smith { 484416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 485416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 486416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 487416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 488416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 489e2f28af5SBarry Smith else if (op == ROWS_SORTED || 490e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 491e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 492e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 493df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 494df719601SLois Curfman McInnes else if (op == COLUMN_ORIENTED) 4954d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:COLUMN_ORIENTED");} 496df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4974d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 498e2f28af5SBarry Smith else 4994d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 50017ab2063SBarry Smith return 0; 50117ab2063SBarry Smith } 50217ab2063SBarry Smith 503416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 50417ab2063SBarry Smith { 505416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 506416022c9SBarry Smith int i,j, n,shift = a->indexshift; 50717ab2063SBarry Smith Scalar *x, zero = 0.0; 50817ab2063SBarry Smith 509416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix"); 51017ab2063SBarry Smith VecSet(&zero,v); 51117ab2063SBarry Smith VecGetArray(v,&x); VecGetLocalSize(v,&n); 512416022c9SBarry Smith if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 513416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 514416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 515416022c9SBarry Smith if (a->j[j]+shift == i) { 516416022c9SBarry Smith x[i] = a->a[j]; 51717ab2063SBarry Smith break; 51817ab2063SBarry Smith } 51917ab2063SBarry Smith } 52017ab2063SBarry Smith } 52117ab2063SBarry Smith return 0; 52217ab2063SBarry Smith } 52317ab2063SBarry Smith 52417ab2063SBarry Smith /* -------------------------------------------------------*/ 52517ab2063SBarry Smith /* Should check that shapes of vectors and matrices match */ 52617ab2063SBarry Smith /* -------------------------------------------------------*/ 527416022c9SBarry Smith static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 52817ab2063SBarry Smith { 529416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 53017ab2063SBarry Smith Scalar *x, *y, *v, alpha; 531416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift; 53217ab2063SBarry Smith 533416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix"); 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 554416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix"); 55517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 55617ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 55717ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 55817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 559416022c9SBarry Smith idx = a->j + a->i[i] + shift; 560416022c9SBarry Smith v = a->a + a->i[i] + shift; 561416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 56217ab2063SBarry Smith alpha = x[i]; 56317ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 56417ab2063SBarry Smith } 56517ab2063SBarry Smith return 0; 56617ab2063SBarry Smith } 56717ab2063SBarry Smith 568416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 56917ab2063SBarry Smith { 570416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 57117ab2063SBarry Smith Scalar *x, *y, *v, sum; 572416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 57317ab2063SBarry Smith 574416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix"); 57517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 57617ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 577416022c9SBarry Smith idx = a->j; 578416022c9SBarry Smith v = a->a; 579416022c9SBarry Smith ii = a->i; 58017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 581416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 58217ab2063SBarry Smith sum = 0.0; 58317ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 58417ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 585416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 58617ab2063SBarry Smith y[i] = sum; 58717ab2063SBarry Smith } 588416022c9SBarry Smith PLogFlops(2*a->nz - m); 58917ab2063SBarry Smith return 0; 59017ab2063SBarry Smith } 59117ab2063SBarry Smith 592416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 59317ab2063SBarry Smith { 594416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 59517ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 596cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 59717ab2063SBarry Smith 59848d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix"); 59917ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 60017ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 601cddf8d76SBarry Smith idx = a->j; 602cddf8d76SBarry Smith v = a->a; 603cddf8d76SBarry Smith ii = a->i; 60417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 605cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 60617ab2063SBarry Smith sum = y[i]; 607cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 608cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 60917ab2063SBarry Smith z[i] = sum; 61017ab2063SBarry Smith } 611416022c9SBarry Smith PLogFlops(2*a->nz); 61217ab2063SBarry Smith return 0; 61317ab2063SBarry Smith } 61417ab2063SBarry Smith 61517ab2063SBarry Smith /* 61617ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 61717ab2063SBarry Smith */ 61817ab2063SBarry Smith 619416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 62017ab2063SBarry Smith { 621416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 622416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 62317ab2063SBarry Smith 624416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:unassembled matrix"); 6250452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 626416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 627416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 628416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 629416022c9SBarry Smith if (a->j[j]+shift == i) { 63017ab2063SBarry Smith diag[i] = j - shift; 63117ab2063SBarry Smith break; 63217ab2063SBarry Smith } 63317ab2063SBarry Smith } 63417ab2063SBarry Smith } 635416022c9SBarry Smith a->diag = diag; 63617ab2063SBarry Smith return 0; 63717ab2063SBarry Smith } 63817ab2063SBarry Smith 639416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 64017ab2063SBarry Smith double fshift,int its,Vec xx) 64117ab2063SBarry Smith { 642416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 643416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 644d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 64517ab2063SBarry Smith 64617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 647416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 648416022c9SBarry Smith diag = a->diag; 649416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 65017ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 65117ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 65217ab2063SBarry Smith bs = b + shift; 65317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 654416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 655416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 656416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 657416022c9SBarry Smith v = a->a + diag[i] + (!shift); 65817ab2063SBarry Smith sum = b[i]*d/omega; 65917ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 66017ab2063SBarry Smith x[i] = sum; 66117ab2063SBarry Smith } 66217ab2063SBarry Smith return 0; 66317ab2063SBarry Smith } 66417ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 66517ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 66617ab2063SBarry Smith } 667416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 66817ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 66917ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 67017ab2063SBarry Smith 67117ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 67217ab2063SBarry Smith 67317ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 67417ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 67517ab2063SBarry Smith is the relaxation factor. 67617ab2063SBarry Smith */ 6770452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 67817ab2063SBarry Smith scale = (2.0/omega) - 1.0; 67917ab2063SBarry Smith 68017ab2063SBarry Smith /* x = (E + U)^{-1} b */ 68117ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 682416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 683416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 684416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 685416022c9SBarry Smith v = a->a + diag[i] + (!shift); 68617ab2063SBarry Smith sum = b[i]; 68717ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 68817ab2063SBarry Smith x[i] = omega*(sum/d); 68917ab2063SBarry Smith } 69017ab2063SBarry Smith 69117ab2063SBarry Smith /* t = b - (2*E - D)x */ 692416022c9SBarry Smith v = a->a; 69317ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 69417ab2063SBarry Smith 69517ab2063SBarry Smith /* t = (E + L)^{-1}t */ 696416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 697416022c9SBarry Smith diag = a->diag; 69817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 699416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 700416022c9SBarry Smith n = diag[i] - a->i[i]; 701416022c9SBarry Smith idx = a->j + a->i[i] + shift; 702416022c9SBarry Smith v = a->a + a->i[i] + shift; 70317ab2063SBarry Smith sum = t[i]; 70417ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 70517ab2063SBarry Smith t[i] = omega*(sum/d); 70617ab2063SBarry Smith } 70717ab2063SBarry Smith 70817ab2063SBarry Smith /* x = x + t */ 70917ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7100452661fSBarry Smith PetscFree(t); 71117ab2063SBarry Smith return 0; 71217ab2063SBarry Smith } 71317ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 71417ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 71517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 716416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 717416022c9SBarry Smith n = diag[i] - a->i[i]; 718416022c9SBarry Smith idx = a->j + a->i[i] + shift; 719416022c9SBarry Smith v = a->a + a->i[i] + shift; 72017ab2063SBarry Smith sum = b[i]; 72117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 72217ab2063SBarry Smith x[i] = omega*(sum/d); 72317ab2063SBarry Smith } 72417ab2063SBarry Smith xb = x; 72517ab2063SBarry Smith } 72617ab2063SBarry Smith else xb = b; 72717ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 72817ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 72917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 730416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 73117ab2063SBarry Smith } 73217ab2063SBarry Smith } 73317ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 73417ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 735416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 736416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 737416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 738416022c9SBarry Smith v = a->a + diag[i] + (!shift); 73917ab2063SBarry Smith sum = xb[i]; 74017ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 74117ab2063SBarry Smith x[i] = omega*(sum/d); 74217ab2063SBarry Smith } 74317ab2063SBarry Smith } 74417ab2063SBarry Smith its--; 74517ab2063SBarry Smith } 74617ab2063SBarry Smith while (its--) { 74717ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 74817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 749416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 750416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 751416022c9SBarry Smith idx = a->j + a->i[i] + shift; 752416022c9SBarry Smith v = a->a + a->i[i] + shift; 75317ab2063SBarry Smith sum = b[i]; 75417ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 75517ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 75617ab2063SBarry Smith } 75717ab2063SBarry Smith } 75817ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 75917ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 760416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 761416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 762416022c9SBarry Smith idx = a->j + a->i[i] + shift; 763416022c9SBarry Smith v = a->a + a->i[i] + shift; 76417ab2063SBarry Smith sum = b[i]; 76517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 76617ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 76717ab2063SBarry Smith } 76817ab2063SBarry Smith } 76917ab2063SBarry Smith } 77017ab2063SBarry Smith return 0; 77117ab2063SBarry Smith } 77217ab2063SBarry Smith 773d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 77417ab2063SBarry Smith { 775416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 776416022c9SBarry Smith *nz = a->nz; 777416022c9SBarry Smith *nzalloc = a->maxnz; 778416022c9SBarry Smith *mem = (int)A->mem; 77917ab2063SBarry Smith return 0; 78017ab2063SBarry Smith } 78117ab2063SBarry Smith 78217ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 78317ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 78417ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 78517ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 78617ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78717ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 78817ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78917ab2063SBarry Smith 79017ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 79117ab2063SBarry Smith { 792416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 793416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 79417ab2063SBarry Smith 79517ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 79617ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 79717ab2063SBarry Smith if (diag) { 79817ab2063SBarry Smith for ( i=0; i<N; i++ ) { 799416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 800416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 801416022c9SBarry Smith a->ilen[rows[i]] = 1; 802416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 803416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 80417ab2063SBarry Smith } 80517ab2063SBarry Smith else { 80617ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 80717ab2063SBarry Smith CHKERRQ(ierr); 80817ab2063SBarry Smith } 80917ab2063SBarry Smith } 81017ab2063SBarry Smith } 81117ab2063SBarry Smith else { 81217ab2063SBarry Smith for ( i=0; i<N; i++ ) { 813416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 814416022c9SBarry Smith a->ilen[rows[i]] = 0; 81517ab2063SBarry Smith } 81617ab2063SBarry Smith } 81717ab2063SBarry Smith ISRestoreIndices(is,&rows); 81817ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 81917ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 82017ab2063SBarry Smith return 0; 82117ab2063SBarry Smith } 82217ab2063SBarry Smith 823416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 82417ab2063SBarry Smith { 825416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 826416022c9SBarry Smith *m = a->m; *n = a->n; 82717ab2063SBarry Smith return 0; 82817ab2063SBarry Smith } 82917ab2063SBarry Smith 830416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 83117ab2063SBarry Smith { 832416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 833416022c9SBarry Smith *m = 0; *n = a->m; 83417ab2063SBarry Smith return 0; 83517ab2063SBarry Smith } 836416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 83717ab2063SBarry Smith { 838416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 839416022c9SBarry Smith int *itmp,i,ierr,shift = a->indexshift; 84017ab2063SBarry Smith 841416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 84217ab2063SBarry Smith 843416022c9SBarry Smith if (!a->assembled) { 844416022c9SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 845416022c9SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 84617ab2063SBarry Smith } 847416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 848416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 84917ab2063SBarry Smith if (idx) { 85017ab2063SBarry Smith if (*nz) { 851416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8520452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 85317ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 85417ab2063SBarry Smith } 85517ab2063SBarry Smith else *idx = 0; 85617ab2063SBarry Smith } 85717ab2063SBarry Smith return 0; 85817ab2063SBarry Smith } 85917ab2063SBarry Smith 860416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 86117ab2063SBarry Smith { 8620452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 86317ab2063SBarry Smith return 0; 86417ab2063SBarry Smith } 86517ab2063SBarry Smith 866cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 86717ab2063SBarry Smith { 868416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 869416022c9SBarry Smith Scalar *v = a->a; 87017ab2063SBarry Smith double sum = 0.0; 871416022c9SBarry Smith int i, j,shift = a->indexshift; 87217ab2063SBarry Smith 873416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix"); 87417ab2063SBarry Smith if (type == NORM_FROBENIUS) { 875416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 87617ab2063SBarry Smith #if defined(PETSC_COMPLEX) 87717ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 87817ab2063SBarry Smith #else 87917ab2063SBarry Smith sum += (*v)*(*v); v++; 88017ab2063SBarry Smith #endif 88117ab2063SBarry Smith } 88217ab2063SBarry Smith *norm = sqrt(sum); 88317ab2063SBarry Smith } 88417ab2063SBarry Smith else if (type == NORM_1) { 88517ab2063SBarry Smith double *tmp; 886416022c9SBarry Smith int *jj = a->j; 8870452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 888cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 88917ab2063SBarry Smith *norm = 0.0; 890416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 891a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 89217ab2063SBarry Smith } 893416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 89417ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 89517ab2063SBarry Smith } 8960452661fSBarry Smith PetscFree(tmp); 89717ab2063SBarry Smith } 89817ab2063SBarry Smith else if (type == NORM_INFINITY) { 89917ab2063SBarry Smith *norm = 0.0; 900416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 901416022c9SBarry Smith v = a->a + a->i[j] + shift; 90217ab2063SBarry Smith sum = 0.0; 903416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 904cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 90517ab2063SBarry Smith } 90617ab2063SBarry Smith if (sum > *norm) *norm = sum; 90717ab2063SBarry Smith } 90817ab2063SBarry Smith } 90917ab2063SBarry Smith else { 91048d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 91117ab2063SBarry Smith } 91217ab2063SBarry Smith return 0; 91317ab2063SBarry Smith } 91417ab2063SBarry Smith 915416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 91617ab2063SBarry Smith { 917416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 918416022c9SBarry Smith Mat C; 919416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 920416022c9SBarry Smith int shift = a->indexshift; 921d5d45c9bSBarry Smith Scalar *array = a->a; 92217ab2063SBarry Smith 923416022c9SBarry Smith if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 9240452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 925cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 92617ab2063SBarry Smith if (shift) { 92717ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 92817ab2063SBarry Smith } 92917ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 930416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 9310452661fSBarry Smith PetscFree(col); 93217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 93317ab2063SBarry Smith len = ai[i+1]-ai[i]; 934416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 93517ab2063SBarry Smith array += len; aj += len; 93617ab2063SBarry Smith } 93717ab2063SBarry Smith if (shift) { 93817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 93917ab2063SBarry Smith } 94017ab2063SBarry Smith 941416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 942416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 94317ab2063SBarry Smith 944416022c9SBarry Smith if (B) { 945416022c9SBarry Smith *B = C; 94617ab2063SBarry Smith } else { 947416022c9SBarry Smith /* This isn't really an in-place transpose */ 9480452661fSBarry Smith PetscFree(a->a); 9490452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9500452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9510452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9520452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9530452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9540452661fSBarry Smith PetscFree(a); 955416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9560452661fSBarry Smith PetscHeaderDestroy(C); 95717ab2063SBarry Smith } 95817ab2063SBarry Smith return 0; 95917ab2063SBarry Smith } 96017ab2063SBarry Smith 961416022c9SBarry Smith static int MatScale_SeqAIJ(Mat A,Vec ll,Vec rr) 96217ab2063SBarry Smith { 963416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 96417ab2063SBarry Smith Scalar *l,*r,x,*v; 965d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 96617ab2063SBarry Smith 96748d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix"); 96817ab2063SBarry Smith if (ll) { 96917ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 970416022c9SBarry Smith if (m != a->m) SETERRQ(1,"MatScale_SeqAIJ:Left scaling vector wrong length"); 971416022c9SBarry Smith v = a->a; 97217ab2063SBarry Smith for ( i=0; i<m; i++ ) { 97317ab2063SBarry Smith x = l[i]; 974416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 97517ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 97617ab2063SBarry Smith } 97717ab2063SBarry Smith } 97817ab2063SBarry Smith if (rr) { 97917ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 980416022c9SBarry Smith if (n != a->n) SETERRQ(1,"MatScale_SeqAIJ:Right scaling vector wrong length"); 981416022c9SBarry Smith v = a->a; jj = a->j; 98217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 98317ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 98417ab2063SBarry Smith } 98517ab2063SBarry Smith } 98617ab2063SBarry Smith return 0; 98717ab2063SBarry Smith } 98817ab2063SBarry Smith 989cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 99017ab2063SBarry Smith { 991db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 99202834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 993a2744918SBarry Smith register int sum,lensi; 99402834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 995a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 996db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 997416022c9SBarry Smith Mat C; 99817ab2063SBarry Smith 999416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix"); 100017ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 100117ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 100217ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 100317ab2063SBarry Smith 100402834360SBarry Smith if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 100502834360SBarry Smith /* special case of contiguous rows */ 10060452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 100702834360SBarry Smith starts = lens + ncols; 100802834360SBarry Smith /* loop over new rows determining lens and starting points */ 100902834360SBarry Smith for (i=0; i<nrows; i++) { 1010a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1011a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 101202834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1013a2744918SBarry Smith if (aj[k] >= first) { 101402834360SBarry Smith starts[i] = k; 101502834360SBarry Smith break; 101602834360SBarry Smith } 101702834360SBarry Smith } 1018a2744918SBarry Smith sum = 0; 101902834360SBarry Smith while (k < kend) { 1020a2744918SBarry Smith if (aj[k++] >= first+ncols) break; 1021a2744918SBarry Smith sum++; 102202834360SBarry Smith } 1023a2744918SBarry Smith lens[i] = sum; 102402834360SBarry Smith } 102502834360SBarry Smith /* create submatrix */ 1026cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 102708480c60SBarry Smith int n_cols,n_rows; 102808480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 102908480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 103008480c60SBarry Smith MatZeroEntries(*B); 103108480c60SBarry Smith C = *B; 103208480c60SBarry Smith } 103308480c60SBarry Smith else { 103402834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 103508480c60SBarry Smith } 1036db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1037db02288aSLois Curfman McInnes 103802834360SBarry Smith /* loop over rows inserting into submatrix */ 1039db02288aSLois Curfman McInnes a_new = c->a; 1040db02288aSLois Curfman McInnes j_new = c->j; 1041db02288aSLois Curfman McInnes i_new = c->i; 1042db02288aSLois Curfman McInnes i_new[0] = -shift; 104302834360SBarry Smith for (i=0; i<nrows; i++) { 1044a2744918SBarry Smith ii = starts[i]; 1045a2744918SBarry Smith lensi = lens[i]; 1046a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1047a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 104802834360SBarry Smith } 1049a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1050a2744918SBarry Smith a_new += lensi; 1051a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1052a2744918SBarry Smith c->ilen[i] = lensi; 105302834360SBarry Smith } 10540452661fSBarry Smith PetscFree(lens); 105502834360SBarry Smith } 105602834360SBarry Smith else { 105702834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 10580452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 105902834360SBarry Smith ssmap = smap + shift; 1060*7eb43aa7SLois Curfman McInnes cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 106102834360SBarry Smith lens = cwork + ncols; 10620452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1063cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 106417ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 106502834360SBarry Smith /* determine lens of each row */ 106602834360SBarry Smith for (i=0; i<nrows; i++) { 106702834360SBarry Smith kstart = a->i[irow[i]]+shift; 106802834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 106902834360SBarry Smith lens[i] = 0; 107002834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 107102834360SBarry Smith if (ssmap[a->j[k]]) { 107202834360SBarry Smith lens[i]++; 107302834360SBarry Smith } 107402834360SBarry Smith } 107502834360SBarry Smith } 107617ab2063SBarry Smith /* Create and fill new matrix */ 1077a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 107808480c60SBarry Smith int n_cols,n_rows; 107908480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 108008480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 108108480c60SBarry Smith MatZeroEntries(*B); 108208480c60SBarry Smith C = *B; 108308480c60SBarry Smith } 108408480c60SBarry Smith else { 108502834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 108608480c60SBarry Smith } 108717ab2063SBarry Smith for (i=0; i<nrows; i++) { 108817ab2063SBarry Smith nznew = 0; 1089416022c9SBarry Smith kstart = a->i[irow[i]]+shift; 1090416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 109117ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 109202834360SBarry Smith if (ssmap[a->j[k]]) { 109302834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1094416022c9SBarry Smith vwork[nznew++] = a->a[k]; 109517ab2063SBarry Smith } 109617ab2063SBarry Smith } 1097416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 109817ab2063SBarry Smith } 109902834360SBarry Smith /* Free work space */ 110002834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 11010452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 110202834360SBarry Smith } 1103416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1104416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 110517ab2063SBarry Smith 110617ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1107416022c9SBarry Smith *B = C; 110817ab2063SBarry Smith return 0; 110917ab2063SBarry Smith } 111017ab2063SBarry Smith 1111a871dcd8SBarry Smith /* 111263b91edcSBarry Smith note: This can only work for identity for row and col. It would 111363b91edcSBarry Smith be good to check this and otherwise generate an error. 1114a871dcd8SBarry Smith */ 111563b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1116a871dcd8SBarry Smith { 111763b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 111808480c60SBarry Smith int ierr; 111963b91edcSBarry Smith Mat outA; 112063b91edcSBarry Smith 1121a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1122a871dcd8SBarry Smith 112363b91edcSBarry Smith outA = inA; 112463b91edcSBarry Smith inA->factor = FACTOR_LU; 112563b91edcSBarry Smith a->row = row; 112663b91edcSBarry Smith a->col = col; 112763b91edcSBarry Smith 11280452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 112963b91edcSBarry Smith 113008480c60SBarry Smith if (!a->diag) { 113108480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 113263b91edcSBarry Smith } 113363b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1134a871dcd8SBarry Smith return 0; 1135a871dcd8SBarry Smith } 1136a871dcd8SBarry Smith 1137cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1138cddf8d76SBarry Smith Mat **B) 1139cddf8d76SBarry Smith { 1140cddf8d76SBarry Smith int ierr,i; 1141cddf8d76SBarry Smith 1142cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 11430452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1144cddf8d76SBarry Smith } 1145cddf8d76SBarry Smith 1146cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1147cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1148cddf8d76SBarry Smith } 1149cddf8d76SBarry Smith return 0; 1150cddf8d76SBarry Smith } 1151cddf8d76SBarry Smith 11524dcbc457SBarry Smith static int MatIncreaseOverlap_SeqAIJ(Mat A, int n, IS *is, int ov) 11534dcbc457SBarry Smith { 11544dcbc457SBarry Smith int i,m,*idx,ierr; 11554dcbc457SBarry Smith 11564dcbc457SBarry Smith for ( i=0; i<n; i++ ) { 11574dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 11584dcbc457SBarry Smith ISGetLocalSize(is[i],&m); 11594dcbc457SBarry Smith } 11604dcbc457SBarry Smith SETERRQ(1,"MatIncreaseOverlap_SeqAIJ:Not implemented"); 11614dcbc457SBarry Smith } 116217ab2063SBarry Smith /* -------------------------------------------------------------------*/ 116317ab2063SBarry Smith 116417ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 116517ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1166416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1167416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 116817ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 116917ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 117017ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 117117ab2063SBarry Smith MatRelax_SeqAIJ, 117217ab2063SBarry Smith MatTranspose_SeqAIJ, 117317ab2063SBarry Smith MatGetInfo_SeqAIJ,0, 117417ab2063SBarry Smith MatGetDiagonal_SeqAIJ,MatScale_SeqAIJ,MatNorm_SeqAIJ, 117517ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 117617ab2063SBarry Smith MatCompress_SeqAIJ, 117717ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 117817ab2063SBarry Smith MatGetReordering_SeqAIJ, 117917ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 118017ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 118117ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 118217ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1183416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 1184a871dcd8SBarry Smith MatCopyPrivate_SeqAIJ,0,0, 1185cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 1186*7eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1187*7eb43aa7SLois Curfman McInnes MatGetValues_SeqAIJ}; 118817ab2063SBarry Smith 118917ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 119017ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 119117ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 119217ab2063SBarry Smith 119317ab2063SBarry Smith /*@C 119417ab2063SBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ format 119517ab2063SBarry Smith (the default uniprocessor PETSc format). 119617ab2063SBarry Smith 119717ab2063SBarry Smith Input Parameters: 119817ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 119917ab2063SBarry Smith . m - number of rows 120017ab2063SBarry Smith . n - number of columns 120117ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 120217ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 120317ab2063SBarry Smith 120417ab2063SBarry Smith Output Parameter: 1205416022c9SBarry Smith . A - the matrix 120617ab2063SBarry Smith 120717ab2063SBarry Smith Notes: 120817ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 120917ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 12100002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 12110002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 121217ab2063SBarry Smith 121317ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 121417ab2063SBarry Smith Set both nz and nnz to zero for PETSc to control dynamic memory 121517ab2063SBarry Smith allocation. 121617ab2063SBarry Smith 121717ab2063SBarry Smith .keywords: matrix, aij, compressed row, sparse 121817ab2063SBarry Smith 121917ab2063SBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 122017ab2063SBarry Smith @*/ 1221416022c9SBarry Smith int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 122217ab2063SBarry Smith { 1223416022c9SBarry Smith Mat B; 1224416022c9SBarry Smith Mat_SeqAIJ *b; 122517ab2063SBarry Smith int i,len,ierr; 1226d5d45c9bSBarry Smith 1227416022c9SBarry Smith *A = 0; 12280452661fSBarry Smith PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1229416022c9SBarry Smith PLogObjectCreate(B); 12300452661fSBarry Smith B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1231416022c9SBarry Smith PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1232416022c9SBarry Smith B->destroy = MatDestroy_SeqAIJ; 1233416022c9SBarry Smith B->view = MatView_SeqAIJ; 1234416022c9SBarry Smith B->factor = 0; 1235416022c9SBarry Smith B->lupivotthreshold = 1.0; 1236416022c9SBarry Smith OptionsGetDouble(0,"-mat_lu_pivotthreshold",&B->lupivotthreshold); 1237416022c9SBarry Smith b->row = 0; 1238416022c9SBarry Smith b->col = 0; 1239416022c9SBarry Smith b->indexshift = 0; 1240416022c9SBarry Smith if (OptionsHasName(0,"-mat_aij_oneindex")) b->indexshift = -1; 124117ab2063SBarry Smith 1242416022c9SBarry Smith b->m = m; 1243416022c9SBarry Smith b->n = n; 12440452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 124517ab2063SBarry Smith if (!nnz) { 124617ab2063SBarry Smith if (nz <= 0) nz = 1; 1247416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 124817ab2063SBarry Smith nz = nz*m; 124917ab2063SBarry Smith } 125017ab2063SBarry Smith else { 125117ab2063SBarry Smith nz = 0; 1252416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 125317ab2063SBarry Smith } 125417ab2063SBarry Smith 125517ab2063SBarry Smith /* allocate the matrix space */ 1256416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 12570452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1258416022c9SBarry Smith b->j = (int *) (b->a + nz); 1259cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1260416022c9SBarry Smith b->i = b->j + nz; 1261416022c9SBarry Smith b->singlemalloc = 1; 126217ab2063SBarry Smith 1263416022c9SBarry Smith b->i[0] = -b->indexshift; 126417ab2063SBarry Smith for (i=1; i<m+1; i++) { 1265416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 126617ab2063SBarry Smith } 126717ab2063SBarry Smith 1268416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 12690452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1270416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1271416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 127217ab2063SBarry Smith 1273416022c9SBarry Smith b->nz = 0; 1274416022c9SBarry Smith b->maxnz = nz; 1275416022c9SBarry Smith b->sorted = 0; 1276416022c9SBarry Smith b->roworiented = 1; 1277416022c9SBarry Smith b->nonew = 0; 1278416022c9SBarry Smith b->diag = 0; 1279416022c9SBarry Smith b->assembled = 0; 1280416022c9SBarry Smith b->solve_work = 0; 128171bd300dSLois Curfman McInnes b->spptr = 0; 1282754ec7b1SSatish Balay b->inode.node_count = 0; 1283754ec7b1SSatish Balay b->inode.size = 0; 128417ab2063SBarry Smith 1285416022c9SBarry Smith *A = B; 128617ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_superlu")) { 1287416022c9SBarry Smith ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); 128817ab2063SBarry Smith } 128917ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_essl")) { 1290416022c9SBarry Smith ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); 129117ab2063SBarry Smith } 129217ab2063SBarry Smith if (OptionsHasName(0,"-mat_aij_dxml")) { 1293416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1294416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 129517ab2063SBarry Smith } 129617ab2063SBarry Smith 129717ab2063SBarry Smith return 0; 129817ab2063SBarry Smith } 129917ab2063SBarry Smith 130008480c60SBarry Smith int MatCopyPrivate_SeqAIJ(Mat A,Mat *B,int cpvalues) 130117ab2063SBarry Smith { 1302416022c9SBarry Smith Mat C; 1303416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 130408480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 130517ab2063SBarry Smith 13064043dd9cSLois Curfman McInnes *B = 0; 1307416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatCopyPrivate_SeqAIJ:Cannot copy unassembled matrix"); 13080452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1309416022c9SBarry Smith PLogObjectCreate(C); 13100452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 131141c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1312416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1313416022c9SBarry Smith C->view = MatView_SeqAIJ; 1314416022c9SBarry Smith C->factor = A->factor; 1315416022c9SBarry Smith c->row = 0; 1316416022c9SBarry Smith c->col = 0; 1317416022c9SBarry Smith c->indexshift = shift; 131817ab2063SBarry Smith 1319416022c9SBarry Smith c->m = a->m; 1320416022c9SBarry Smith c->n = a->n; 132117ab2063SBarry Smith 13220452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 13230452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 132417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1325416022c9SBarry Smith c->imax[i] = a->imax[i]; 1326416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 132717ab2063SBarry Smith } 132817ab2063SBarry Smith 132917ab2063SBarry Smith /* allocate the matrix space */ 1330416022c9SBarry Smith c->singlemalloc = 1; 1331416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 13320452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1333416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1334416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1335416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 133617ab2063SBarry Smith if (m > 0) { 1337416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 133808480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1339416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 134017ab2063SBarry Smith } 134108480c60SBarry Smith } 134217ab2063SBarry Smith 1343416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1344416022c9SBarry Smith c->sorted = a->sorted; 1345416022c9SBarry Smith c->roworiented = a->roworiented; 1346416022c9SBarry Smith c->nonew = a->nonew; 134717ab2063SBarry Smith 1348416022c9SBarry Smith if (a->diag) { 13490452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1350416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 135117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1352416022c9SBarry Smith c->diag[i] = a->diag[i]; 135317ab2063SBarry Smith } 135417ab2063SBarry Smith } 1355416022c9SBarry Smith else c->diag = 0; 1356754ec7b1SSatish Balay if( a->inode.size){ 1357754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1358754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1359754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1360754ec7b1SSatish Balay } else { 1361754ec7b1SSatish Balay c->inode.size = 0; 1362754ec7b1SSatish Balay c->inode.node_count = 0; 1363754ec7b1SSatish Balay } 1364416022c9SBarry Smith c->assembled = 1; 1365416022c9SBarry Smith c->nz = a->nz; 1366416022c9SBarry Smith c->maxnz = a->maxnz; 1367416022c9SBarry Smith c->solve_work = 0; 136876dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1369754ec7b1SSatish Balay 1370416022c9SBarry Smith *B = C; 137117ab2063SBarry Smith return 0; 137217ab2063SBarry Smith } 137317ab2063SBarry Smith 1374416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 137517ab2063SBarry Smith { 1376416022c9SBarry Smith Mat_SeqAIJ *a; 1377416022c9SBarry Smith Mat B; 137817699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 137917ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 138017ab2063SBarry Smith MPI_Comm comm = vobj->comm; 138117ab2063SBarry Smith 138217699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 138317699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 138417ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1385416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 138648d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 138717ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 138817ab2063SBarry Smith 138917ab2063SBarry Smith /* read in row lengths */ 13900452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1391416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 139217ab2063SBarry Smith 139317ab2063SBarry Smith /* create our matrix */ 1394416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1395416022c9SBarry Smith B = *A; 1396416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1397416022c9SBarry Smith shift = a->indexshift; 139817ab2063SBarry Smith 139917ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1400416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 140117ab2063SBarry Smith if (shift) { 140217ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1403416022c9SBarry Smith a->j[i] += 1; 140417ab2063SBarry Smith } 140517ab2063SBarry Smith } 140617ab2063SBarry Smith 140717ab2063SBarry Smith /* read in nonzero values */ 1408416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 140917ab2063SBarry Smith 141017ab2063SBarry Smith /* set matrix "i" values */ 1411416022c9SBarry Smith a->i[0] = -shift; 141217ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1413416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1414416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 141517ab2063SBarry Smith } 14160452661fSBarry Smith PetscFree(rowlengths); 141717ab2063SBarry Smith 1418416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1419416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 142017ab2063SBarry Smith return 0; 142117ab2063SBarry Smith } 142217ab2063SBarry Smith 142317ab2063SBarry Smith 142417ab2063SBarry Smith 1425