117ab2063SBarry Smith #ifndef lint 2*04a348a9SBarry Smith static char vcid[] = "$Id: aij.c,v 1.137 1996/01/19 14:56:18 balay Exp bsmith $"; 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" 1317ab2063SBarry Smith 1417ab2063SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(Mat_SeqAIJ*,int**,int**); 1517ab2063SBarry Smith 16416022c9SBarry Smith static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 1717ab2063SBarry Smith { 18416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 19a2744918SBarry Smith int ierr, *ia, *ja,n,*idx,i; 203d54f168SSatish Balay /*Viewer V1, V2;*/ 2117ab2063SBarry Smith 22416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetReordering_SeqAIJ:Not for unassembled matrix"); 2317ab2063SBarry Smith 24a2744918SBarry Smith /* 25a2744918SBarry Smith this is tacky: In the future when we have written special factorization 26a2744918SBarry Smith and solve routines for the identity permutation we should use a 27a2744918SBarry Smith stride index set instead of the general one. 28a2744918SBarry Smith */ 29a2744918SBarry Smith if (type == ORDER_NATURAL) { 30a2744918SBarry Smith n = a->n; 31a2744918SBarry Smith idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 32a2744918SBarry Smith for ( i=0; i<n; i++ ) idx[i] = i; 33a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 34a2744918SBarry Smith ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 35a2744918SBarry Smith PetscFree(idx); 36a2744918SBarry Smith ISSetPermutation(*rperm); 37a2744918SBarry Smith ISSetPermutation(*cperm); 38a2744918SBarry Smith ISSetIdentity(*rperm); 39a2744918SBarry Smith ISSetIdentity(*cperm); 40a2744918SBarry Smith return 0; 41a2744918SBarry Smith } 42a2744918SBarry Smith 43416022c9SBarry Smith ierr = MatToSymmetricIJ_SeqAIJ( a, &ia, &ja ); CHKERRQ(ierr); 44416022c9SBarry Smith ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 45f31bba2fSSatish Balay 460452661fSBarry Smith PetscFree(ia); PetscFree(ja); 4717ab2063SBarry Smith return 0; 4817ab2063SBarry Smith } 4917ab2063SBarry Smith 50416022c9SBarry Smith #define CHUNKSIZE 10 5117ab2063SBarry Smith 5217ab2063SBarry Smith /* This version has row oriented v */ 53416022c9SBarry Smith static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 5417ab2063SBarry Smith { 55416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 56416022c9SBarry Smith int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 574b0e389bSBarry Smith int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 58d5d45c9bSBarry Smith int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 59416022c9SBarry Smith Scalar *ap,value, *aa = a->a; 6017ab2063SBarry Smith 6117ab2063SBarry Smith for ( k=0; k<m; k++ ) { /* loop over added rows */ 62416022c9SBarry Smith row = im[k]; 6317ab2063SBarry Smith if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 64416022c9SBarry Smith if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 6517ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 6617ab2063SBarry Smith rmax = imax[row]; nrow = ailen[row]; 67416022c9SBarry Smith low = 0; 6817ab2063SBarry Smith for ( l=0; l<n; l++ ) { /* loop over added columns */ 69416022c9SBarry Smith if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 70416022c9SBarry Smith if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 714b0e389bSBarry Smith col = in[l] - shift; 724b0e389bSBarry Smith if (roworiented) { 734b0e389bSBarry Smith value = *v++; 744b0e389bSBarry Smith } 754b0e389bSBarry Smith else { 764b0e389bSBarry Smith value = v[k + l*m]; 774b0e389bSBarry Smith } 78416022c9SBarry Smith if (!sorted) low = 0; high = nrow; 79416022c9SBarry Smith while (high-low > 5) { 80416022c9SBarry Smith t = (low+high)/2; 81416022c9SBarry Smith if (rp[t] > col) high = t; 82416022c9SBarry Smith else low = t; 8317ab2063SBarry Smith } 84416022c9SBarry Smith for ( i=low; i<high; i++ ) { 8517ab2063SBarry Smith if (rp[i] > col) break; 8617ab2063SBarry Smith if (rp[i] == col) { 87416022c9SBarry Smith if (is == ADD_VALUES) ap[i] += value; 8817ab2063SBarry Smith else ap[i] = value; 8917ab2063SBarry Smith goto noinsert; 9017ab2063SBarry Smith } 9117ab2063SBarry Smith } 9217ab2063SBarry Smith if (nonew) goto noinsert; 9317ab2063SBarry Smith if (nrow >= rmax) { 9417ab2063SBarry Smith /* there is no extra room in row, therefore enlarge */ 95416022c9SBarry Smith int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 9617ab2063SBarry Smith Scalar *new_a; 9717ab2063SBarry Smith 9817ab2063SBarry Smith /* malloc new storage space */ 99416022c9SBarry Smith len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 1000452661fSBarry Smith new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 10117ab2063SBarry Smith new_j = (int *) (new_a + new_nz); 10217ab2063SBarry Smith new_i = new_j + new_nz; 10317ab2063SBarry Smith 10417ab2063SBarry Smith /* copy over old data into new slots */ 10517ab2063SBarry Smith for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 106416022c9SBarry Smith for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 107416022c9SBarry Smith PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 108416022c9SBarry Smith len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 109416022c9SBarry Smith PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 11017ab2063SBarry Smith len*sizeof(int)); 111416022c9SBarry Smith PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 112416022c9SBarry Smith PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 11317ab2063SBarry Smith len*sizeof(Scalar)); 11417ab2063SBarry Smith /* free up old matrix storage */ 1150452661fSBarry Smith PetscFree(a->a); 1160452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 117416022c9SBarry Smith aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 118416022c9SBarry Smith a->singlemalloc = 1; 11917ab2063SBarry Smith 12017ab2063SBarry Smith rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 121416022c9SBarry Smith rmax = imax[row] = imax[row] + CHUNKSIZE; 122416022c9SBarry Smith PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 123416022c9SBarry Smith a->maxnz += CHUNKSIZE; 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 if (!a->assembled) SETERRQ(1,"MatGetValues_SeqAIJ:Not for unassembled matrix"); 1497eb43aa7SLois Curfman McInnes for ( k=0; k<m; k++ ) { /* loop over rows */ 1507eb43aa7SLois Curfman McInnes row = im[k]; 1517eb43aa7SLois Curfman McInnes if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 1527eb43aa7SLois Curfman McInnes if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 1537eb43aa7SLois Curfman McInnes rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 1547eb43aa7SLois Curfman McInnes nrow = ailen[row]; 1557eb43aa7SLois Curfman McInnes for ( l=0; l<n; l++ ) { /* loop over columns */ 1567eb43aa7SLois Curfman McInnes if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 1577eb43aa7SLois Curfman McInnes if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 1587eb43aa7SLois Curfman McInnes col = in[l] - shift; 1597eb43aa7SLois Curfman McInnes high = nrow; low = 0; /* assume unsorted */ 1607eb43aa7SLois Curfman McInnes while (high-low > 5) { 1617eb43aa7SLois Curfman McInnes t = (low+high)/2; 1627eb43aa7SLois Curfman McInnes if (rp[t] > col) high = t; 1637eb43aa7SLois Curfman McInnes else low = t; 1647eb43aa7SLois Curfman McInnes } 1657eb43aa7SLois Curfman McInnes for ( i=low; i<high; i++ ) { 1667eb43aa7SLois Curfman McInnes if (rp[i] > col) break; 1677eb43aa7SLois Curfman McInnes if (rp[i] == col) { 168b49de8d1SLois Curfman McInnes *v++ = ap[i]; 1697eb43aa7SLois Curfman McInnes goto finished; 1707eb43aa7SLois Curfman McInnes } 1717eb43aa7SLois Curfman McInnes } 172b49de8d1SLois Curfman McInnes *v++ = zero; 1737eb43aa7SLois Curfman McInnes finished:; 1747eb43aa7SLois Curfman McInnes } 1757eb43aa7SLois Curfman McInnes } 1767eb43aa7SLois Curfman McInnes return 0; 1777eb43aa7SLois Curfman McInnes } 1787eb43aa7SLois Curfman McInnes 17917ab2063SBarry Smith #include "draw.h" 18017ab2063SBarry Smith #include "pinclude/pviewer.h" 181416022c9SBarry Smith #include "sysio.h" 18217ab2063SBarry Smith 183416022c9SBarry Smith static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 18417ab2063SBarry Smith { 185416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 186416022c9SBarry Smith int i, fd, *col_lens, ierr; 18717ab2063SBarry Smith 188416022c9SBarry Smith ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 1890452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 190416022c9SBarry Smith col_lens[0] = MAT_COOKIE; 191416022c9SBarry Smith col_lens[1] = a->m; 192416022c9SBarry Smith col_lens[2] = a->n; 193416022c9SBarry Smith col_lens[3] = a->nz; 194416022c9SBarry Smith 195416022c9SBarry Smith /* store lengths of each row and write (including header) to file */ 196416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 197416022c9SBarry Smith col_lens[4+i] = a->i[i+1] - a->i[i]; 19817ab2063SBarry Smith } 199416022c9SBarry Smith ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 2000452661fSBarry Smith PetscFree(col_lens); 201416022c9SBarry Smith 202416022c9SBarry Smith /* store column indices (zero start index) */ 203416022c9SBarry Smith if (a->indexshift) { 204416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]--; 20517ab2063SBarry Smith } 206416022c9SBarry Smith ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 207416022c9SBarry Smith if (a->indexshift) { 208416022c9SBarry Smith for ( i=0; i<a->nz; i++ ) a->j[i]++; 20917ab2063SBarry Smith } 210416022c9SBarry Smith 211416022c9SBarry Smith /* store nonzero values */ 212416022c9SBarry Smith ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 21317ab2063SBarry Smith return 0; 21417ab2063SBarry Smith } 215416022c9SBarry Smith 216416022c9SBarry Smith static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 217416022c9SBarry Smith { 218416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 219416022c9SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,format; 22017ab2063SBarry Smith FILE *fd; 22117ab2063SBarry Smith char *outputname; 22217ab2063SBarry Smith 223416022c9SBarry Smith ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 224416022c9SBarry Smith ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 225416022c9SBarry Smith ierr = ViewerFileGetFormat_Private(viewer,&format); 22617ab2063SBarry Smith if (format == FILE_FORMAT_INFO) { 22708480c60SBarry Smith /* no need to print additional information */ ; 22817ab2063SBarry Smith } 22917ab2063SBarry Smith else if (format == FILE_FORMAT_MATLAB) { 23017ab2063SBarry Smith int nz, nzalloc, mem; 231416022c9SBarry Smith MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 232416022c9SBarry Smith fprintf(fd,"%% Size = %d %d \n",m,a->n); 23317ab2063SBarry Smith fprintf(fd,"%% Nonzeros = %d \n",nz); 23417ab2063SBarry Smith fprintf(fd,"zzz = zeros(%d,3);\n",nz); 23517ab2063SBarry Smith fprintf(fd,"zzz = [\n"); 23617ab2063SBarry Smith 23717ab2063SBarry Smith for (i=0; i<m; i++) { 238416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 23917ab2063SBarry Smith #if defined(PETSC_COMPLEX) 240416022c9SBarry Smith fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j],real(a->a[j]), 241416022c9SBarry Smith imag(a->a[j])); 24217ab2063SBarry Smith #else 243416022c9SBarry Smith fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j], a->a[j]); 24417ab2063SBarry Smith #endif 24517ab2063SBarry Smith } 24617ab2063SBarry Smith } 24717ab2063SBarry Smith fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 24817ab2063SBarry Smith } 24917ab2063SBarry Smith else { 25017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 25117ab2063SBarry Smith fprintf(fd,"row %d:",i); 252416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 25317ab2063SBarry Smith #if defined(PETSC_COMPLEX) 254416022c9SBarry Smith if (imag(a->a[j]) != 0.0) { 255416022c9SBarry Smith fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 25617ab2063SBarry Smith } 25717ab2063SBarry Smith else { 258416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 25917ab2063SBarry Smith } 26017ab2063SBarry Smith #else 261416022c9SBarry Smith fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 26217ab2063SBarry Smith #endif 26317ab2063SBarry Smith } 26417ab2063SBarry Smith fprintf(fd,"\n"); 26517ab2063SBarry Smith } 26617ab2063SBarry Smith } 26717ab2063SBarry Smith fflush(fd); 268416022c9SBarry Smith return 0; 269416022c9SBarry Smith } 270416022c9SBarry Smith 271416022c9SBarry Smith static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 272416022c9SBarry Smith { 273416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 274cddf8d76SBarry Smith int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 275cddf8d76SBarry Smith double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 276d7e8b826SBarry Smith Draw draw = (Draw) viewer; 277cddf8d76SBarry Smith DrawButton button; 278cddf8d76SBarry Smith 279416022c9SBarry Smith xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 280416022c9SBarry Smith xr += w; yr += h; xl = -w; yl = -h; 281416022c9SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 282416022c9SBarry Smith /* loop over matrix elements drawing boxes */ 283cddf8d76SBarry Smith color = DRAW_BLUE; 284416022c9SBarry Smith for ( i=0; i<m; i++ ) { 285cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 286416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 287cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 288cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 289cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 290cddf8d76SBarry Smith #else 291cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 292cddf8d76SBarry Smith #endif 293cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 294cddf8d76SBarry Smith } 295cddf8d76SBarry Smith } 296cddf8d76SBarry Smith color = DRAW_CYAN; 297cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 298cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 299cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 300cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 301cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 302cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 303cddf8d76SBarry Smith } 304cddf8d76SBarry Smith } 305cddf8d76SBarry Smith color = DRAW_RED; 306cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 307cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 308cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 309cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 310cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 311cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 312cddf8d76SBarry Smith #else 313cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 314cddf8d76SBarry Smith #endif 315cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 316416022c9SBarry Smith } 317416022c9SBarry Smith } 318416022c9SBarry Smith DrawFlush(draw); 319cddf8d76SBarry Smith DrawGetPause(draw,&pause); 320cddf8d76SBarry Smith if (pause >= 0) { PetscSleep(pause); return 0;} 321cddf8d76SBarry Smith 322cddf8d76SBarry Smith /* allow the matrix to zoom or shrink */ 323cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 324cddf8d76SBarry Smith while (button != BUTTON_RIGHT) { 325cddf8d76SBarry Smith DrawClear(draw); 326cddf8d76SBarry Smith if (button == BUTTON_LEFT) scale = .5; 327cddf8d76SBarry Smith else if (button == BUTTON_CENTER) scale = 2.; 328cddf8d76SBarry Smith xl = scale*(xl + w - xc) + xc - w*scale; 329cddf8d76SBarry Smith xr = scale*(xr - w - xc) + xc + w*scale; 330cddf8d76SBarry Smith yl = scale*(yl + h - yc) + yc - h*scale; 331cddf8d76SBarry Smith yr = scale*(yr - h - yc) + yc + h*scale; 332cddf8d76SBarry Smith w *= scale; h *= scale; 333cddf8d76SBarry Smith ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 334cddf8d76SBarry Smith color = DRAW_BLUE; 335cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 336cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 337cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 338cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 339cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 340cddf8d76SBarry Smith if (real(a->a[j]) >= 0.) continue; 341cddf8d76SBarry Smith #else 342cddf8d76SBarry Smith if (a->a[j] >= 0.) continue; 343cddf8d76SBarry Smith #endif 344cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 345cddf8d76SBarry Smith } 346cddf8d76SBarry Smith } 347cddf8d76SBarry Smith color = DRAW_CYAN; 348cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 349cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 350cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 351cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 352cddf8d76SBarry Smith if (a->a[j] != 0.) continue; 353cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 354cddf8d76SBarry Smith } 355cddf8d76SBarry Smith } 356cddf8d76SBarry Smith color = DRAW_RED; 357cddf8d76SBarry Smith for ( i=0; i<m; i++ ) { 358cddf8d76SBarry Smith y_l = m - i - 1.0; y_r = y_l + 1.0; 359cddf8d76SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 360cddf8d76SBarry Smith x_l = a->j[j] + shift; x_r = x_l + 1.0; 361cddf8d76SBarry Smith #if defined(PETSC_COMPLEX) 362cddf8d76SBarry Smith if (real(a->a[j]) <= 0.) continue; 363cddf8d76SBarry Smith #else 364cddf8d76SBarry Smith if (a->a[j] <= 0.) continue; 365cddf8d76SBarry Smith #endif 366cddf8d76SBarry Smith DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 367cddf8d76SBarry Smith } 368cddf8d76SBarry Smith } 369cddf8d76SBarry Smith ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 370cddf8d76SBarry Smith } 371416022c9SBarry Smith return 0; 372416022c9SBarry Smith } 373416022c9SBarry Smith 374416022c9SBarry Smith static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 375416022c9SBarry Smith { 376416022c9SBarry Smith Mat A = (Mat) obj; 377416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 378416022c9SBarry Smith PetscObject vobj = (PetscObject) viewer; 379416022c9SBarry Smith 380416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatView_SeqAIJ:Not for unassembled matrix"); 381416022c9SBarry Smith if (!viewer) { 382416022c9SBarry Smith viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 383416022c9SBarry Smith } 384416022c9SBarry Smith if (vobj->cookie == VIEWER_COOKIE) { 385416022c9SBarry Smith if (vobj->type == MATLAB_VIEWER) { 386416022c9SBarry Smith return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 387416022c9SBarry Smith } 388416022c9SBarry Smith else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 389416022c9SBarry Smith return MatView_SeqAIJ_ASCII(A,viewer); 390416022c9SBarry Smith } 391416022c9SBarry Smith else if (vobj->type == BINARY_FILE_VIEWER) { 392416022c9SBarry Smith return MatView_SeqAIJ_Binary(A,viewer); 393416022c9SBarry Smith } 394416022c9SBarry Smith } 395416022c9SBarry Smith else if (vobj->cookie == DRAW_COOKIE) { 396416022c9SBarry Smith if (vobj->type == NULLWINDOW) return 0; 397416022c9SBarry Smith else return MatView_SeqAIJ_Draw(A,viewer); 39817ab2063SBarry Smith } 39917ab2063SBarry Smith return 0; 40017ab2063SBarry Smith } 40141c01911SSatish Balay int Mat_AIJ_CheckInode(Mat); 402416022c9SBarry Smith static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 40317ab2063SBarry Smith { 404416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 40541c01911SSatish Balay int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 406416022c9SBarry Smith int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 407416022c9SBarry Smith Scalar *aa = a->a, *ap; 40817ab2063SBarry Smith 40917ab2063SBarry Smith if (mode == FLUSH_ASSEMBLY) return 0; 41017ab2063SBarry Smith 41117ab2063SBarry Smith for ( i=1; i<m; i++ ) { 412416022c9SBarry Smith /* move each row back by the amount of empty slots (fshift) before it*/ 41317ab2063SBarry Smith fshift += imax[i-1] - ailen[i-1]; 41417ab2063SBarry Smith if (fshift) { 415416022c9SBarry Smith ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 41617ab2063SBarry Smith N = ailen[i]; 41717ab2063SBarry Smith for ( j=0; j<N; j++ ) { 41817ab2063SBarry Smith ip[j-fshift] = ip[j]; 41917ab2063SBarry Smith ap[j-fshift] = ap[j]; 42017ab2063SBarry Smith } 42117ab2063SBarry Smith } 42217ab2063SBarry Smith ai[i] = ai[i-1] + ailen[i-1]; 42317ab2063SBarry Smith } 42417ab2063SBarry Smith if (m) { 42517ab2063SBarry Smith fshift += imax[m-1] - ailen[m-1]; 42617ab2063SBarry Smith ai[m] = ai[m-1] + ailen[m-1]; 42717ab2063SBarry Smith } 42817ab2063SBarry Smith /* reset ilen and imax for each row */ 42917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 43017ab2063SBarry Smith ailen[i] = imax[i] = ai[i+1] - ai[i]; 43117ab2063SBarry Smith } 432416022c9SBarry Smith a->nz = ai[m] + shift; 43317ab2063SBarry Smith 43417ab2063SBarry Smith /* diagonals may have moved, so kill the diagonal pointers */ 435416022c9SBarry Smith if (fshift && a->diag) { 4360452661fSBarry Smith PetscFree(a->diag); 437416022c9SBarry Smith PLogObjectMemory(A,-(m+1)*sizeof(int)); 438416022c9SBarry Smith a->diag = 0; 43917ab2063SBarry Smith } 44076dd722bSSatish Balay /* check out for identical nodes. If found, use inode functions */ 44141c01911SSatish Balay ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 442416022c9SBarry Smith a->assembled = 1; 44317ab2063SBarry Smith return 0; 44417ab2063SBarry Smith } 44517ab2063SBarry Smith 446416022c9SBarry Smith static int MatZeroEntries_SeqAIJ(Mat A) 44717ab2063SBarry Smith { 448416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 449cddf8d76SBarry Smith PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 45017ab2063SBarry Smith return 0; 45117ab2063SBarry Smith } 452416022c9SBarry Smith 45317ab2063SBarry Smith int MatDestroy_SeqAIJ(PetscObject obj) 45417ab2063SBarry Smith { 455416022c9SBarry Smith Mat A = (Mat) obj; 456416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 457d5d45c9bSBarry Smith 45817ab2063SBarry Smith #if defined(PETSC_LOG) 459416022c9SBarry Smith PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 46017ab2063SBarry Smith #endif 4610452661fSBarry Smith PetscFree(a->a); 4620452661fSBarry Smith if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 4630452661fSBarry Smith if (a->diag) PetscFree(a->diag); 4640452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 4650452661fSBarry Smith if (a->imax) PetscFree(a->imax); 4660452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 46776dd722bSSatish Balay if (a->inode.size) PetscFree(a->inode.size); 4680452661fSBarry Smith PetscFree(a); 469416022c9SBarry Smith PLogObjectDestroy(A); 4700452661fSBarry Smith PetscHeaderDestroy(A); 47117ab2063SBarry Smith return 0; 47217ab2063SBarry Smith } 47317ab2063SBarry Smith 474416022c9SBarry Smith static int MatCompress_SeqAIJ(Mat A) 47517ab2063SBarry Smith { 47617ab2063SBarry Smith return 0; 47717ab2063SBarry Smith } 47817ab2063SBarry Smith 479416022c9SBarry Smith static int MatSetOption_SeqAIJ(Mat A,MatOption op) 48017ab2063SBarry Smith { 481416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 482416022c9SBarry Smith if (op == ROW_ORIENTED) a->roworiented = 1; 4834b0e389bSBarry Smith else if (op == COLUMN_ORIENTED) a->roworiented = 0; 484416022c9SBarry Smith else if (op == COLUMNS_SORTED) a->sorted = 1; 485416022c9SBarry Smith else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 486416022c9SBarry Smith else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 487e2f28af5SBarry Smith else if (op == ROWS_SORTED || 488e2f28af5SBarry Smith op == SYMMETRIC_MATRIX || 489e2f28af5SBarry Smith op == STRUCTURALLY_SYMMETRIC_MATRIX || 490e2f28af5SBarry Smith op == YES_NEW_DIAGONALS) 491df719601SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 492df719601SLois Curfman McInnes else if (op == NO_NEW_DIAGONALS) 4934d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 4946c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_1) a->inode.limit = 1; 4956c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_2) a->inode.limit = 2; 4966c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_3) a->inode.limit = 3; 4976c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_4) a->inode.limit = 4; 4986c7ebb05SLois Curfman McInnes else if (op == INODE_LIMIT_5) a->inode.limit = 5; 499e2f28af5SBarry Smith else 5004d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 50117ab2063SBarry Smith return 0; 50217ab2063SBarry Smith } 50317ab2063SBarry Smith 504416022c9SBarry Smith static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 50517ab2063SBarry Smith { 506416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 507416022c9SBarry Smith int i,j, n,shift = a->indexshift; 50817ab2063SBarry Smith Scalar *x, zero = 0.0; 50917ab2063SBarry Smith 510416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Not for unassembled matrix"); 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 534416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTrans_SeqAIJ:Not for unassembled matrix"); 53517ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 536cddf8d76SBarry Smith PetscMemzero(y,a->n*sizeof(Scalar)); 53717ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 53817ab2063SBarry Smith for ( i=0; i<m; i++ ) { 539416022c9SBarry Smith idx = a->j + a->i[i] + shift; 540416022c9SBarry Smith v = a->a + a->i[i] + shift; 541416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 54217ab2063SBarry Smith alpha = x[i]; 54317ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 54417ab2063SBarry Smith } 545416022c9SBarry Smith PLogFlops(2*a->nz - a->n); 54617ab2063SBarry Smith return 0; 54717ab2063SBarry Smith } 548d5d45c9bSBarry Smith 549416022c9SBarry Smith static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 55017ab2063SBarry Smith { 551416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 55217ab2063SBarry Smith Scalar *x, *y, *v, alpha; 553416022c9SBarry Smith int m = a->m, n, i, *idx,shift = a->indexshift; 55417ab2063SBarry Smith 555416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultTransAdd_SeqAIJ:Not for unassembled matrix"); 55617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 55717ab2063SBarry Smith if (zz != yy) VecCopy(zz,yy); 55817ab2063SBarry Smith y = y + shift; /* shift for Fortran start by 1 indexing */ 55917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 560416022c9SBarry Smith idx = a->j + a->i[i] + shift; 561416022c9SBarry Smith v = a->a + a->i[i] + shift; 562416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 56317ab2063SBarry Smith alpha = x[i]; 56417ab2063SBarry Smith while (n-->0) {y[*idx++] += alpha * *v++;} 56517ab2063SBarry Smith } 56617ab2063SBarry Smith return 0; 56717ab2063SBarry Smith } 56817ab2063SBarry Smith 569416022c9SBarry Smith static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 57017ab2063SBarry Smith { 571416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 57217ab2063SBarry Smith Scalar *x, *y, *v, sum; 573416022c9SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 57417ab2063SBarry Smith 575416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatMult_SeqAIJ:Not for unassembled matrix"); 57617ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 57717ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 578416022c9SBarry Smith idx = a->j; 579416022c9SBarry Smith v = a->a; 580416022c9SBarry Smith ii = a->i; 58117ab2063SBarry Smith for ( i=0; i<m; i++ ) { 582416022c9SBarry Smith n = ii[1] - ii[0]; ii++; 58317ab2063SBarry Smith sum = 0.0; 58417ab2063SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 58517ab2063SBarry Smith /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 586416022c9SBarry Smith while (n--) sum += *v++ * x[*idx++]; 58717ab2063SBarry Smith y[i] = sum; 58817ab2063SBarry Smith } 589416022c9SBarry Smith PLogFlops(2*a->nz - m); 59017ab2063SBarry Smith return 0; 59117ab2063SBarry Smith } 59217ab2063SBarry Smith 593416022c9SBarry Smith static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 59417ab2063SBarry Smith { 595416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 59617ab2063SBarry Smith Scalar *x, *y, *z, *v, sum; 597cddf8d76SBarry Smith int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 59817ab2063SBarry Smith 59948d91487SBarry Smith if (!a->assembled) SETERRQ(1,"MatMultAdd_SeqAIJ:Not for unassembled matrix"); 60017ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 60117ab2063SBarry Smith x = x + shift; /* shift for Fortran start by 1 indexing */ 602cddf8d76SBarry Smith idx = a->j; 603cddf8d76SBarry Smith v = a->a; 604cddf8d76SBarry Smith ii = a->i; 60517ab2063SBarry Smith for ( i=0; i<m; i++ ) { 606cddf8d76SBarry Smith n = ii[1] - ii[0]; ii++; 60717ab2063SBarry Smith sum = y[i]; 608cddf8d76SBarry Smith /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 609cddf8d76SBarry Smith while (n--) sum += *v++ * x[*idx++]; 61017ab2063SBarry Smith z[i] = sum; 61117ab2063SBarry Smith } 612416022c9SBarry Smith PLogFlops(2*a->nz); 61317ab2063SBarry Smith return 0; 61417ab2063SBarry Smith } 61517ab2063SBarry Smith 61617ab2063SBarry Smith /* 61717ab2063SBarry Smith Adds diagonal pointers to sparse matrix structure. 61817ab2063SBarry Smith */ 61917ab2063SBarry Smith 620416022c9SBarry Smith int MatMarkDiag_SeqAIJ(Mat A) 62117ab2063SBarry Smith { 622416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 623416022c9SBarry Smith int i,j, *diag, m = a->m,shift = a->indexshift; 62417ab2063SBarry Smith 625f0b747eeSBarry Smith if (!a->assembled) SETERRQ(1,"MatMarkDiag_SeqAIJ:Not for unassembled matrix"); 6260452661fSBarry Smith diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 627416022c9SBarry Smith PLogObjectMemory(A,(m+1)*sizeof(int)); 628416022c9SBarry Smith for ( i=0; i<a->m; i++ ) { 629416022c9SBarry Smith for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 630416022c9SBarry Smith if (a->j[j]+shift == i) { 63117ab2063SBarry Smith diag[i] = j - shift; 63217ab2063SBarry Smith break; 63317ab2063SBarry Smith } 63417ab2063SBarry Smith } 63517ab2063SBarry Smith } 636416022c9SBarry Smith a->diag = diag; 63717ab2063SBarry Smith return 0; 63817ab2063SBarry Smith } 63917ab2063SBarry Smith 640416022c9SBarry Smith static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 64117ab2063SBarry Smith double fshift,int its,Vec xx) 64217ab2063SBarry Smith { 643416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 644416022c9SBarry Smith Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 645d5d45c9bSBarry Smith int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 64617ab2063SBarry Smith 64717ab2063SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 648416022c9SBarry Smith if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 649416022c9SBarry Smith diag = a->diag; 650416022c9SBarry Smith xs = x + shift; /* shifted by one for index start of a or a->j*/ 65117ab2063SBarry Smith if (flag == SOR_APPLY_UPPER) { 65217ab2063SBarry Smith /* apply ( U + D/omega) to the vector */ 65317ab2063SBarry Smith bs = b + shift; 65417ab2063SBarry Smith for ( i=0; i<m; i++ ) { 655416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 656416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 657416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 658416022c9SBarry Smith v = a->a + diag[i] + (!shift); 65917ab2063SBarry Smith sum = b[i]*d/omega; 66017ab2063SBarry Smith SPARSEDENSEDOT(sum,bs,v,idx,n); 66117ab2063SBarry Smith x[i] = sum; 66217ab2063SBarry Smith } 66317ab2063SBarry Smith return 0; 66417ab2063SBarry Smith } 66517ab2063SBarry Smith if (flag == SOR_APPLY_LOWER) { 66617ab2063SBarry Smith SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 66717ab2063SBarry Smith } 668416022c9SBarry Smith else if (flag & SOR_EISENSTAT) { 66917ab2063SBarry Smith /* Let A = L + U + D; where L is lower trianglar, 67017ab2063SBarry Smith U is upper triangular, E is diagonal; This routine applies 67117ab2063SBarry Smith 67217ab2063SBarry Smith (L + E)^{-1} A (U + E)^{-1} 67317ab2063SBarry Smith 67417ab2063SBarry Smith to a vector efficiently using Eisenstat's trick. This is for 67517ab2063SBarry Smith the case of SSOR preconditioner, so E is D/omega where omega 67617ab2063SBarry Smith is the relaxation factor. 67717ab2063SBarry Smith */ 6780452661fSBarry Smith t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 67917ab2063SBarry Smith scale = (2.0/omega) - 1.0; 68017ab2063SBarry Smith 68117ab2063SBarry Smith /* x = (E + U)^{-1} b */ 68217ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 683416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 684416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 685416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 686416022c9SBarry Smith v = a->a + diag[i] + (!shift); 68717ab2063SBarry Smith sum = b[i]; 68817ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 68917ab2063SBarry Smith x[i] = omega*(sum/d); 69017ab2063SBarry Smith } 69117ab2063SBarry Smith 69217ab2063SBarry Smith /* t = b - (2*E - D)x */ 693416022c9SBarry Smith v = a->a; 69417ab2063SBarry Smith for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 69517ab2063SBarry Smith 69617ab2063SBarry Smith /* t = (E + L)^{-1}t */ 697416022c9SBarry Smith ts = t + shift; /* shifted by one for index start of a or a->j*/ 698416022c9SBarry Smith diag = a->diag; 69917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 700416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 701416022c9SBarry Smith n = diag[i] - a->i[i]; 702416022c9SBarry Smith idx = a->j + a->i[i] + shift; 703416022c9SBarry Smith v = a->a + a->i[i] + shift; 70417ab2063SBarry Smith sum = t[i]; 70517ab2063SBarry Smith SPARSEDENSEMDOT(sum,ts,v,idx,n); 70617ab2063SBarry Smith t[i] = omega*(sum/d); 70717ab2063SBarry Smith } 70817ab2063SBarry Smith 70917ab2063SBarry Smith /* x = x + t */ 71017ab2063SBarry Smith for ( i=0; i<m; i++ ) { x[i] += t[i]; } 7110452661fSBarry Smith PetscFree(t); 71217ab2063SBarry Smith return 0; 71317ab2063SBarry Smith } 71417ab2063SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 71517ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 71617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 717416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 718416022c9SBarry Smith n = diag[i] - a->i[i]; 719416022c9SBarry Smith idx = a->j + a->i[i] + shift; 720416022c9SBarry Smith v = a->a + a->i[i] + shift; 72117ab2063SBarry Smith sum = b[i]; 72217ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 72317ab2063SBarry Smith x[i] = omega*(sum/d); 72417ab2063SBarry Smith } 72517ab2063SBarry Smith xb = x; 72617ab2063SBarry Smith } 72717ab2063SBarry Smith else xb = b; 72817ab2063SBarry Smith if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 72917ab2063SBarry Smith (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 73017ab2063SBarry Smith for ( i=0; i<m; i++ ) { 731416022c9SBarry Smith x[i] *= a->a[diag[i]+shift]; 73217ab2063SBarry Smith } 73317ab2063SBarry Smith } 73417ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 73517ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 736416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 737416022c9SBarry Smith n = a->i[i+1] - diag[i] - 1; 738416022c9SBarry Smith idx = a->j + diag[i] + (!shift); 739416022c9SBarry Smith v = a->a + diag[i] + (!shift); 74017ab2063SBarry Smith sum = xb[i]; 74117ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 74217ab2063SBarry Smith x[i] = omega*(sum/d); 74317ab2063SBarry Smith } 74417ab2063SBarry Smith } 74517ab2063SBarry Smith its--; 74617ab2063SBarry Smith } 74717ab2063SBarry Smith while (its--) { 74817ab2063SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 74917ab2063SBarry Smith for ( i=0; i<m; i++ ) { 750416022c9SBarry Smith d = fshift + a->a[diag[i]+shift]; 751416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 752416022c9SBarry Smith idx = a->j + a->i[i] + shift; 753416022c9SBarry Smith v = a->a + a->i[i] + shift; 75417ab2063SBarry Smith sum = b[i]; 75517ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 75617ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 75717ab2063SBarry Smith } 75817ab2063SBarry Smith } 75917ab2063SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 76017ab2063SBarry Smith for ( i=m-1; i>=0; i-- ) { 761416022c9SBarry Smith d = fshift + a->a[diag[i] + shift]; 762416022c9SBarry Smith n = a->i[i+1] - a->i[i]; 763416022c9SBarry Smith idx = a->j + a->i[i] + shift; 764416022c9SBarry Smith v = a->a + a->i[i] + shift; 76517ab2063SBarry Smith sum = b[i]; 76617ab2063SBarry Smith SPARSEDENSEMDOT(sum,xs,v,idx,n); 76717ab2063SBarry Smith x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 76817ab2063SBarry Smith } 76917ab2063SBarry Smith } 77017ab2063SBarry Smith } 77117ab2063SBarry Smith return 0; 77217ab2063SBarry Smith } 77317ab2063SBarry Smith 774d5d45c9bSBarry Smith static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 77517ab2063SBarry Smith { 776416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 777416022c9SBarry Smith *nz = a->nz; 778416022c9SBarry Smith *nzalloc = a->maxnz; 779416022c9SBarry Smith *mem = (int)A->mem; 78017ab2063SBarry Smith return 0; 78117ab2063SBarry Smith } 78217ab2063SBarry Smith 78317ab2063SBarry Smith extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 78417ab2063SBarry Smith extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 78517ab2063SBarry Smith extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 78617ab2063SBarry Smith extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 78717ab2063SBarry Smith extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 78817ab2063SBarry Smith extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 78917ab2063SBarry Smith extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 79017ab2063SBarry Smith 79117ab2063SBarry Smith static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 79217ab2063SBarry Smith { 793416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 794416022c9SBarry Smith int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 79517ab2063SBarry Smith 79617ab2063SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 79717ab2063SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 79817ab2063SBarry Smith if (diag) { 79917ab2063SBarry Smith for ( i=0; i<N; i++ ) { 800416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 801416022c9SBarry Smith if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 802416022c9SBarry Smith a->ilen[rows[i]] = 1; 803416022c9SBarry Smith a->a[a->i[rows[i]]+shift] = *diag; 804416022c9SBarry Smith a->j[a->i[rows[i]]+shift] = rows[i]+shift; 80517ab2063SBarry Smith } 80617ab2063SBarry Smith else { 80717ab2063SBarry Smith ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 80817ab2063SBarry Smith CHKERRQ(ierr); 80917ab2063SBarry Smith } 81017ab2063SBarry Smith } 81117ab2063SBarry Smith } 81217ab2063SBarry Smith else { 81317ab2063SBarry Smith for ( i=0; i<N; i++ ) { 814416022c9SBarry Smith if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 815416022c9SBarry Smith a->ilen[rows[i]] = 0; 81617ab2063SBarry Smith } 81717ab2063SBarry Smith } 81817ab2063SBarry Smith ISRestoreIndices(is,&rows); 81917ab2063SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 82017ab2063SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 82117ab2063SBarry Smith return 0; 82217ab2063SBarry Smith } 82317ab2063SBarry Smith 824416022c9SBarry Smith static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 82517ab2063SBarry Smith { 826416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 827416022c9SBarry Smith *m = a->m; *n = a->n; 82817ab2063SBarry Smith return 0; 82917ab2063SBarry Smith } 83017ab2063SBarry Smith 831416022c9SBarry Smith static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 83217ab2063SBarry Smith { 833416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 834416022c9SBarry Smith *m = 0; *n = a->m; 83517ab2063SBarry Smith return 0; 83617ab2063SBarry Smith } 837416022c9SBarry Smith static int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 83817ab2063SBarry Smith { 839416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 840416022c9SBarry Smith int *itmp,i,ierr,shift = a->indexshift; 84117ab2063SBarry Smith 842416022c9SBarry Smith if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 84317ab2063SBarry Smith 844416022c9SBarry Smith if (!a->assembled) { 845416022c9SBarry Smith ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 846416022c9SBarry Smith ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 84717ab2063SBarry Smith } 848416022c9SBarry Smith *nz = a->i[row+1] - a->i[row]; 849416022c9SBarry Smith if (v) *v = a->a + a->i[row] + shift; 85017ab2063SBarry Smith if (idx) { 85117ab2063SBarry Smith if (*nz) { 852416022c9SBarry Smith itmp = a->j + a->i[row] + shift; 8530452661fSBarry Smith *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 85417ab2063SBarry Smith for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 85517ab2063SBarry Smith } 85617ab2063SBarry Smith else *idx = 0; 85717ab2063SBarry Smith } 85817ab2063SBarry Smith return 0; 85917ab2063SBarry Smith } 86017ab2063SBarry Smith 861416022c9SBarry Smith static int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 86217ab2063SBarry Smith { 8630452661fSBarry Smith if (idx) {if (*idx) PetscFree(*idx);} 86417ab2063SBarry Smith return 0; 86517ab2063SBarry Smith } 86617ab2063SBarry Smith 867cddf8d76SBarry Smith static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 86817ab2063SBarry Smith { 869416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 870416022c9SBarry Smith Scalar *v = a->a; 87117ab2063SBarry Smith double sum = 0.0; 872416022c9SBarry Smith int i, j,shift = a->indexshift; 87317ab2063SBarry Smith 874416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatNorm_SeqAIJ:Not for unassembled matrix"); 87517ab2063SBarry Smith if (type == NORM_FROBENIUS) { 876416022c9SBarry Smith for (i=0; i<a->nz; i++ ) { 87717ab2063SBarry Smith #if defined(PETSC_COMPLEX) 87817ab2063SBarry Smith sum += real(conj(*v)*(*v)); v++; 87917ab2063SBarry Smith #else 88017ab2063SBarry Smith sum += (*v)*(*v); v++; 88117ab2063SBarry Smith #endif 88217ab2063SBarry Smith } 88317ab2063SBarry Smith *norm = sqrt(sum); 88417ab2063SBarry Smith } 88517ab2063SBarry Smith else if (type == NORM_1) { 88617ab2063SBarry Smith double *tmp; 887416022c9SBarry Smith int *jj = a->j; 8880452661fSBarry Smith tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 889cddf8d76SBarry Smith PetscMemzero(tmp,a->n*sizeof(double)); 89017ab2063SBarry Smith *norm = 0.0; 891416022c9SBarry Smith for ( j=0; j<a->nz; j++ ) { 892a2744918SBarry Smith tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 89317ab2063SBarry Smith } 894416022c9SBarry Smith for ( j=0; j<a->n; j++ ) { 89517ab2063SBarry Smith if (tmp[j] > *norm) *norm = tmp[j]; 89617ab2063SBarry Smith } 8970452661fSBarry Smith PetscFree(tmp); 89817ab2063SBarry Smith } 89917ab2063SBarry Smith else if (type == NORM_INFINITY) { 90017ab2063SBarry Smith *norm = 0.0; 901416022c9SBarry Smith for ( j=0; j<a->m; j++ ) { 902416022c9SBarry Smith v = a->a + a->i[j] + shift; 90317ab2063SBarry Smith sum = 0.0; 904416022c9SBarry Smith for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 905cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v++; 90617ab2063SBarry Smith } 90717ab2063SBarry Smith if (sum > *norm) *norm = sum; 90817ab2063SBarry Smith } 90917ab2063SBarry Smith } 91017ab2063SBarry Smith else { 91148d91487SBarry Smith SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 91217ab2063SBarry Smith } 91317ab2063SBarry Smith return 0; 91417ab2063SBarry Smith } 91517ab2063SBarry Smith 916416022c9SBarry Smith static int MatTranspose_SeqAIJ(Mat A,Mat *B) 91717ab2063SBarry Smith { 918416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 919416022c9SBarry Smith Mat C; 920416022c9SBarry Smith int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 921416022c9SBarry Smith int shift = a->indexshift; 922d5d45c9bSBarry Smith Scalar *array = a->a; 92317ab2063SBarry Smith 924416022c9SBarry Smith if (!B && m != a->n) SETERRQ(1,"MatTranspose_SeqAIJ:Not for rectangular mat in place"); 9250452661fSBarry Smith col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 926cddf8d76SBarry Smith PetscMemzero(col,(1+a->n)*sizeof(int)); 92717ab2063SBarry Smith if (shift) { 92817ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 92917ab2063SBarry Smith } 93017ab2063SBarry Smith for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 931416022c9SBarry Smith ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 9320452661fSBarry Smith PetscFree(col); 93317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 93417ab2063SBarry Smith len = ai[i+1]-ai[i]; 935416022c9SBarry Smith ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 93617ab2063SBarry Smith array += len; aj += len; 93717ab2063SBarry Smith } 93817ab2063SBarry Smith if (shift) { 93917ab2063SBarry Smith for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 94017ab2063SBarry Smith } 94117ab2063SBarry Smith 942416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 943416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 94417ab2063SBarry Smith 945416022c9SBarry Smith if (B) { 946416022c9SBarry Smith *B = C; 94717ab2063SBarry Smith } else { 948416022c9SBarry Smith /* This isn't really an in-place transpose */ 9490452661fSBarry Smith PetscFree(a->a); 9500452661fSBarry Smith if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 9510452661fSBarry Smith if (a->diag) PetscFree(a->diag); 9520452661fSBarry Smith if (a->ilen) PetscFree(a->ilen); 9530452661fSBarry Smith if (a->imax) PetscFree(a->imax); 9540452661fSBarry Smith if (a->solve_work) PetscFree(a->solve_work); 9550452661fSBarry Smith PetscFree(a); 956416022c9SBarry Smith PetscMemcpy(A,C,sizeof(struct _Mat)); 9570452661fSBarry Smith PetscHeaderDestroy(C); 95817ab2063SBarry Smith } 95917ab2063SBarry Smith return 0; 96017ab2063SBarry Smith } 96117ab2063SBarry Smith 962f0b747eeSBarry Smith static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 96317ab2063SBarry Smith { 964416022c9SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 96517ab2063SBarry Smith Scalar *l,*r,x,*v; 966d5d45c9bSBarry Smith int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 96717ab2063SBarry Smith 968f0b747eeSBarry Smith if (!a->assembled) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Not for unassembled matrix"); 96917ab2063SBarry Smith if (ll) { 97017ab2063SBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 971f0b747eeSBarry Smith if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 972416022c9SBarry Smith v = a->a; 97317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 97417ab2063SBarry Smith x = l[i]; 975416022c9SBarry Smith M = a->i[i+1] - a->i[i]; 97617ab2063SBarry Smith for ( j=0; j<M; j++ ) { (*v++) *= x;} 97717ab2063SBarry Smith } 97817ab2063SBarry Smith } 97917ab2063SBarry Smith if (rr) { 98017ab2063SBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 981f0b747eeSBarry Smith if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 982416022c9SBarry Smith v = a->a; jj = a->j; 98317ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 98417ab2063SBarry Smith (*v++) *= r[*jj++ + shift]; 98517ab2063SBarry Smith } 98617ab2063SBarry Smith } 98717ab2063SBarry Smith return 0; 98817ab2063SBarry Smith } 98917ab2063SBarry Smith 990cddf8d76SBarry Smith static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 99117ab2063SBarry Smith { 992db02288aSLois Curfman McInnes Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 99302834360SBarry Smith int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 994a2744918SBarry Smith register int sum,lensi; 99502834360SBarry Smith int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 996a2744918SBarry Smith int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 997db02288aSLois Curfman McInnes Scalar *vwork,*a_new; 998416022c9SBarry Smith Mat C; 99917ab2063SBarry Smith 1000416022c9SBarry Smith if (!a->assembled) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:Not for unassembled matrix"); 100117ab2063SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 100217ab2063SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 100317ab2063SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 100417ab2063SBarry Smith 100502834360SBarry Smith if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { 100602834360SBarry Smith /* special case of contiguous rows */ 10070452661fSBarry Smith lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 100802834360SBarry Smith starts = lens + ncols; 100902834360SBarry Smith /* loop over new rows determining lens and starting points */ 101002834360SBarry Smith for (i=0; i<nrows; i++) { 1011a2744918SBarry Smith kstart = ai[irow[i]]+shift; 1012a2744918SBarry Smith kend = kstart + ailen[irow[i]]; 101302834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1014d8ced48eSBarry Smith if (aj[k]+shift >= first) { 101502834360SBarry Smith starts[i] = k; 101602834360SBarry Smith break; 101702834360SBarry Smith } 101802834360SBarry Smith } 1019a2744918SBarry Smith sum = 0; 102002834360SBarry Smith while (k < kend) { 1021d8ced48eSBarry Smith if (aj[k++]+shift >= first+ncols) break; 1022a2744918SBarry Smith sum++; 102302834360SBarry Smith } 1024a2744918SBarry Smith lens[i] = sum; 102502834360SBarry Smith } 102602834360SBarry Smith /* create submatrix */ 1027cddf8d76SBarry Smith if (scall == MAT_REUSE_MATRIX) { 102808480c60SBarry Smith int n_cols,n_rows; 102908480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 103008480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1031d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 103208480c60SBarry Smith C = *B; 103308480c60SBarry Smith } 103408480c60SBarry Smith else { 103502834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 103608480c60SBarry Smith } 1037db02288aSLois Curfman McInnes c = (Mat_SeqAIJ*) C->data; 1038db02288aSLois Curfman McInnes 103902834360SBarry Smith /* loop over rows inserting into submatrix */ 1040db02288aSLois Curfman McInnes a_new = c->a; 1041db02288aSLois Curfman McInnes j_new = c->j; 1042db02288aSLois Curfman McInnes i_new = c->i; 1043db02288aSLois Curfman McInnes i_new[0] = -shift; 104402834360SBarry Smith for (i=0; i<nrows; i++) { 1045a2744918SBarry Smith ii = starts[i]; 1046a2744918SBarry Smith lensi = lens[i]; 1047a2744918SBarry Smith for ( k=0; k<lensi; k++ ) { 1048a2744918SBarry Smith *j_new++ = aj[ii+k] - first; 104902834360SBarry Smith } 1050a2744918SBarry Smith PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1051a2744918SBarry Smith a_new += lensi; 1052a2744918SBarry Smith i_new[i+1] = i_new[i] + lensi; 1053a2744918SBarry Smith c->ilen[i] = lensi; 105402834360SBarry Smith } 10550452661fSBarry Smith PetscFree(lens); 105602834360SBarry Smith } 105702834360SBarry Smith else { 105802834360SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 10590452661fSBarry Smith smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 106002834360SBarry Smith ssmap = smap + shift; 10617eb43aa7SLois Curfman McInnes cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 106202834360SBarry Smith lens = cwork + ncols; 10630452661fSBarry Smith vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1064cddf8d76SBarry Smith PetscMemzero(smap,oldcols*sizeof(int)); 106517ab2063SBarry Smith for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 106602834360SBarry Smith /* determine lens of each row */ 106702834360SBarry Smith for (i=0; i<nrows; i++) { 1068d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 106902834360SBarry Smith kend = kstart + a->ilen[irow[i]]; 107002834360SBarry Smith lens[i] = 0; 107102834360SBarry Smith for ( k=kstart; k<kend; k++ ) { 1072d8ced48eSBarry Smith if (ssmap[aj[k]]) { 107302834360SBarry Smith lens[i]++; 107402834360SBarry Smith } 107502834360SBarry Smith } 107602834360SBarry Smith } 107717ab2063SBarry Smith /* Create and fill new matrix */ 1078a2744918SBarry Smith if (scall == MAT_REUSE_MATRIX) { 107908480c60SBarry Smith int n_cols,n_rows; 108008480c60SBarry Smith ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 108108480c60SBarry Smith if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1082d8ced48eSBarry Smith ierr = MatZeroEntries(*B); CHKERRQ(ierr); 108308480c60SBarry Smith C = *B; 108408480c60SBarry Smith } 108508480c60SBarry Smith else { 108602834360SBarry Smith ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 108708480c60SBarry Smith } 108817ab2063SBarry Smith for (i=0; i<nrows; i++) { 108917ab2063SBarry Smith nznew = 0; 1090d8ced48eSBarry Smith kstart = ai[irow[i]]+shift; 1091416022c9SBarry Smith kend = kstart + a->ilen[irow[i]]; 109217ab2063SBarry Smith for ( k=kstart; k<kend; k++ ) { 109302834360SBarry Smith if (ssmap[a->j[k]]) { 109402834360SBarry Smith cwork[nznew] = ssmap[a->j[k]] - 1; 1095416022c9SBarry Smith vwork[nznew++] = a->a[k]; 109617ab2063SBarry Smith } 109717ab2063SBarry Smith } 1098416022c9SBarry Smith ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 109917ab2063SBarry Smith } 110002834360SBarry Smith /* Free work space */ 110102834360SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 11020452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 110302834360SBarry Smith } 1104416022c9SBarry Smith ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1105416022c9SBarry Smith ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 110617ab2063SBarry Smith 110717ab2063SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1108416022c9SBarry Smith *B = C; 110917ab2063SBarry Smith return 0; 111017ab2063SBarry Smith } 111117ab2063SBarry Smith 1112a871dcd8SBarry Smith /* 111363b91edcSBarry Smith note: This can only work for identity for row and col. It would 111463b91edcSBarry Smith be good to check this and otherwise generate an error. 1115a871dcd8SBarry Smith */ 111663b91edcSBarry Smith static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1117a871dcd8SBarry Smith { 111863b91edcSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 111908480c60SBarry Smith int ierr; 112063b91edcSBarry Smith Mat outA; 112163b91edcSBarry Smith 1122a871dcd8SBarry Smith if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1123a871dcd8SBarry Smith 112463b91edcSBarry Smith outA = inA; 112563b91edcSBarry Smith inA->factor = FACTOR_LU; 112663b91edcSBarry Smith a->row = row; 112763b91edcSBarry Smith a->col = col; 112863b91edcSBarry Smith 11290452661fSBarry Smith a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 113063b91edcSBarry Smith 113108480c60SBarry Smith if (!a->diag) { 113208480c60SBarry Smith ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 113363b91edcSBarry Smith } 113463b91edcSBarry Smith ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1135a871dcd8SBarry Smith return 0; 1136a871dcd8SBarry Smith } 1137a871dcd8SBarry Smith 1138f0b747eeSBarry Smith #include "pinclude/plapack.h" 1139f0b747eeSBarry Smith static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1140f0b747eeSBarry Smith { 1141f0b747eeSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1142f0b747eeSBarry Smith int one = 1; 1143f0b747eeSBarry Smith if (!a->assembled) SETERRQ(1,"MatScale_SeqAIJ:Not for unassembled matrix"); 1144f0b747eeSBarry Smith BLscal_( &a->nz, alpha, a->a, &one ); 1145f0b747eeSBarry Smith PLogFlops(a->nz); 1146f0b747eeSBarry Smith return 0; 1147f0b747eeSBarry Smith } 1148f0b747eeSBarry Smith 1149cddf8d76SBarry Smith static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1150cddf8d76SBarry Smith Mat **B) 1151cddf8d76SBarry Smith { 1152cddf8d76SBarry Smith int ierr,i; 1153cddf8d76SBarry Smith 1154cddf8d76SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 11550452661fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1156cddf8d76SBarry Smith } 1157cddf8d76SBarry Smith 1158cddf8d76SBarry Smith for ( i=0; i<n; i++ ) { 1159cddf8d76SBarry Smith ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1160cddf8d76SBarry Smith } 1161cddf8d76SBarry Smith return 0; 1162cddf8d76SBarry Smith } 1163cddf8d76SBarry Smith 1164e4d965acSSatish Balay static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 11654dcbc457SBarry Smith { 1166e4d965acSSatish Balay Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 11678a047759SSatish Balay int shift, row, i,j,k,l,m,n, *idx,ierr, *table, *nidx, isz, val; 11688a047759SSatish Balay int start, end, *ai, *aj; 11694dcbc457SBarry Smith 11708a047759SSatish Balay shift = a->indexshift; 1171e4d965acSSatish Balay m = a->m; 1172e4d965acSSatish Balay ai = a->i; 11738a047759SSatish Balay aj = a->j+shift; 11748a047759SSatish Balay 1175*04a348a9SBarry Smith table = (int *) PetscMalloc(2*m*sizeof(int)); CHKPTRQ(table); nidx = table + m; 11768a047759SSatish Balay 1177e4d965acSSatish Balay if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 1178e4d965acSSatish Balay for ( i=0; i<is_max; i++ ) { 1179e4d965acSSatish Balay /* Initialise the two local arrays */ 1180e4d965acSSatish Balay isz = 0; 1181e4d965acSSatish Balay PetscMemzero(table,m*sizeof(int)); 1182e4d965acSSatish Balay 1183e4d965acSSatish Balay /* Extract the indices, assume there can be duplicate entries */ 11844dcbc457SBarry Smith ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1185e4d965acSSatish Balay ierr = ISGetLocalSize(is[i],&n); CHKERRQ(ierr); 1186e4d965acSSatish Balay 1187e4d965acSSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1188e4d965acSSatish Balay for ( j=0; j<n ; ++j){ 1189e4d965acSSatish Balay if(!table[idx[j]]++) { nidx[isz++] = idx[j];} 11904dcbc457SBarry Smith } 1191e4d965acSSatish Balay 1192*04a348a9SBarry Smith k = 0; 1193*04a348a9SBarry Smith for ( j=0; j<ov; j++){ /* for each overlap*/ 1194*04a348a9SBarry Smith n = isz; 1195*04a348a9SBarry Smith for ( ; k<n ; k++){ /* do only those rows in nidx[k] not yet done */ 1196e4d965acSSatish Balay row = nidx[k]; 1197e4d965acSSatish Balay start = ai[row]; 1198e4d965acSSatish Balay end = ai[row+1]; 1199*04a348a9SBarry Smith for ( l = start; l<end ; l++){ 12008a047759SSatish Balay val = aj[l] + shift; 1201e4d965acSSatish Balay if (!table[val]++) {nidx[isz++] = val;} 1202e4d965acSSatish Balay } 1203e4d965acSSatish Balay } 1204e4d965acSSatish Balay } 1205e4d965acSSatish Balay ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 1206e4d965acSSatish Balay ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1207e4d965acSSatish Balay ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1208e4d965acSSatish Balay } 1209*04a348a9SBarry Smith PetscFree(table); 1210e4d965acSSatish Balay return 0; 12114dcbc457SBarry Smith } 121217ab2063SBarry Smith 1213682d7d0cSBarry Smith int MatPrintHelp_SeqAIJ(Mat A) 1214682d7d0cSBarry Smith { 1215682d7d0cSBarry Smith static int called = 0; 1216682d7d0cSBarry Smith MPI_Comm comm = A->comm; 1217682d7d0cSBarry Smith 1218682d7d0cSBarry Smith if (called) return 0; else called = 1; 12192a7368beSLois Curfman McInnes MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 1220682d7d0cSBarry Smith MPIU_printf(comm," -mat_lu_pivotthreshold <threshold>\n"); 12212a7368beSLois Curfman McInnes MPIU_printf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 1222682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_no_inode - Do not use inodes\n"); 1223682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1224682d7d0cSBarry Smith #if defined(HAVE_ESSL) 1225682d7d0cSBarry Smith MPIU_printf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1226682d7d0cSBarry Smith #endif 1227682d7d0cSBarry Smith return 0; 1228682d7d0cSBarry Smith } 1229682d7d0cSBarry Smith /* -------------------------------------------------------------------*/ 123017ab2063SBarry Smith static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 123117ab2063SBarry Smith MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1232416022c9SBarry Smith MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1233416022c9SBarry Smith MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 123417ab2063SBarry Smith MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 123517ab2063SBarry Smith MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 123617ab2063SBarry Smith MatLUFactor_SeqAIJ,0, 123717ab2063SBarry Smith MatRelax_SeqAIJ, 123817ab2063SBarry Smith MatTranspose_SeqAIJ, 123917ab2063SBarry Smith MatGetInfo_SeqAIJ,0, 1240f0b747eeSBarry Smith MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 124117ab2063SBarry Smith 0,MatAssemblyEnd_SeqAIJ, 124217ab2063SBarry Smith MatCompress_SeqAIJ, 124317ab2063SBarry Smith MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 124417ab2063SBarry Smith MatGetReordering_SeqAIJ, 124517ab2063SBarry Smith MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 124617ab2063SBarry Smith MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 124717ab2063SBarry Smith MatILUFactorSymbolic_SeqAIJ,0, 124817ab2063SBarry Smith 0,0,MatConvert_SeqAIJ, 1249416022c9SBarry Smith MatGetSubMatrix_SeqAIJ,0, 12503d1612f7SBarry Smith MatConvertSameType_SeqAIJ,0,0, 1251cddf8d76SBarry Smith MatILUFactor_SeqAIJ,0,0, 12527eb43aa7SLois Curfman McInnes MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1253682d7d0cSBarry Smith MatGetValues_SeqAIJ,0, 1254f0b747eeSBarry Smith MatPrintHelp_SeqAIJ, 1255f0b747eeSBarry Smith MatScale_SeqAIJ}; 125617ab2063SBarry Smith 125717ab2063SBarry Smith extern int MatUseSuperLU_SeqAIJ(Mat); 125817ab2063SBarry Smith extern int MatUseEssl_SeqAIJ(Mat); 125917ab2063SBarry Smith extern int MatUseDXML_SeqAIJ(Mat); 126017ab2063SBarry Smith 126117ab2063SBarry Smith /*@C 1262682d7d0cSBarry Smith MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 126317ab2063SBarry Smith (the default uniprocessor PETSc format). 126417ab2063SBarry Smith 126517ab2063SBarry Smith Input Parameters: 126617ab2063SBarry Smith . comm - MPI communicator, set to MPI_COMM_SELF 126717ab2063SBarry Smith . m - number of rows 126817ab2063SBarry Smith . n - number of columns 126917ab2063SBarry Smith . nz - number of nonzeros per row (same for all rows) 127017ab2063SBarry Smith . nzz - number of nonzeros per row or null (possibly different for each row) 127117ab2063SBarry Smith 127217ab2063SBarry Smith Output Parameter: 1273416022c9SBarry Smith . A - the matrix 127417ab2063SBarry Smith 127517ab2063SBarry Smith Notes: 127617ab2063SBarry Smith The AIJ format (also called the Yale sparse matrix format or 127717ab2063SBarry Smith compressed row storage), is fully compatible with standard Fortran 77 12780002213bSLois Curfman McInnes storage. That is, the stored row and column indices can begin at 12790002213bSLois Curfman McInnes either one (as in Fortran) or zero. See the users manual for details. 128017ab2063SBarry Smith 128117ab2063SBarry Smith Specify the preallocated storage with either nz or nnz (not both). 1282a40aa06bSLois Curfman McInnes Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 128317ab2063SBarry Smith allocation. 128417ab2063SBarry Smith 1285682d7d0cSBarry Smith By default, this format uses inodes (identical nodes) when possible, to 1286682d7d0cSBarry Smith improve numerical efficiency of Matrix vector products and solves. We 1287682d7d0cSBarry Smith search for consecutive rows with the same nonzero structure, thereby 12886c7ebb05SLois Curfman McInnes reusing matrix information to achieve increased efficiency. 12896c7ebb05SLois Curfman McInnes 12906c7ebb05SLois Curfman McInnes Options Database Keys: 12916c7ebb05SLois Curfman McInnes $ -mat_aij_no_inode - Do not use inodes 12926c7ebb05SLois Curfman McInnes $ -mat_aij_inode_limit <limit> - Set inode limit (max limit=5) 12934b0e389bSBarry Smith $ -mat_aij_oneindex - Internally using indexing starting at 1 rather then zero. 12944b0e389bSBarry Smith . Note: You still 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; 131369957df2SSatish Balay ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, \ 131469957df2SSatish Balay &flg); CHKERRQ(ierr); 1315416022c9SBarry Smith b->row = 0; 1316416022c9SBarry Smith b->col = 0; 1317416022c9SBarry Smith b->indexshift = 0; 131869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 131969957df2SSatish Balay if (flg) b->indexshift = -1; 132017ab2063SBarry Smith 1321416022c9SBarry Smith b->m = m; 1322416022c9SBarry Smith b->n = n; 13230452661fSBarry Smith b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1324b4fd4287SBarry Smith if (nnz == PETSC_NULL) { 13257b8455f0SLois Curfman McInnes if (nz == PETSC_DEFAULT) nz = 10; 13267b8455f0SLois Curfman McInnes else if (nz <= 0) nz = 1; 1327416022c9SBarry Smith for ( i=0; i<m; i++ ) b->imax[i] = nz; 132817ab2063SBarry Smith nz = nz*m; 132917ab2063SBarry Smith } 133017ab2063SBarry Smith else { 133117ab2063SBarry Smith nz = 0; 1332416022c9SBarry Smith for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 133317ab2063SBarry Smith } 133417ab2063SBarry Smith 133517ab2063SBarry Smith /* allocate the matrix space */ 1336416022c9SBarry Smith len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 13370452661fSBarry Smith b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1338416022c9SBarry Smith b->j = (int *) (b->a + nz); 1339cddf8d76SBarry Smith PetscMemzero(b->j,nz*sizeof(int)); 1340416022c9SBarry Smith b->i = b->j + nz; 1341416022c9SBarry Smith b->singlemalloc = 1; 134217ab2063SBarry Smith 1343416022c9SBarry Smith b->i[0] = -b->indexshift; 134417ab2063SBarry Smith for (i=1; i<m+1; i++) { 1345416022c9SBarry Smith b->i[i] = b->i[i-1] + b->imax[i-1]; 134617ab2063SBarry Smith } 134717ab2063SBarry Smith 1348416022c9SBarry Smith /* b->ilen will count nonzeros in each row so far. */ 13490452661fSBarry Smith b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1350416022c9SBarry Smith PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1351416022c9SBarry Smith for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 135217ab2063SBarry Smith 1353416022c9SBarry Smith b->nz = 0; 1354416022c9SBarry Smith b->maxnz = nz; 1355416022c9SBarry Smith b->sorted = 0; 1356416022c9SBarry Smith b->roworiented = 1; 1357416022c9SBarry Smith b->nonew = 0; 1358416022c9SBarry Smith b->diag = 0; 1359416022c9SBarry Smith b->assembled = 0; 1360416022c9SBarry Smith b->solve_work = 0; 136171bd300dSLois Curfman McInnes b->spptr = 0; 1362754ec7b1SSatish Balay b->inode.node_count = 0; 1363754ec7b1SSatish Balay b->inode.size = 0; 13646c7ebb05SLois Curfman McInnes b->inode.limit = 5; 13656c7ebb05SLois Curfman McInnes b->inode.max_limit = 5; 136617ab2063SBarry Smith 1367416022c9SBarry Smith *A = B; 136869957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 136969957df2SSatish Balay if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 137069957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 137169957df2SSatish Balay if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 137269957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 137369957df2SSatish Balay if (flg) { 1374416022c9SBarry Smith if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1375416022c9SBarry Smith ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 137617ab2063SBarry Smith } 137769957df2SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 137869957df2SSatish Balay if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 137917ab2063SBarry Smith return 0; 138017ab2063SBarry Smith } 138117ab2063SBarry Smith 13823d1612f7SBarry Smith int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 138317ab2063SBarry Smith { 1384416022c9SBarry Smith Mat C; 1385416022c9SBarry Smith Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 138608480c60SBarry Smith int i,len, m = a->m,shift = a->indexshift; 138717ab2063SBarry Smith 13884043dd9cSLois Curfman McInnes *B = 0; 1389f0b747eeSBarry Smith if (!a->assembled) SETERRQ(1,"MatConvertSameType_SeqAIJ:Not for unassembled matrix"); 13900452661fSBarry Smith PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1391416022c9SBarry Smith PLogObjectCreate(C); 13920452661fSBarry Smith C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 139341c01911SSatish Balay PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1394416022c9SBarry Smith C->destroy = MatDestroy_SeqAIJ; 1395416022c9SBarry Smith C->view = MatView_SeqAIJ; 1396416022c9SBarry Smith C->factor = A->factor; 1397416022c9SBarry Smith c->row = 0; 1398416022c9SBarry Smith c->col = 0; 1399416022c9SBarry Smith c->indexshift = shift; 140017ab2063SBarry Smith 1401416022c9SBarry Smith c->m = a->m; 1402416022c9SBarry Smith c->n = a->n; 140317ab2063SBarry Smith 14040452661fSBarry Smith c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 14050452661fSBarry Smith c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 140617ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1407416022c9SBarry Smith c->imax[i] = a->imax[i]; 1408416022c9SBarry Smith c->ilen[i] = a->ilen[i]; 140917ab2063SBarry Smith } 141017ab2063SBarry Smith 141117ab2063SBarry Smith /* allocate the matrix space */ 1412416022c9SBarry Smith c->singlemalloc = 1; 1413416022c9SBarry Smith len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 14140452661fSBarry Smith c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1415416022c9SBarry Smith c->j = (int *) (c->a + a->i[m] + shift); 1416416022c9SBarry Smith c->i = c->j + a->i[m] + shift; 1417416022c9SBarry Smith PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 141817ab2063SBarry Smith if (m > 0) { 1419416022c9SBarry Smith PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 142008480c60SBarry Smith if (cpvalues == COPY_VALUES) { 1421416022c9SBarry Smith PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 142217ab2063SBarry Smith } 142308480c60SBarry Smith } 142417ab2063SBarry Smith 1425416022c9SBarry Smith PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1426416022c9SBarry Smith c->sorted = a->sorted; 1427416022c9SBarry Smith c->roworiented = a->roworiented; 1428416022c9SBarry Smith c->nonew = a->nonew; 142917ab2063SBarry Smith 1430416022c9SBarry Smith if (a->diag) { 14310452661fSBarry Smith c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1432416022c9SBarry Smith PLogObjectMemory(C,(m+1)*sizeof(int)); 143317ab2063SBarry Smith for ( i=0; i<m; i++ ) { 1434416022c9SBarry Smith c->diag[i] = a->diag[i]; 143517ab2063SBarry Smith } 143617ab2063SBarry Smith } 1437416022c9SBarry Smith else c->diag = 0; 14386c7ebb05SLois Curfman McInnes c->inode.limit = a->inode.limit; 14396c7ebb05SLois Curfman McInnes c->inode.max_limit = a->inode.max_limit; 1440754ec7b1SSatish Balay if (a->inode.size){ 1441754ec7b1SSatish Balay c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1442754ec7b1SSatish Balay c->inode.node_count = a->inode.node_count; 1443754ec7b1SSatish Balay PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1444754ec7b1SSatish Balay } else { 1445754ec7b1SSatish Balay c->inode.size = 0; 1446754ec7b1SSatish Balay c->inode.node_count = 0; 1447754ec7b1SSatish Balay } 1448416022c9SBarry Smith c->assembled = 1; 1449416022c9SBarry Smith c->nz = a->nz; 1450416022c9SBarry Smith c->maxnz = a->maxnz; 1451416022c9SBarry Smith c->solve_work = 0; 145276dd722bSSatish Balay c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1453754ec7b1SSatish Balay 1454416022c9SBarry Smith *B = C; 145517ab2063SBarry Smith return 0; 145617ab2063SBarry Smith } 145717ab2063SBarry Smith 1458416022c9SBarry Smith int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 145917ab2063SBarry Smith { 1460416022c9SBarry Smith Mat_SeqAIJ *a; 1461416022c9SBarry Smith Mat B; 146217699dbbSLois Curfman McInnes int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 146317ab2063SBarry Smith PetscObject vobj = (PetscObject) bview; 146417ab2063SBarry Smith MPI_Comm comm = vobj->comm; 146517ab2063SBarry Smith 146617699dbbSLois Curfman McInnes MPI_Comm_size(comm,&size); 146717699dbbSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 146817ab2063SBarry Smith ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1469416022c9SBarry Smith ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 147048d91487SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 147117ab2063SBarry Smith M = header[1]; N = header[2]; nz = header[3]; 147217ab2063SBarry Smith 147317ab2063SBarry Smith /* read in row lengths */ 14740452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1475416022c9SBarry Smith ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 147617ab2063SBarry Smith 147717ab2063SBarry Smith /* create our matrix */ 1478416022c9SBarry Smith ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1479416022c9SBarry Smith B = *A; 1480416022c9SBarry Smith a = (Mat_SeqAIJ *) B->data; 1481416022c9SBarry Smith shift = a->indexshift; 148217ab2063SBarry Smith 148317ab2063SBarry Smith /* read in column indices and adjust for Fortran indexing*/ 1484416022c9SBarry Smith ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 148517ab2063SBarry Smith if (shift) { 148617ab2063SBarry Smith for ( i=0; i<nz; i++ ) { 1487416022c9SBarry Smith a->j[i] += 1; 148817ab2063SBarry Smith } 148917ab2063SBarry Smith } 149017ab2063SBarry Smith 149117ab2063SBarry Smith /* read in nonzero values */ 1492416022c9SBarry Smith ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 149317ab2063SBarry Smith 149417ab2063SBarry Smith /* set matrix "i" values */ 1495416022c9SBarry Smith a->i[0] = -shift; 149617ab2063SBarry Smith for ( i=1; i<= M; i++ ) { 1497416022c9SBarry Smith a->i[i] = a->i[i-1] + rowlengths[i-1]; 1498416022c9SBarry Smith a->ilen[i-1] = rowlengths[i-1]; 149917ab2063SBarry Smith } 15000452661fSBarry Smith PetscFree(rowlengths); 150117ab2063SBarry Smith 1502416022c9SBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1503416022c9SBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 150417ab2063SBarry Smith return 0; 150517ab2063SBarry Smith } 150617ab2063SBarry Smith 150717ab2063SBarry Smith 150817ab2063SBarry Smith 1509