xref: /petsc/src/mat/impls/baij/seq/baij.c (revision af5da2bf8042f6414d2176c885c7f2028a1539d5)
12593348eSBarry Smith #ifndef lint
2*af5da2bfSBarry Smith static char vcid[] = "$Id: baij.c,v 1.69 1996/11/07 15:09:41 bsmith Exp bsmith $";
32593348eSBarry Smith #endif
42593348eSBarry Smith 
52593348eSBarry Smith /*
6b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
72593348eSBarry Smith   matrix storage format.
82593348eSBarry Smith */
970f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
101a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
111a6a6d98SBarry Smith #include "src/inline/spops.h"
122593348eSBarry Smith #include "petsc.h"
133b547af2SSatish Balay 
142593348eSBarry Smith 
15de6a44a3SBarry Smith /*
16de6a44a3SBarry Smith      Adds diagonal pointers to sparse matrix structure.
17de6a44a3SBarry Smith */
18de6a44a3SBarry Smith 
19de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A)
20de6a44a3SBarry Smith {
21de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
227fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
23de6a44a3SBarry Smith 
24de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
25de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
267fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
27de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
28de6a44a3SBarry Smith       if (a->j[j] == i) {
29de6a44a3SBarry Smith         diag[i] = j;
30de6a44a3SBarry Smith         break;
31de6a44a3SBarry Smith       }
32de6a44a3SBarry Smith     }
33de6a44a3SBarry Smith   }
34de6a44a3SBarry Smith   a->diag = diag;
35de6a44a3SBarry Smith   return 0;
36de6a44a3SBarry Smith }
372593348eSBarry Smith 
382593348eSBarry Smith #include "draw.h"
392593348eSBarry Smith #include "pinclude/pviewer.h"
4077c4ece6SBarry Smith #include "sys.h"
412593348eSBarry Smith 
423b2fbd54SBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
433b2fbd54SBarry Smith 
443b2fbd54SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
453b2fbd54SBarry Smith                             PetscTruth *done)
463b2fbd54SBarry Smith {
473b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
483b2fbd54SBarry Smith   int         ierr, n = a->mbs,i;
493b2fbd54SBarry Smith 
503b2fbd54SBarry Smith   *nn = n;
513b2fbd54SBarry Smith   if (!ia) return 0;
523b2fbd54SBarry Smith   if (symmetric) {
533b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
543b2fbd54SBarry Smith   } else if (oshift == 1) {
553b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
563b2fbd54SBarry Smith     int nz = a->i[n] + 1;
573b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]++;
583b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]++;
593b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
603b2fbd54SBarry Smith   } else {
613b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
623b2fbd54SBarry Smith   }
633b2fbd54SBarry Smith 
643b2fbd54SBarry Smith   return 0;
653b2fbd54SBarry Smith }
663b2fbd54SBarry Smith 
673b2fbd54SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
683b2fbd54SBarry Smith                                 PetscTruth *done)
693b2fbd54SBarry Smith {
703b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
713b2fbd54SBarry Smith   int         i,n = a->mbs;
723b2fbd54SBarry Smith 
733b2fbd54SBarry Smith   if (!ia) return 0;
743b2fbd54SBarry Smith   if (symmetric) {
753b2fbd54SBarry Smith     PetscFree(*ia);
763b2fbd54SBarry Smith     PetscFree(*ja);
77*af5da2bfSBarry Smith   } else if (oshift == 1) {
783b2fbd54SBarry Smith     int nz = a->i[n];
793b2fbd54SBarry Smith     for ( i=0; i<nz; i++ ) a->j[i]--;
803b2fbd54SBarry Smith     for ( i=0; i<n+1; i++ ) a->i[i]--;
813b2fbd54SBarry Smith   }
823b2fbd54SBarry Smith   return 0;
833b2fbd54SBarry Smith }
843b2fbd54SBarry Smith 
853b2fbd54SBarry Smith 
86b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
872593348eSBarry Smith {
88b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
899df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
90b6490206SBarry Smith   Scalar      *aa;
912593348eSBarry Smith 
9290ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
932593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
942593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
953b2fbd54SBarry Smith 
962593348eSBarry Smith   col_lens[1] = a->m;
972593348eSBarry Smith   col_lens[2] = a->n;
987e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
992593348eSBarry Smith 
1002593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
101b6490206SBarry Smith   count = 0;
102b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
103b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
104b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
105b6490206SBarry Smith     }
1062593348eSBarry Smith   }
10777c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
1082593348eSBarry Smith   PetscFree(col_lens);
1092593348eSBarry Smith 
1102593348eSBarry Smith   /* store column indices (zero start index) */
1117e67e3f9SSatish Balay   jj = (int *) PetscMalloc( a->nz*bs2*sizeof(int) ); CHKPTRQ(jj);
112b6490206SBarry Smith   count = 0;
113b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
114b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
115b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
116b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
117b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
1182593348eSBarry Smith         }
1192593348eSBarry Smith       }
120b6490206SBarry Smith     }
121b6490206SBarry Smith   }
1227e67e3f9SSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,BINARY_INT,0); CHKERRQ(ierr);
123b6490206SBarry Smith   PetscFree(jj);
1242593348eSBarry Smith 
1252593348eSBarry Smith   /* store nonzero values */
1267e67e3f9SSatish Balay   aa = (Scalar *) PetscMalloc(a->nz*bs2*sizeof(Scalar)); CHKPTRQ(aa);
127b6490206SBarry Smith   count = 0;
128b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
129b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
130b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
131b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
1327e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
133b6490206SBarry Smith         }
134b6490206SBarry Smith       }
135b6490206SBarry Smith     }
136b6490206SBarry Smith   }
1377e67e3f9SSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
138b6490206SBarry Smith   PetscFree(aa);
1392593348eSBarry Smith   return 0;
1402593348eSBarry Smith }
1412593348eSBarry Smith 
142b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
1432593348eSBarry Smith {
144b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1459df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
1462593348eSBarry Smith   FILE        *fd;
1472593348eSBarry Smith   char        *outputname;
1482593348eSBarry Smith 
14990ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
1502593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
15190ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
152639f9d9dSBarry Smith   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
1537ddc982cSLois Curfman McInnes     fprintf(fd,"  block size is %d\n",bs);
1542593348eSBarry Smith   }
155639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
156b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported");
1572593348eSBarry Smith   }
158639f9d9dSBarry Smith   else if (format == VIEWER_FORMAT_ASCII_COMMON) {
15944cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
16044cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
16144cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
16244cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
16344cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
16444cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
16544cd7ae7SLois Curfman McInnes           if (imag(a->a[bs2*k + l*bs + j]) != 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
16644cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
16744cd7ae7SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
16844cd7ae7SLois Curfman McInnes           else if (real(a->a[bs2*k + l*bs + j]) != 0.0)
16944cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
17044cd7ae7SLois Curfman McInnes #else
17144cd7ae7SLois Curfman McInnes           if (a->a[bs2*k + l*bs + j] != 0.0)
17244cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
17344cd7ae7SLois Curfman McInnes #endif
17444cd7ae7SLois Curfman McInnes           }
17544cd7ae7SLois Curfman McInnes         }
17644cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
17744cd7ae7SLois Curfman McInnes       }
17844cd7ae7SLois Curfman McInnes     }
17944cd7ae7SLois Curfman McInnes   }
1802593348eSBarry Smith   else {
181b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
182b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
183b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
184b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
185b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
18688685aaeSLois Curfman McInnes #if defined(PETSC_COMPLEX)
1877e67e3f9SSatish Balay           if (imag(a->a[bs2*k + l*bs + j]) != 0.0) {
18888685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
1897e67e3f9SSatish Balay               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
19088685aaeSLois Curfman McInnes           }
19188685aaeSLois Curfman McInnes           else {
1927e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
19388685aaeSLois Curfman McInnes           }
19488685aaeSLois Curfman McInnes #else
1957e67e3f9SSatish Balay             fprintf(fd," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
19688685aaeSLois Curfman McInnes #endif
1972593348eSBarry Smith           }
1982593348eSBarry Smith         }
1992593348eSBarry Smith         fprintf(fd,"\n");
2002593348eSBarry Smith       }
2012593348eSBarry Smith     }
202b6490206SBarry Smith   }
2032593348eSBarry Smith   fflush(fd);
2042593348eSBarry Smith   return 0;
2052593348eSBarry Smith }
2062593348eSBarry Smith 
2073270192aSSatish Balay static int MatView_SeqBAIJ_Draw(Mat A,Viewer viewer)
2083270192aSSatish Balay {
2093270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ *) A->data;
2103270192aSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,pause,color,bs=a->bs,bs2=a->bs2;
2113270192aSSatish Balay   double       xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r;
2123270192aSSatish Balay   Scalar       *aa;
2133270192aSSatish Balay   Draw         draw;
2143270192aSSatish Balay   DrawButton   button;
2153270192aSSatish Balay   PetscTruth   isnull;
2163270192aSSatish Balay 
2173270192aSSatish Balay   ViewerDrawGetDraw(viewer,&draw);
2183270192aSSatish Balay   ierr = DrawIsNull(draw,&isnull); CHKERRQ(ierr); if (isnull) return 0;
2193270192aSSatish Balay 
2203270192aSSatish Balay   xr  = a->n; yr = a->m; h = yr/10.0; w = xr/10.0;
2213270192aSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
2223270192aSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
2233270192aSSatish Balay   /* loop over matrix elements drawing boxes */
2243270192aSSatish Balay   color = DRAW_BLUE;
2253270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
2263270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
2273270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
2283270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
2293270192aSSatish Balay       aa = a->a + j*bs2;
2303270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
2313270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
2323270192aSSatish Balay           if (PetscReal(*aa++) >=  0.) continue;
2333270192aSSatish Balay           DrawRectangle(draw,x_l+k,y_l+l,x_r+k,y_r+l,color,color,color,color);
2343270192aSSatish Balay         }
2353270192aSSatish Balay       }
2363270192aSSatish Balay     }
2373270192aSSatish Balay   }
2383270192aSSatish Balay   color = DRAW_CYAN;
2393270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
2403270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
2413270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
2423270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
2433270192aSSatish Balay       aa = a->a + j*bs2;
2443270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
2453270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
2463270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
2473270192aSSatish Balay           DrawRectangle(draw,x_l+k,y_l+l,x_r+k,y_r+l,color,color,color,color);
2483270192aSSatish Balay         }
2493270192aSSatish Balay       }
2503270192aSSatish Balay     }
2513270192aSSatish Balay   }
2523270192aSSatish Balay 
2533270192aSSatish Balay   color = DRAW_RED;
2543270192aSSatish Balay   for ( i=0,row=0; i<mbs; i++,row+=bs ) {
2553270192aSSatish Balay     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
2563270192aSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
2573270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
2583270192aSSatish Balay       aa = a->a + j*bs2;
2593270192aSSatish Balay       for ( k=0; k<bs; k++ ) {
2603270192aSSatish Balay         for ( l=0; l<bs; l++ ) {
2613270192aSSatish Balay           if (PetscReal(*aa++) <= 0.) continue;
2623270192aSSatish Balay           DrawRectangle(draw,x_l+k,y_l+l,x_r+k,y_r+l,color,color,color,color);
2633270192aSSatish Balay         }
2643270192aSSatish Balay       }
2653270192aSSatish Balay     }
2663270192aSSatish Balay   }
2673270192aSSatish Balay 
2683270192aSSatish Balay   DrawFlush(draw);
2693270192aSSatish Balay   DrawGetPause(draw,&pause);
2703270192aSSatish Balay   if (pause >= 0) { PetscSleep(pause); return 0;}
2713270192aSSatish Balay 
2723270192aSSatish Balay   /* allow the matrix to zoom or shrink */
2733270192aSSatish Balay   ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
2743270192aSSatish Balay   while (button != BUTTON_RIGHT) {
2753270192aSSatish Balay     DrawClear(draw);
2763270192aSSatish Balay     if (button == BUTTON_LEFT) scale = .5;
2773270192aSSatish Balay     else if (button == BUTTON_CENTER) scale = 2.;
2783270192aSSatish Balay     xl = scale*(xl + w - xc) + xc - w*scale;
2793270192aSSatish Balay     xr = scale*(xr - w - xc) + xc + w*scale;
2803270192aSSatish Balay     yl = scale*(yl + h - yc) + yc - h*scale;
2813270192aSSatish Balay     yr = scale*(yr - h - yc) + yc + h*scale;
2823270192aSSatish Balay     w *= scale; h *= scale;
2833270192aSSatish Balay     ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr);
2843270192aSSatish Balay     color = DRAW_BLUE;
2853270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
2863270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
2873270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
2883270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
2893270192aSSatish Balay         aa = a->a + j*bs2;
2903270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
2913270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
2923270192aSSatish Balay             if (PetscReal(*aa++) >=  0.) continue;
2933270192aSSatish Balay             DrawRectangle(draw,x_l+k,y_l+l,x_r+k,y_r+l,color,color,color,color);
2943270192aSSatish Balay           }
2953270192aSSatish Balay         }
2963270192aSSatish Balay       }
2973270192aSSatish Balay     }
2983270192aSSatish Balay     color = DRAW_CYAN;
2993270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
3003270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
3013270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
3023270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
3033270192aSSatish Balay         aa = a->a + j*bs2;
3043270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
3053270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
3063270192aSSatish Balay           if (PetscReal(*aa++) != 0.) continue;
3073270192aSSatish Balay           DrawRectangle(draw,x_l+k,y_l+l,x_r+k,y_r+l,color,color,color,color);
3083270192aSSatish Balay           }
3093270192aSSatish Balay         }
3103270192aSSatish Balay       }
3113270192aSSatish Balay     }
3123270192aSSatish Balay 
3133270192aSSatish Balay     color = DRAW_RED;
3143270192aSSatish Balay     for ( i=0,row=0; i<mbs; i++,row+=bs ) {
3153270192aSSatish Balay       for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
3163270192aSSatish Balay         y_l = a->m - row - 1.0; y_r = y_l + 1.0;
3173270192aSSatish Balay         x_l = a->j[j]*bs; x_r = x_l + 1.0;
3183270192aSSatish Balay         aa = a->a + j*bs2;
3193270192aSSatish Balay         for ( k=0; k<bs; k++ ) {
3203270192aSSatish Balay           for ( l=0; l<bs; l++ ) {
3213270192aSSatish Balay             if (PetscReal(*aa++) <= 0.) continue;
3223270192aSSatish Balay             DrawRectangle(draw,x_l+k,y_l+l,x_r+k,y_r+l,color,color,color,color);
3233270192aSSatish Balay           }
3243270192aSSatish Balay         }
3253270192aSSatish Balay       }
3263270192aSSatish Balay     }
3273270192aSSatish Balay     ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0);
3283270192aSSatish Balay   }
3293270192aSSatish Balay   return 0;
3303270192aSSatish Balay }
3313270192aSSatish Balay 
332b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
3332593348eSBarry Smith {
3342593348eSBarry Smith   Mat         A = (Mat) obj;
33519bcc07fSBarry Smith   ViewerType  vtype;
33619bcc07fSBarry Smith   int         ierr;
3372593348eSBarry Smith 
33819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
33919bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
340b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported");
3412593348eSBarry Smith   }
34219bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
343b6490206SBarry Smith     return MatView_SeqBAIJ_ASCII(A,viewer);
3442593348eSBarry Smith   }
34519bcc07fSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
346b6490206SBarry Smith     return MatView_SeqBAIJ_Binary(A,viewer);
3472593348eSBarry Smith   }
34819bcc07fSBarry Smith   else if (vtype == DRAW_VIEWER) {
3493270192aSSatish Balay     return MatView_SeqBAIJ_Draw(A,viewer);
3502593348eSBarry Smith   }
3512593348eSBarry Smith   return 0;
3522593348eSBarry Smith }
353b6490206SBarry Smith 
354119a934aSSatish Balay #define CHUNKSIZE  10
355cd0e1443SSatish Balay 
356cd0e1443SSatish Balay /* This version has row oriented v  */
357639f9d9dSBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
358cd0e1443SSatish Balay {
359cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
360cd0e1443SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
361cd0e1443SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
362a7c10996SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
3639df24120SSatish Balay   int          ridx,cidx,bs2=a->bs2;
364cd0e1443SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
365cd0e1443SSatish Balay 
366cd0e1443SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
367cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
3683b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
369cd0e1443SSatish Balay     if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row");
370cd0e1443SSatish Balay     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large");
3713b2fbd54SBarry Smith #endif
372a7c10996SSatish Balay     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
373cd0e1443SSatish Balay     rmax = imax[brow]; nrow = ailen[brow];
374cd0e1443SSatish Balay     low = 0;
375cd0e1443SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
3763b2fbd54SBarry Smith #if defined(PETSC_BOPT_g)
377cd0e1443SSatish Balay       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column");
378cd0e1443SSatish Balay       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large");
3793b2fbd54SBarry Smith #endif
380a7c10996SSatish Balay       col = in[l]; bcol = col/bs;
381cd0e1443SSatish Balay       ridx = row % bs; cidx = col % bs;
382cd0e1443SSatish Balay       if (roworiented) {
383cd0e1443SSatish Balay         value = *v++;
384cd0e1443SSatish Balay       }
385cd0e1443SSatish Balay       else {
386cd0e1443SSatish Balay         value = v[k + l*m];
387cd0e1443SSatish Balay       }
388cd0e1443SSatish Balay       if (!sorted) low = 0; high = nrow;
389cd0e1443SSatish Balay       while (high-low > 5) {
390cd0e1443SSatish Balay         t = (low+high)/2;
391cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
392cd0e1443SSatish Balay         else              low  = t;
393cd0e1443SSatish Balay       }
394cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
395cd0e1443SSatish Balay         if (rp[i] > bcol) break;
396cd0e1443SSatish Balay         if (rp[i] == bcol) {
3977e67e3f9SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
398cd0e1443SSatish Balay           if (is == ADD_VALUES) *bap += value;
399cd0e1443SSatish Balay           else                  *bap  = value;
400cd0e1443SSatish Balay           goto noinsert;
401cd0e1443SSatish Balay         }
402cd0e1443SSatish Balay       }
403cd0e1443SSatish Balay       if (nonew) goto noinsert;
404cd0e1443SSatish Balay       if (nrow >= rmax) {
405cd0e1443SSatish Balay         /* there is no extra room in row, therefore enlarge */
406cd0e1443SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
407cd0e1443SSatish Balay         Scalar *new_a;
408cd0e1443SSatish Balay 
409cd0e1443SSatish Balay         /* malloc new storage space */
4107e67e3f9SSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
411cd0e1443SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
4127e67e3f9SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
413cd0e1443SSatish Balay         new_i   = new_j + new_nz;
414cd0e1443SSatish Balay 
415cd0e1443SSatish Balay         /* copy over old data into new slots */
416cd0e1443SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
4177e67e3f9SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
418a7c10996SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
419a7c10996SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
420a7c10996SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
421cd0e1443SSatish Balay                                                            len*sizeof(int));
422a7c10996SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar));
423a7c10996SSatish Balay         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
424a7c10996SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
425a7c10996SSatish Balay                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar));
426cd0e1443SSatish Balay         /* free up old matrix storage */
427cd0e1443SSatish Balay         PetscFree(a->a);
428cd0e1443SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
429cd0e1443SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
430cd0e1443SSatish Balay         a->singlemalloc = 1;
431cd0e1443SSatish Balay 
432a7c10996SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
433cd0e1443SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
4347e67e3f9SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
43518e7b8e4SLois Curfman McInnes         a->maxnz += bs2*CHUNKSIZE;
436cd0e1443SSatish Balay         a->reallocs++;
437119a934aSSatish Balay         a->nz++;
438cd0e1443SSatish Balay       }
4397e67e3f9SSatish Balay       N = nrow++ - 1;
440cd0e1443SSatish Balay       /* shift up all the later entries in this row */
441cd0e1443SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
442cd0e1443SSatish Balay         rp[ii+1] = rp[ii];
4437e67e3f9SSatish Balay         PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
444cd0e1443SSatish Balay       }
4457e67e3f9SSatish Balay       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
446cd0e1443SSatish Balay       rp[i]                      = bcol;
4477e67e3f9SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
448cd0e1443SSatish Balay       noinsert:;
449cd0e1443SSatish Balay       low = i;
450cd0e1443SSatish Balay     }
451cd0e1443SSatish Balay     ailen[brow] = nrow;
452cd0e1443SSatish Balay   }
453cd0e1443SSatish Balay   return 0;
454cd0e1443SSatish Balay }
455cd0e1443SSatish Balay 
4560b824a48SSatish Balay static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
4570b824a48SSatish Balay {
4580b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4590b824a48SSatish Balay   *m = a->m; *n = a->n;
4600b824a48SSatish Balay   return 0;
4610b824a48SSatish Balay }
4620b824a48SSatish Balay 
4630b824a48SSatish Balay static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
4640b824a48SSatish Balay {
4650b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4660b824a48SSatish Balay   *m = 0; *n = a->m;
4670b824a48SSatish Balay   return 0;
4680b824a48SSatish Balay }
4690b824a48SSatish Balay 
4709867e207SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
4719867e207SSatish Balay {
4729867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4737e67e3f9SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
4749867e207SSatish Balay   Scalar      *aa,*v_i,*aa_i;
4759867e207SSatish Balay 
4769867e207SSatish Balay   bs  = a->bs;
4779867e207SSatish Balay   ai  = a->i;
4789867e207SSatish Balay   aj  = a->j;
4799867e207SSatish Balay   aa  = a->a;
4809df24120SSatish Balay   bs2 = a->bs2;
4819867e207SSatish Balay 
4829867e207SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range");
4839867e207SSatish Balay 
4849867e207SSatish Balay   bn  = row/bs;   /* Block number */
4859867e207SSatish Balay   bp  = row % bs; /* Block Position */
4869867e207SSatish Balay   M   = ai[bn+1] - ai[bn];
4879867e207SSatish Balay   *nz = bs*M;
4889867e207SSatish Balay 
4899867e207SSatish Balay   if (v) {
4909867e207SSatish Balay     *v = 0;
4919867e207SSatish Balay     if (*nz) {
4929867e207SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
4939867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
4949867e207SSatish Balay         v_i  = *v + i*bs;
4957e67e3f9SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
4967e67e3f9SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
4979867e207SSatish Balay       }
4989867e207SSatish Balay     }
4999867e207SSatish Balay   }
5009867e207SSatish Balay 
5019867e207SSatish Balay   if (idx) {
5029867e207SSatish Balay     *idx = 0;
5039867e207SSatish Balay     if (*nz) {
5049867e207SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
5059867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
5069867e207SSatish Balay         idx_i = *idx + i*bs;
5079867e207SSatish Balay         itmp  = bs*aj[ai[bn] + i];
5089867e207SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
5099867e207SSatish Balay       }
5109867e207SSatish Balay     }
5119867e207SSatish Balay   }
5129867e207SSatish Balay   return 0;
5139867e207SSatish Balay }
5149867e207SSatish Balay 
5159867e207SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
5169867e207SSatish Balay {
5179867e207SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
5189867e207SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
5199867e207SSatish Balay   return 0;
5209867e207SSatish Balay }
521b6490206SBarry Smith 
522f2501298SSatish Balay static int MatTranspose_SeqBAIJ(Mat A,Mat *B)
523f2501298SSatish Balay {
524f2501298SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
525f2501298SSatish Balay   Mat         C;
526f2501298SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
5279df24120SSatish Balay   int         *rows,*cols,bs2=a->bs2;
528f2501298SSatish Balay   Scalar      *array=a->a;
529f2501298SSatish Balay 
530f2501298SSatish Balay   if (B==PETSC_NULL && mbs!=nbs)
531f2501298SSatish Balay     SETERRQ(1,"MatTranspose_SeqBAIJ:Square matrix only for in-place");
532f2501298SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
533f2501298SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
534a7c10996SSatish Balay 
535a7c10996SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
536f2501298SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
537f2501298SSatish Balay   PetscFree(col);
538f2501298SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
539f2501298SSatish Balay   cols = rows + bs;
540f2501298SSatish Balay   for ( i=0; i<mbs; i++ ) {
541f2501298SSatish Balay     cols[0] = i*bs;
542f2501298SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
543f2501298SSatish Balay     len = ai[i+1] - ai[i];
544f2501298SSatish Balay     for ( j=0; j<len; j++ ) {
545f2501298SSatish Balay       rows[0] = (*aj++)*bs;
546f2501298SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
547f2501298SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
548f2501298SSatish Balay       array += bs2;
549f2501298SSatish Balay     }
550f2501298SSatish Balay   }
5511073c447SSatish Balay   PetscFree(rows);
552f2501298SSatish Balay 
5536d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
5546d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
555f2501298SSatish Balay 
556f2501298SSatish Balay   if (B != PETSC_NULL) {
557f2501298SSatish Balay     *B = C;
558f2501298SSatish Balay   } else {
559f2501298SSatish Balay     /* This isn't really an in-place transpose */
560f2501298SSatish Balay     PetscFree(a->a);
561f2501298SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
562f2501298SSatish Balay     if (a->diag) PetscFree(a->diag);
563f2501298SSatish Balay     if (a->ilen) PetscFree(a->ilen);
564f2501298SSatish Balay     if (a->imax) PetscFree(a->imax);
565f2501298SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
566f2501298SSatish Balay     PetscFree(a);
567f2501298SSatish Balay     PetscMemcpy(A,C,sizeof(struct _Mat));
568f2501298SSatish Balay     PetscHeaderDestroy(C);
569f2501298SSatish Balay   }
570f2501298SSatish Balay   return 0;
571f2501298SSatish Balay }
572f2501298SSatish Balay 
573f2501298SSatish Balay 
574584200bdSSatish Balay static int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
575584200bdSSatish Balay {
576584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
577584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
578a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
5799df24120SSatish Balay   int        mbs = a->mbs, bs2 = a->bs2;
580584200bdSSatish Balay   Scalar     *aa = a->a, *ap;
581584200bdSSatish Balay 
5826d4a8577SBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) return 0;
583584200bdSSatish Balay 
584584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
585584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
586584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
587584200bdSSatish Balay     if (fshift) {
588a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
589584200bdSSatish Balay       N = ailen[i];
590584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
591584200bdSSatish Balay         ip[j-fshift] = ip[j];
5927e67e3f9SSatish Balay         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
593584200bdSSatish Balay       }
594584200bdSSatish Balay     }
595584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
596584200bdSSatish Balay   }
597584200bdSSatish Balay   if (mbs) {
598584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
599584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
600584200bdSSatish Balay   }
601584200bdSSatish Balay   /* reset ilen and imax for each row */
602584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
603584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
604584200bdSSatish Balay   }
605a7c10996SSatish Balay   a->nz = ai[mbs];
606584200bdSSatish Balay 
607584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
608584200bdSSatish Balay   if (fshift && a->diag) {
609584200bdSSatish Balay     PetscFree(a->diag);
610584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
611584200bdSSatish Balay     a->diag = 0;
612584200bdSSatish Balay   }
613537820f0SBarry Smith   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Unneed storage space %d used %d, rows %d, block size %d\n",
614537820f0SBarry Smith            fshift*bs2,a->nz*bs2,m,a->bs);
615584200bdSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues %d\n",
616584200bdSSatish Balay            a->reallocs);
6174e220ebcSLois Curfman McInnes   A->info.nz_unneeded  = (double)fshift*bs2;
6184e220ebcSLois Curfman McInnes 
619584200bdSSatish Balay   return 0;
620584200bdSSatish Balay }
621584200bdSSatish Balay 
622b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A)
6232593348eSBarry Smith {
624b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6259df24120SSatish Balay   PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(Scalar));
6262593348eSBarry Smith   return 0;
6272593348eSBarry Smith }
6282593348eSBarry Smith 
629b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj)
6302593348eSBarry Smith {
6312593348eSBarry Smith   Mat         A  = (Mat) obj;
632b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6332593348eSBarry Smith 
6342593348eSBarry Smith #if defined(PETSC_LOG)
6352593348eSBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
6362593348eSBarry Smith #endif
6372593348eSBarry Smith   PetscFree(a->a);
6382593348eSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
6392593348eSBarry Smith   if (a->diag) PetscFree(a->diag);
6402593348eSBarry Smith   if (a->ilen) PetscFree(a->ilen);
6412593348eSBarry Smith   if (a->imax) PetscFree(a->imax);
6422593348eSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
643de6a44a3SBarry Smith   if (a->mult_work) PetscFree(a->mult_work);
6442593348eSBarry Smith   PetscFree(a);
645f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
646f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
6472593348eSBarry Smith   return 0;
6482593348eSBarry Smith }
6492593348eSBarry Smith 
650b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
6512593348eSBarry Smith {
652b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6536d4a8577SBarry Smith   if      (op == MAT_ROW_ORIENTED)              a->roworiented = 1;
6546d4a8577SBarry Smith   else if (op == MAT_COLUMN_ORIENTED)           a->roworiented = 0;
6556d4a8577SBarry Smith   else if (op == MAT_COLUMNS_SORTED)            a->sorted      = 1;
6566d4a8577SBarry Smith   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
6576d4a8577SBarry Smith   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
6586d4a8577SBarry Smith   else if (op == MAT_ROWS_SORTED ||
6596d4a8577SBarry Smith            op == MAT_SYMMETRIC ||
6606d4a8577SBarry Smith            op == MAT_STRUCTURALLY_SYMMETRIC ||
6616d4a8577SBarry Smith            op == MAT_YES_NEW_DIAGONALS)
66294a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
6636d4a8577SBarry Smith   else if (op == MAT_NO_NEW_DIAGONALS)
6646d4a8577SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:MAT_NO_NEW_DIAGONALS");}
6652593348eSBarry Smith   else
666b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
6672593348eSBarry Smith   return 0;
6682593348eSBarry Smith }
6692593348eSBarry Smith 
6702593348eSBarry Smith 
6712593348eSBarry Smith /* -------------------------------------------------------*/
6722593348eSBarry Smith /* Should check that shapes of vectors and matrices match */
6732593348eSBarry Smith /* -------------------------------------------------------*/
674b6490206SBarry Smith #include "pinclude/plapack.h"
675b6490206SBarry Smith 
67639b95ed1SBarry Smith static int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
6772593348eSBarry Smith {
678b6490206SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
67939b95ed1SBarry Smith   register Scalar *x,*z,*v,sum;
680c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,n;
6812593348eSBarry Smith 
682c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
683c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
684b6490206SBarry Smith 
685119a934aSSatish Balay   idx   = a->j;
686119a934aSSatish Balay   v     = a->a;
687119a934aSSatish Balay   ii    = a->i;
688119a934aSSatish Balay 
689119a934aSSatish Balay   for ( i=0; i<mbs; i++ ) {
690119a934aSSatish Balay     n    = ii[1] - ii[0]; ii++;
691119a934aSSatish Balay     sum  = 0.0;
692119a934aSSatish Balay     while (n--) sum += *v++ * x[*idx++];
6931a6a6d98SBarry Smith     z[i] = sum;
694119a934aSSatish Balay   }
695c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
696c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
69739b95ed1SBarry Smith   PLogFlops(2*a->nz - a->m);
69839b95ed1SBarry Smith   return 0;
69939b95ed1SBarry Smith }
70039b95ed1SBarry Smith 
70139b95ed1SBarry Smith static int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
70239b95ed1SBarry Smith {
70339b95ed1SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
70439b95ed1SBarry Smith   register Scalar *x,*z,*v,*xb,sum1,sum2;
70539b95ed1SBarry Smith   register Scalar x1,x2;
70639b95ed1SBarry Smith   int             mbs=a->mbs,i,*idx,*ii;
707c16cb8f2SBarry Smith   int             j,n;
70839b95ed1SBarry Smith 
709c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
710c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
71139b95ed1SBarry Smith 
71239b95ed1SBarry Smith   idx   = a->j;
71339b95ed1SBarry Smith   v     = a->a;
71439b95ed1SBarry Smith   ii    = a->i;
71539b95ed1SBarry Smith 
716119a934aSSatish Balay   for ( i=0; i<mbs; i++ ) {
717119a934aSSatish Balay     n  = ii[1] - ii[0]; ii++;
718119a934aSSatish Balay     sum1 = 0.0; sum2 = 0.0;
719119a934aSSatish Balay     for ( j=0; j<n; j++ ) {
720119a934aSSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
721119a934aSSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
722119a934aSSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
723119a934aSSatish Balay       v += 4;
724119a934aSSatish Balay     }
7251a6a6d98SBarry Smith     z[0] = sum1; z[1] = sum2;
726119a934aSSatish Balay     z += 2;
727119a934aSSatish Balay   }
728c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
729c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
73039b95ed1SBarry Smith   PLogFlops(4*a->nz - a->m);
73139b95ed1SBarry Smith   return 0;
73239b95ed1SBarry Smith }
73339b95ed1SBarry Smith 
73439b95ed1SBarry Smith static int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
73539b95ed1SBarry Smith {
73639b95ed1SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
73739b95ed1SBarry Smith   register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,x1,x2,x3;
738c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,j,n;
73939b95ed1SBarry Smith 
740c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
741c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
74239b95ed1SBarry Smith 
74339b95ed1SBarry Smith   idx   = a->j;
74439b95ed1SBarry Smith   v     = a->a;
74539b95ed1SBarry Smith   ii    = a->i;
74639b95ed1SBarry Smith 
747119a934aSSatish Balay   for ( i=0; i<mbs; i++ ) {
748119a934aSSatish Balay     n  = ii[1] - ii[0]; ii++;
749119a934aSSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
750119a934aSSatish Balay     for ( j=0; j<n; j++ ) {
751119a934aSSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
752119a934aSSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
753119a934aSSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
754119a934aSSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
755119a934aSSatish Balay       v += 9;
756119a934aSSatish Balay     }
7571a6a6d98SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3;
758119a934aSSatish Balay     z += 3;
759119a934aSSatish Balay   }
760c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
761c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
76239b95ed1SBarry Smith   PLogFlops(18*a->nz - a->m);
76339b95ed1SBarry Smith   return 0;
76439b95ed1SBarry Smith }
76539b95ed1SBarry Smith 
76639b95ed1SBarry Smith static int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
76739b95ed1SBarry Smith {
76839b95ed1SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
76939b95ed1SBarry Smith   register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,sum4;
77039b95ed1SBarry Smith   register Scalar x1,x2,x3,x4;
77139b95ed1SBarry Smith   int             mbs=a->mbs,i,*idx,*ii;
772c16cb8f2SBarry Smith   int             j,n;
77339b95ed1SBarry Smith 
774c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
775c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
77639b95ed1SBarry Smith 
77739b95ed1SBarry Smith   idx   = a->j;
77839b95ed1SBarry Smith   v     = a->a;
77939b95ed1SBarry Smith   ii    = a->i;
78039b95ed1SBarry Smith 
781119a934aSSatish Balay   for ( i=0; i<mbs; i++ ) {
782119a934aSSatish Balay     n  = ii[1] - ii[0]; ii++;
783119a934aSSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
784119a934aSSatish Balay     for ( j=0; j<n; j++ ) {
785119a934aSSatish Balay       xb = x + 4*(*idx++);
786119a934aSSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
787119a934aSSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
788119a934aSSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
789119a934aSSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
790119a934aSSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
791119a934aSSatish Balay       v += 16;
792119a934aSSatish Balay     }
7931a6a6d98SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
794119a934aSSatish Balay     z += 4;
795119a934aSSatish Balay   }
796c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
797c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
79839b95ed1SBarry Smith   PLogFlops(32*a->nz - a->m);
79939b95ed1SBarry Smith   return 0;
80039b95ed1SBarry Smith }
80139b95ed1SBarry Smith 
80239b95ed1SBarry Smith static int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
80339b95ed1SBarry Smith {
80439b95ed1SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
80539b95ed1SBarry Smith   register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,sum4,sum5;
80639b95ed1SBarry Smith   register Scalar x1,x2,x3,x4,x5;
807c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,j,n;
80839b95ed1SBarry Smith 
809c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
810c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
81139b95ed1SBarry Smith 
81239b95ed1SBarry Smith   idx   = a->j;
81339b95ed1SBarry Smith   v     = a->a;
81439b95ed1SBarry Smith   ii    = a->i;
81539b95ed1SBarry Smith 
816119a934aSSatish Balay   for ( i=0; i<mbs; i++ ) {
817119a934aSSatish Balay     n  = ii[1] - ii[0]; ii++;
818119a934aSSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
819119a934aSSatish Balay     for ( j=0; j<n; j++ ) {
820119a934aSSatish Balay       xb = x + 5*(*idx++);
821119a934aSSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
822119a934aSSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
823119a934aSSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
824119a934aSSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
825119a934aSSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
826119a934aSSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
827119a934aSSatish Balay       v += 25;
828119a934aSSatish Balay     }
8291a6a6d98SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
830119a934aSSatish Balay     z += 5;
831119a934aSSatish Balay   }
832c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
833c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
83439b95ed1SBarry Smith   PLogFlops(50*a->nz - a->m);
83539b95ed1SBarry Smith   return 0;
83639b95ed1SBarry Smith }
83739b95ed1SBarry Smith 
83839b95ed1SBarry Smith static int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
83939b95ed1SBarry Smith {
84039b95ed1SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
84139b95ed1SBarry Smith   register Scalar *x,*z,*v,*xb;
842c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
843c16cb8f2SBarry Smith   int             _One = 1,ncols,k;
844c16cb8f2SBarry Smith   Scalar          _DOne = 1.0,*work,*workt,_DZero = 0.0;
84539b95ed1SBarry Smith 
846c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
847c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
84839b95ed1SBarry Smith 
84939b95ed1SBarry Smith   idx   = a->j;
85039b95ed1SBarry Smith   v     = a->a;
85139b95ed1SBarry Smith   ii    = a->i;
85239b95ed1SBarry Smith 
85339b95ed1SBarry Smith 
854119a934aSSatish Balay   if (!a->mult_work) {
8553b547af2SSatish Balay     k = PetscMax(a->m,a->n);
8562b3076dcSSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
857119a934aSSatish Balay   }
858119a934aSSatish Balay   work = a->mult_work;
859119a934aSSatish Balay   for ( i=0; i<mbs; i++ ) {
860119a934aSSatish Balay     n     = ii[1] - ii[0]; ii++;
861119a934aSSatish Balay     ncols = n*bs;
862119a934aSSatish Balay     workt = work;
863119a934aSSatish Balay     for ( j=0; j<n; j++ ) {
864119a934aSSatish Balay       xb = x + bs*(*idx++);
865119a934aSSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
866119a934aSSatish Balay       workt += bs;
867119a934aSSatish Balay     }
8681a6a6d98SBarry Smith     LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One);
869119a934aSSatish Balay     v += n*bs2;
870119a934aSSatish Balay     z += bs;
871119a934aSSatish Balay   }
872c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
873c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
8741a6a6d98SBarry Smith   PLogFlops(2*a->nz*bs2 - a->m);
875bb42667fSSatish Balay   return 0;
876bb42667fSSatish Balay }
877bb42667fSSatish Balay 
878f44d3a6dSSatish Balay static int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
879f44d3a6dSSatish Balay {
880f44d3a6dSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
881f44d3a6dSSatish Balay   register Scalar *x,*y,*z,*v,sum;
882c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,n;
883f44d3a6dSSatish Balay 
884c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
885c16cb8f2SBarry Smith   VecGetArray_Fast(yy,y);
886c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
887f44d3a6dSSatish Balay 
888f44d3a6dSSatish Balay   idx   = a->j;
889f44d3a6dSSatish Balay   v     = a->a;
890f44d3a6dSSatish Balay   ii    = a->i;
891f44d3a6dSSatish Balay 
892f44d3a6dSSatish Balay   for ( i=0; i<mbs; i++ ) {
893f44d3a6dSSatish Balay     n    = ii[1] - ii[0]; ii++;
894f44d3a6dSSatish Balay     sum  = y[i];
895f44d3a6dSSatish Balay     while (n--) sum += *v++ * x[*idx++];
896f44d3a6dSSatish Balay     z[i] = sum;
897f44d3a6dSSatish Balay   }
898c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
899c16cb8f2SBarry Smith   VecRestoreArray_Fast(yy,y);
900c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
901f44d3a6dSSatish Balay   PLogFlops(2*a->nz);
902f44d3a6dSSatish Balay   return 0;
903f44d3a6dSSatish Balay }
904f44d3a6dSSatish Balay 
905f44d3a6dSSatish Balay static int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
906f44d3a6dSSatish Balay {
907f44d3a6dSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
908f44d3a6dSSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2;
909f44d3a6dSSatish Balay   register Scalar x1,x2;
910f44d3a6dSSatish Balay   int             mbs=a->mbs,i,*idx,*ii;
911c16cb8f2SBarry Smith   int             j,n;
912f44d3a6dSSatish Balay 
913c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
914c16cb8f2SBarry Smith   VecGetArray_Fast(yy,y);
915c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
916f44d3a6dSSatish Balay 
917f44d3a6dSSatish Balay   idx   = a->j;
918f44d3a6dSSatish Balay   v     = a->a;
919f44d3a6dSSatish Balay   ii    = a->i;
920f44d3a6dSSatish Balay 
921f44d3a6dSSatish Balay   for ( i=0; i<mbs; i++ ) {
922f44d3a6dSSatish Balay     n  = ii[1] - ii[0]; ii++;
923f44d3a6dSSatish Balay     sum1 = y[0]; sum2 = y[1];
924f44d3a6dSSatish Balay     for ( j=0; j<n; j++ ) {
925f44d3a6dSSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
926f44d3a6dSSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
927f44d3a6dSSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
928f44d3a6dSSatish Balay       v += 4;
929f44d3a6dSSatish Balay     }
930f44d3a6dSSatish Balay     z[0] = sum1; z[1] = sum2;
931f44d3a6dSSatish Balay     z += 2; y += 2;
932f44d3a6dSSatish Balay   }
933c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
934c16cb8f2SBarry Smith   VecRestoreArray_Fast(yy,y);
935c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
936f44d3a6dSSatish Balay   PLogFlops(4*a->nz);
937f44d3a6dSSatish Balay   return 0;
938f44d3a6dSSatish Balay }
939f44d3a6dSSatish Balay 
940f44d3a6dSSatish Balay static int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
941f44d3a6dSSatish Balay {
942f44d3a6dSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
943f44d3a6dSSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,x1,x2,x3;
944c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,j,n;
945f44d3a6dSSatish Balay 
946c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
947c16cb8f2SBarry Smith   VecGetArray_Fast(yy,y);
948c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
949f44d3a6dSSatish Balay 
950f44d3a6dSSatish Balay   idx   = a->j;
951f44d3a6dSSatish Balay   v     = a->a;
952f44d3a6dSSatish Balay   ii    = a->i;
953f44d3a6dSSatish Balay 
954f44d3a6dSSatish Balay   for ( i=0; i<mbs; i++ ) {
955f44d3a6dSSatish Balay     n  = ii[1] - ii[0]; ii++;
956f44d3a6dSSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
957f44d3a6dSSatish Balay     for ( j=0; j<n; j++ ) {
958f44d3a6dSSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
959f44d3a6dSSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
960f44d3a6dSSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
961f44d3a6dSSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
962f44d3a6dSSatish Balay       v += 9;
963f44d3a6dSSatish Balay     }
964f44d3a6dSSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
965f44d3a6dSSatish Balay     z += 3; y += 3;
966f44d3a6dSSatish Balay   }
967c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
968c16cb8f2SBarry Smith   VecRestoreArray_Fast(yy,y);
969c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
970f44d3a6dSSatish Balay   PLogFlops(18*a->nz);
971f44d3a6dSSatish Balay   return 0;
972f44d3a6dSSatish Balay }
973f44d3a6dSSatish Balay 
974f44d3a6dSSatish Balay static int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
975f44d3a6dSSatish Balay {
976f44d3a6dSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
977f44d3a6dSSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,sum4;
978f44d3a6dSSatish Balay   register Scalar x1,x2,x3,x4;
979f44d3a6dSSatish Balay   int             mbs=a->mbs,i,*idx,*ii;
980c16cb8f2SBarry Smith   int             j,n;
981f44d3a6dSSatish Balay 
982c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
983c16cb8f2SBarry Smith   VecGetArray_Fast(yy,y);
984c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
985f44d3a6dSSatish Balay 
986f44d3a6dSSatish Balay   idx   = a->j;
987f44d3a6dSSatish Balay   v     = a->a;
988f44d3a6dSSatish Balay   ii    = a->i;
989f44d3a6dSSatish Balay 
990f44d3a6dSSatish Balay   for ( i=0; i<mbs; i++ ) {
991f44d3a6dSSatish Balay     n  = ii[1] - ii[0]; ii++;
992f44d3a6dSSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
993f44d3a6dSSatish Balay     for ( j=0; j<n; j++ ) {
994f44d3a6dSSatish Balay       xb = x + 4*(*idx++);
995f44d3a6dSSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
996f44d3a6dSSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
997f44d3a6dSSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
998f44d3a6dSSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
999f44d3a6dSSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
1000f44d3a6dSSatish Balay       v += 16;
1001f44d3a6dSSatish Balay     }
1002f44d3a6dSSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
1003f44d3a6dSSatish Balay     z += 4; y += 4;
1004f44d3a6dSSatish Balay   }
1005c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
1006c16cb8f2SBarry Smith   VecRestoreArray_Fast(yy,y);
1007c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
1008f44d3a6dSSatish Balay   PLogFlops(32*a->nz);
1009f44d3a6dSSatish Balay   return 0;
1010f44d3a6dSSatish Balay }
1011f44d3a6dSSatish Balay 
1012f44d3a6dSSatish Balay static int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
1013f44d3a6dSSatish Balay {
1014f44d3a6dSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
1015f44d3a6dSSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,sum4,sum5;
1016f44d3a6dSSatish Balay   register Scalar x1,x2,x3,x4,x5;
1017c16cb8f2SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,j,n;
1018f44d3a6dSSatish Balay 
1019c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
1020c16cb8f2SBarry Smith   VecGetArray_Fast(yy,y);
1021c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
1022f44d3a6dSSatish Balay 
1023f44d3a6dSSatish Balay   idx   = a->j;
1024f44d3a6dSSatish Balay   v     = a->a;
1025f44d3a6dSSatish Balay   ii    = a->i;
1026f44d3a6dSSatish Balay 
1027f44d3a6dSSatish Balay   for ( i=0; i<mbs; i++ ) {
1028f44d3a6dSSatish Balay     n  = ii[1] - ii[0]; ii++;
1029f44d3a6dSSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
1030f44d3a6dSSatish Balay     for ( j=0; j<n; j++ ) {
1031f44d3a6dSSatish Balay       xb = x + 5*(*idx++);
1032f44d3a6dSSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
1033f44d3a6dSSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
1034f44d3a6dSSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
1035f44d3a6dSSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
1036f44d3a6dSSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
1037f44d3a6dSSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
1038f44d3a6dSSatish Balay       v += 25;
1039f44d3a6dSSatish Balay     }
1040f44d3a6dSSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
1041f44d3a6dSSatish Balay     z += 5; y += 5;
1042f44d3a6dSSatish Balay   }
1043c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
1044c16cb8f2SBarry Smith   VecRestoreArray_Fast(yy,y);
1045c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
1046f44d3a6dSSatish Balay   PLogFlops(50*a->nz);
1047f44d3a6dSSatish Balay   return 0;
1048f44d3a6dSSatish Balay }
1049f44d3a6dSSatish Balay 
1050f44d3a6dSSatish Balay static int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
1051f44d3a6dSSatish Balay {
1052f44d3a6dSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
1053f44d3a6dSSatish Balay   register Scalar *x,*z,*v,*xb;
1054f44d3a6dSSatish Balay   int             mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr;
1055f44d3a6dSSatish Balay   int             _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
1056f44d3a6dSSatish Balay 
1057f44d3a6dSSatish Balay   if ( xx != yy) { ierr = VecCopy(yy,zz); CHKERRQ(ierr); }
1058f44d3a6dSSatish Balay 
1059c16cb8f2SBarry Smith   VecGetArray_Fast(xx,x);
1060c16cb8f2SBarry Smith   VecGetArray_Fast(zz,z);
1061f44d3a6dSSatish Balay 
1062f44d3a6dSSatish Balay   idx   = a->j;
1063f44d3a6dSSatish Balay   v     = a->a;
1064f44d3a6dSSatish Balay   ii    = a->i;
1065f44d3a6dSSatish Balay 
1066f44d3a6dSSatish Balay 
1067f44d3a6dSSatish Balay   if (!a->mult_work) {
1068f44d3a6dSSatish Balay     k = PetscMax(a->m,a->n);
1069f44d3a6dSSatish Balay     a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));CHKPTRQ(a->mult_work);
1070f44d3a6dSSatish Balay   }
1071f44d3a6dSSatish Balay   work = a->mult_work;
1072f44d3a6dSSatish Balay   for ( i=0; i<mbs; i++ ) {
1073f44d3a6dSSatish Balay     n     = ii[1] - ii[0]; ii++;
1074f44d3a6dSSatish Balay     ncols = n*bs;
1075f44d3a6dSSatish Balay     workt = work;
1076f44d3a6dSSatish Balay     for ( j=0; j<n; j++ ) {
1077f44d3a6dSSatish Balay       xb = x + bs*(*idx++);
1078f44d3a6dSSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
1079f44d3a6dSSatish Balay       workt += bs;
1080f44d3a6dSSatish Balay     }
1081f44d3a6dSSatish Balay     LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One);
1082f44d3a6dSSatish Balay     v += n*bs2;
1083f44d3a6dSSatish Balay     z += bs;
1084f44d3a6dSSatish Balay   }
1085c16cb8f2SBarry Smith   VecRestoreArray_Fast(xx,x);
1086c16cb8f2SBarry Smith   VecRestoreArray_Fast(zz,z);
1087f44d3a6dSSatish Balay   PLogFlops(2*a->nz*bs2 );
1088f44d3a6dSSatish Balay   return 0;
1089f44d3a6dSSatish Balay }
1090f44d3a6dSSatish Balay 
10911a6a6d98SBarry Smith static int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz)
1092bb42667fSSatish Balay {
1093bb42667fSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
10941a6a6d98SBarry Smith   Scalar          *xg,*zg,*zb;
1095bb42667fSSatish Balay   register Scalar *x,*z,*v,*xb,x1,x2,x3,x4,x5;
1096bb1453f0SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
10979df24120SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
1098119a934aSSatish Balay 
1099119a934aSSatish Balay 
1100bb42667fSSatish Balay   ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg;
1101bb42667fSSatish Balay   ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg;
1102bb1453f0SSatish Balay   PetscMemzero(z,N*sizeof(Scalar));
1103bb42667fSSatish Balay 
1104119a934aSSatish Balay   idx   = a->j;
1105119a934aSSatish Balay   v     = a->a;
1106119a934aSSatish Balay   ii    = a->i;
1107119a934aSSatish Balay 
1108119a934aSSatish Balay   switch (bs) {
1109119a934aSSatish Balay   case 1:
1110119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
1111119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
1112119a934aSSatish Balay       xb = x + i; x1 = xb[0];
1113119a934aSSatish Balay       ib = idx + ai[i];
1114119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
1115bb1453f0SSatish Balay         rval    = ib[j];
1116bb1453f0SSatish Balay         z[rval] += *v++ * x1;
1117119a934aSSatish Balay       }
1118119a934aSSatish Balay     }
1119119a934aSSatish Balay     break;
1120119a934aSSatish Balay   case 2:
1121119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
1122119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
1123119a934aSSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
1124119a934aSSatish Balay       ib = idx + ai[i];
1125119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
1126119a934aSSatish Balay         rval      = ib[j]*2;
1127bb1453f0SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
1128bb1453f0SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
1129119a934aSSatish Balay         v += 4;
1130119a934aSSatish Balay       }
1131119a934aSSatish Balay     }
1132119a934aSSatish Balay     break;
1133119a934aSSatish Balay   case 3:
1134119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
1135119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
1136119a934aSSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1137119a934aSSatish Balay       ib = idx + ai[i];
1138119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
1139119a934aSSatish Balay         rval      = ib[j]*3;
1140bb1453f0SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
1141bb1453f0SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
1142bb1453f0SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
1143119a934aSSatish Balay         v += 9;
1144119a934aSSatish Balay       }
1145119a934aSSatish Balay     }
1146119a934aSSatish Balay     break;
1147119a934aSSatish Balay   case 4:
1148119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
1149119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
1150119a934aSSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1151119a934aSSatish Balay       ib = idx + ai[i];
1152119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
1153119a934aSSatish Balay         rval      = ib[j]*4;
1154bb1453f0SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
1155bb1453f0SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
1156bb1453f0SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
1157bb1453f0SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
1158119a934aSSatish Balay         v += 16;
1159119a934aSSatish Balay       }
1160119a934aSSatish Balay     }
1161119a934aSSatish Balay     break;
1162119a934aSSatish Balay   case 5:
1163119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
1164119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
1165119a934aSSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1166119a934aSSatish Balay       x4 = xb[3];   x5 = xb[4];
1167119a934aSSatish Balay       ib = idx + ai[i];
1168119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
1169119a934aSSatish Balay         rval      = ib[j]*5;
1170bb1453f0SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
1171bb1453f0SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
1172bb1453f0SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
1173bb1453f0SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
1174bb1453f0SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
1175119a934aSSatish Balay         v += 25;
1176119a934aSSatish Balay       }
1177119a934aSSatish Balay     }
1178119a934aSSatish Balay     break;
1179119a934aSSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
1180119a934aSSatish Balay     default: {
1181119a934aSSatish Balay       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
1182119a934aSSatish Balay       if (!a->mult_work) {
11833b547af2SSatish Balay         k = PetscMax(a->m,a->n);
1184bb42667fSSatish Balay         a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));
1185119a934aSSatish Balay         CHKPTRQ(a->mult_work);
1186119a934aSSatish Balay       }
1187119a934aSSatish Balay       work = a->mult_work;
1188119a934aSSatish Balay       for ( i=0; i<mbs; i++ ) {
1189119a934aSSatish Balay         n     = ii[1] - ii[0]; ii++;
1190119a934aSSatish Balay         ncols = n*bs;
1191119a934aSSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
1192119a934aSSatish Balay         LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One);
1193119a934aSSatish Balay         v += n*bs2;
1194119a934aSSatish Balay         x += bs;
1195119a934aSSatish Balay         workt = work;
1196119a934aSSatish Balay         for ( j=0; j<n; j++ ) {
1197119a934aSSatish Balay           zb = z + bs*(*idx++);
1198bb1453f0SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
1199119a934aSSatish Balay           workt += bs;
1200119a934aSSatish Balay         }
1201119a934aSSatish Balay       }
1202119a934aSSatish Balay     }
1203119a934aSSatish Balay   }
12040419e6cdSSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
12050419e6cdSSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
1206faf6580fSSatish Balay   PLogFlops(2*a->nz*a->bs2 - a->n);
1207faf6580fSSatish Balay   return 0;
1208faf6580fSSatish Balay }
1209faf6580fSSatish Balay static int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1210faf6580fSSatish Balay {
1211faf6580fSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
1212faf6580fSSatish Balay   Scalar          *xg,*zg,*zb;
1213faf6580fSSatish Balay   register Scalar *x,*z,*v,*xb,x1,x2,x3,x4,x5;
1214faf6580fSSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
1215faf6580fSSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
1216faf6580fSSatish Balay 
1217faf6580fSSatish Balay 
1218faf6580fSSatish Balay 
1219faf6580fSSatish Balay   ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg;
1220faf6580fSSatish Balay   ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg;
1221faf6580fSSatish Balay 
1222648c1d95SSatish Balay   if ( yy != zz ) { ierr = VecCopy(yy,zz); CHKERRQ(ierr); }
1223648c1d95SSatish Balay   else PetscMemzero(z,N*sizeof(Scalar));
1224648c1d95SSatish Balay 
1225faf6580fSSatish Balay 
1226faf6580fSSatish Balay   idx   = a->j;
1227faf6580fSSatish Balay   v     = a->a;
1228faf6580fSSatish Balay   ii    = a->i;
1229faf6580fSSatish Balay 
1230faf6580fSSatish Balay   switch (bs) {
1231faf6580fSSatish Balay   case 1:
1232faf6580fSSatish Balay     for ( i=0; i<mbs; i++ ) {
1233faf6580fSSatish Balay       n  = ii[1] - ii[0]; ii++;
1234faf6580fSSatish Balay       xb = x + i; x1 = xb[0];
1235faf6580fSSatish Balay       ib = idx + ai[i];
1236faf6580fSSatish Balay       for ( j=0; j<n; j++ ) {
1237faf6580fSSatish Balay         rval    = ib[j];
1238faf6580fSSatish Balay         z[rval] += *v++ * x1;
1239faf6580fSSatish Balay       }
1240faf6580fSSatish Balay     }
1241faf6580fSSatish Balay     break;
1242faf6580fSSatish Balay   case 2:
1243faf6580fSSatish Balay     for ( i=0; i<mbs; i++ ) {
1244faf6580fSSatish Balay       n  = ii[1] - ii[0]; ii++;
1245faf6580fSSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
1246faf6580fSSatish Balay       ib = idx + ai[i];
1247faf6580fSSatish Balay       for ( j=0; j<n; j++ ) {
1248faf6580fSSatish Balay         rval      = ib[j]*2;
1249faf6580fSSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
1250faf6580fSSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
1251faf6580fSSatish Balay         v += 4;
1252faf6580fSSatish Balay       }
1253faf6580fSSatish Balay     }
1254faf6580fSSatish Balay     break;
1255faf6580fSSatish Balay   case 3:
1256faf6580fSSatish Balay     for ( i=0; i<mbs; i++ ) {
1257faf6580fSSatish Balay       n  = ii[1] - ii[0]; ii++;
1258faf6580fSSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1259faf6580fSSatish Balay       ib = idx + ai[i];
1260faf6580fSSatish Balay       for ( j=0; j<n; j++ ) {
1261faf6580fSSatish Balay         rval      = ib[j]*3;
1262faf6580fSSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
1263faf6580fSSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
1264faf6580fSSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
1265faf6580fSSatish Balay         v += 9;
1266faf6580fSSatish Balay       }
1267faf6580fSSatish Balay     }
1268faf6580fSSatish Balay     break;
1269faf6580fSSatish Balay   case 4:
1270faf6580fSSatish Balay     for ( i=0; i<mbs; i++ ) {
1271faf6580fSSatish Balay       n  = ii[1] - ii[0]; ii++;
1272faf6580fSSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1273faf6580fSSatish Balay       ib = idx + ai[i];
1274faf6580fSSatish Balay       for ( j=0; j<n; j++ ) {
1275faf6580fSSatish Balay         rval      = ib[j]*4;
1276faf6580fSSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
1277faf6580fSSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
1278faf6580fSSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
1279faf6580fSSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
1280faf6580fSSatish Balay         v += 16;
1281faf6580fSSatish Balay       }
1282faf6580fSSatish Balay     }
1283faf6580fSSatish Balay     break;
1284faf6580fSSatish Balay   case 5:
1285faf6580fSSatish Balay     for ( i=0; i<mbs; i++ ) {
1286faf6580fSSatish Balay       n  = ii[1] - ii[0]; ii++;
1287faf6580fSSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1288faf6580fSSatish Balay       x4 = xb[3];   x5 = xb[4];
1289faf6580fSSatish Balay       ib = idx + ai[i];
1290faf6580fSSatish Balay       for ( j=0; j<n; j++ ) {
1291faf6580fSSatish Balay         rval      = ib[j]*5;
1292faf6580fSSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
1293faf6580fSSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
1294faf6580fSSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
1295faf6580fSSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
1296faf6580fSSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
1297faf6580fSSatish Balay         v += 25;
1298faf6580fSSatish Balay       }
1299faf6580fSSatish Balay     }
1300faf6580fSSatish Balay     break;
1301faf6580fSSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
1302faf6580fSSatish Balay     default: {
1303faf6580fSSatish Balay       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
1304faf6580fSSatish Balay       if (!a->mult_work) {
1305faf6580fSSatish Balay         k = PetscMax(a->m,a->n);
1306faf6580fSSatish Balay         a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));
1307faf6580fSSatish Balay         CHKPTRQ(a->mult_work);
1308faf6580fSSatish Balay       }
1309faf6580fSSatish Balay       work = a->mult_work;
1310faf6580fSSatish Balay       for ( i=0; i<mbs; i++ ) {
1311faf6580fSSatish Balay         n     = ii[1] - ii[0]; ii++;
1312faf6580fSSatish Balay         ncols = n*bs;
1313faf6580fSSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
1314faf6580fSSatish Balay         LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One);
1315faf6580fSSatish Balay         v += n*bs2;
1316faf6580fSSatish Balay         x += bs;
1317faf6580fSSatish Balay         workt = work;
1318faf6580fSSatish Balay         for ( j=0; j<n; j++ ) {
1319faf6580fSSatish Balay           zb = z + bs*(*idx++);
1320faf6580fSSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
1321faf6580fSSatish Balay           workt += bs;
1322faf6580fSSatish Balay         }
1323faf6580fSSatish Balay       }
1324faf6580fSSatish Balay     }
1325faf6580fSSatish Balay   }
1326faf6580fSSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
1327faf6580fSSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
1328faf6580fSSatish Balay   PLogFlops(2*a->nz*a->bs2);
13292593348eSBarry Smith   return 0;
13302593348eSBarry Smith }
13312593348eSBarry Smith 
13324e220ebcSLois Curfman McInnes static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
13332593348eSBarry Smith {
1334b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
13354e220ebcSLois Curfman McInnes 
13364e220ebcSLois Curfman McInnes   info->rows_global    = (double)a->m;
13374e220ebcSLois Curfman McInnes   info->columns_global = (double)a->n;
13384e220ebcSLois Curfman McInnes   info->rows_local     = (double)a->m;
13394e220ebcSLois Curfman McInnes   info->columns_local  = (double)a->n;
13404e220ebcSLois Curfman McInnes   info->block_size     = a->bs2;
13414e220ebcSLois Curfman McInnes   info->nz_allocated   = a->maxnz;
13424e220ebcSLois Curfman McInnes   info->nz_used        = a->bs2*a->nz;
13434e220ebcSLois Curfman McInnes   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
13444e220ebcSLois Curfman McInnes   /*  if (info->nz_unneeded != A->info.nz_unneeded)
13454e220ebcSLois Curfman McInnes     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
13464e220ebcSLois Curfman McInnes   info->assemblies   = A->num_ass;
13474e220ebcSLois Curfman McInnes   info->mallocs      = a->reallocs;
13484e220ebcSLois Curfman McInnes   info->memory       = A->mem;
13494e220ebcSLois Curfman McInnes   if (A->factor) {
13504e220ebcSLois Curfman McInnes     info->fill_ratio_given  = A->info.fill_ratio_given;
13514e220ebcSLois Curfman McInnes     info->fill_ratio_needed = A->info.fill_ratio_needed;
13524e220ebcSLois Curfman McInnes     info->factor_mallocs    = A->info.factor_mallocs;
13534e220ebcSLois Curfman McInnes   } else {
13544e220ebcSLois Curfman McInnes     info->fill_ratio_given  = 0;
13554e220ebcSLois Curfman McInnes     info->fill_ratio_needed = 0;
13564e220ebcSLois Curfman McInnes     info->factor_mallocs    = 0;
13574e220ebcSLois Curfman McInnes   }
13582593348eSBarry Smith   return 0;
13592593348eSBarry Smith }
13602593348eSBarry Smith 
136191d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
136291d316f6SSatish Balay {
136391d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
136491d316f6SSatish Balay 
136591d316f6SSatish Balay   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
136691d316f6SSatish Balay 
136791d316f6SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
136891d316f6SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
1369a7c10996SSatish Balay       (a->nz != b->nz)) {
137091d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
137191d316f6SSatish Balay   }
137291d316f6SSatish Balay 
137391d316f6SSatish Balay   /* if the a->i are the same */
137491d316f6SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
137591d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
137691d316f6SSatish Balay   }
137791d316f6SSatish Balay 
137891d316f6SSatish Balay   /* if a->j are the same */
137991d316f6SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
138091d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
138191d316f6SSatish Balay   }
138291d316f6SSatish Balay 
138391d316f6SSatish Balay   /* if a->a are the same */
138491d316f6SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
138591d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
138691d316f6SSatish Balay   }
138791d316f6SSatish Balay   *flg = PETSC_TRUE;
138891d316f6SSatish Balay   return 0;
138991d316f6SSatish Balay 
139091d316f6SSatish Balay }
139191d316f6SSatish Balay 
139291d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
139391d316f6SSatish Balay {
139491d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
13957e67e3f9SSatish Balay   int         i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
139617e48fc4SSatish Balay   Scalar      *x, zero = 0.0,*aa,*aa_j;
139717e48fc4SSatish Balay 
139817e48fc4SSatish Balay   bs   = a->bs;
139917e48fc4SSatish Balay   aa   = a->a;
140017e48fc4SSatish Balay   ai   = a->i;
140117e48fc4SSatish Balay   aj   = a->j;
140217e48fc4SSatish Balay   ambs = a->mbs;
14039df24120SSatish Balay   bs2  = a->bs2;
140491d316f6SSatish Balay 
140591d316f6SSatish Balay   VecSet(&zero,v);
140691d316f6SSatish Balay   VecGetArray(v,&x); VecGetLocalSize(v,&n);
14079867e207SSatish Balay   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector");
140817e48fc4SSatish Balay   for ( i=0; i<ambs; i++ ) {
140917e48fc4SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
141017e48fc4SSatish Balay       if (aj[j] == i) {
141117e48fc4SSatish Balay         row  = i*bs;
14127e67e3f9SSatish Balay         aa_j = aa+j*bs2;
14137e67e3f9SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
141491d316f6SSatish Balay         break;
141591d316f6SSatish Balay       }
141691d316f6SSatish Balay     }
141791d316f6SSatish Balay   }
141891d316f6SSatish Balay   return 0;
141991d316f6SSatish Balay }
142091d316f6SSatish Balay 
14219867e207SSatish Balay static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
14229867e207SSatish Balay {
14239867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
14249867e207SSatish Balay   Scalar      *l,*r,x,*v,*aa,*li,*ri;
14257e67e3f9SSatish Balay   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
14269867e207SSatish Balay 
14279867e207SSatish Balay   ai  = a->i;
14289867e207SSatish Balay   aj  = a->j;
14299867e207SSatish Balay   aa  = a->a;
14309867e207SSatish Balay   m   = a->m;
14319867e207SSatish Balay   n   = a->n;
14329867e207SSatish Balay   bs  = a->bs;
14339867e207SSatish Balay   mbs = a->mbs;
14349df24120SSatish Balay   bs2 = a->bs2;
14359867e207SSatish Balay   if (ll) {
14369867e207SSatish Balay     VecGetArray(ll,&l); VecGetSize(ll,&lm);
14379867e207SSatish Balay     if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length");
14389867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
14399867e207SSatish Balay       M  = ai[i+1] - ai[i];
14409867e207SSatish Balay       li = l + i*bs;
14417e67e3f9SSatish Balay       v  = aa + bs2*ai[i];
14429867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
14437e67e3f9SSatish Balay         for ( k=0; k<bs2; k++ ) {
14449867e207SSatish Balay           (*v++) *= li[k%bs];
14459867e207SSatish Balay         }
14469867e207SSatish Balay       }
14479867e207SSatish Balay     }
14489867e207SSatish Balay   }
14499867e207SSatish Balay 
14509867e207SSatish Balay   if (rr) {
14519867e207SSatish Balay     VecGetArray(rr,&r); VecGetSize(rr,&rn);
14529867e207SSatish Balay     if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length");
14539867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
14549867e207SSatish Balay       M  = ai[i+1] - ai[i];
14557e67e3f9SSatish Balay       v  = aa + bs2*ai[i];
14569867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
14579867e207SSatish Balay         ri = r + bs*aj[ai[i]+j];
14589867e207SSatish Balay         for ( k=0; k<bs; k++ ) {
14599867e207SSatish Balay           x = ri[k];
14609867e207SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
14619867e207SSatish Balay         }
14629867e207SSatish Balay       }
14639867e207SSatish Balay     }
14649867e207SSatish Balay   }
14659867e207SSatish Balay   return 0;
14669867e207SSatish Balay }
14679867e207SSatish Balay 
14689867e207SSatish Balay 
1469b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
1470b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
14712a38ea04SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
1472736121d4SSatish Balay extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,MatGetSubMatrixCall,Mat*);
1473736121d4SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**);
14741a6a6d98SBarry Smith 
14751a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
14761a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
14771a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
14781a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
14791a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
14801a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
14811a6a6d98SBarry Smith 
14827fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
14837fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
14847fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
14857fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
14867fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
14877fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
14882593348eSBarry Smith 
1489b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
14902593348eSBarry Smith {
1491b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
14922593348eSBarry Smith   Scalar      *v = a->a;
14932593348eSBarry Smith   double      sum = 0.0;
14949df24120SSatish Balay   int         i,nz=a->nz,bs2=a->bs2;
14952593348eSBarry Smith 
14962593348eSBarry Smith   if (type == NORM_FROBENIUS) {
14970419e6cdSSatish Balay     for (i=0; i< bs2*nz; i++ ) {
14982593348eSBarry Smith #if defined(PETSC_COMPLEX)
14992593348eSBarry Smith       sum += real(conj(*v)*(*v)); v++;
15002593348eSBarry Smith #else
15012593348eSBarry Smith       sum += (*v)*(*v); v++;
15022593348eSBarry Smith #endif
15032593348eSBarry Smith     }
15042593348eSBarry Smith     *norm = sqrt(sum);
15052593348eSBarry Smith   }
15062593348eSBarry Smith   else {
1507b6490206SBarry Smith     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
15082593348eSBarry Smith   }
15092593348eSBarry Smith   return 0;
15102593348eSBarry Smith }
15112593348eSBarry Smith 
15122593348eSBarry Smith /*
15132593348eSBarry Smith      note: This can only work for identity for row and col. It would
15142593348eSBarry Smith    be good to check this and otherwise generate an error.
15152593348eSBarry Smith */
1516b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
15172593348eSBarry Smith {
1518b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
15192593348eSBarry Smith   Mat         outA;
1520de6a44a3SBarry Smith   int         ierr;
15212593348eSBarry Smith 
1522b6490206SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
15232593348eSBarry Smith 
15242593348eSBarry Smith   outA          = inA;
15252593348eSBarry Smith   inA->factor   = FACTOR_LU;
15262593348eSBarry Smith   a->row        = row;
15272593348eSBarry Smith   a->col        = col;
15282593348eSBarry Smith 
1529de6a44a3SBarry Smith   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
15302593348eSBarry Smith 
15312593348eSBarry Smith   if (!a->diag) {
1532b6490206SBarry Smith     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
15332593348eSBarry Smith   }
15347fc0212eSBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
15352593348eSBarry Smith   return 0;
15362593348eSBarry Smith }
15372593348eSBarry Smith 
1538b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
15392593348eSBarry Smith {
1540b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
15419df24120SSatish Balay   int         one = 1, totalnz = a->bs2*a->nz;
1542b6490206SBarry Smith   BLscal_( &totalnz, alpha, a->a, &one );
1543b6490206SBarry Smith   PLogFlops(totalnz);
15442593348eSBarry Smith   return 0;
15452593348eSBarry Smith }
15462593348eSBarry Smith 
15477e67e3f9SSatish Balay static int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
15487e67e3f9SSatish Balay {
15497e67e3f9SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
15507e67e3f9SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1551a7c10996SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
15529df24120SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
15537e67e3f9SSatish Balay   Scalar     *ap, *aa = a->a, zero = 0.0;
15547e67e3f9SSatish Balay 
15557e67e3f9SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
15567e67e3f9SSatish Balay     row  = im[k]; brow = row/bs;
15577e67e3f9SSatish Balay     if (row < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative row");
15587e67e3f9SSatish Balay     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqBAIJ:Row too large");
1559a7c10996SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
15607e67e3f9SSatish Balay     nrow = ailen[brow];
15617e67e3f9SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
15627e67e3f9SSatish Balay       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative column");
15637e67e3f9SSatish Balay       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqBAIJ:Column too large");
1564a7c10996SSatish Balay       col  = in[l] ;
15657e67e3f9SSatish Balay       bcol = col/bs;
15667e67e3f9SSatish Balay       cidx = col%bs;
15677e67e3f9SSatish Balay       ridx = row%bs;
15687e67e3f9SSatish Balay       high = nrow;
15697e67e3f9SSatish Balay       low  = 0; /* assume unsorted */
15707e67e3f9SSatish Balay       while (high-low > 5) {
15717e67e3f9SSatish Balay         t = (low+high)/2;
15727e67e3f9SSatish Balay         if (rp[t] > bcol) high = t;
15737e67e3f9SSatish Balay         else             low  = t;
15747e67e3f9SSatish Balay       }
15757e67e3f9SSatish Balay       for ( i=low; i<high; i++ ) {
15767e67e3f9SSatish Balay         if (rp[i] > bcol) break;
15777e67e3f9SSatish Balay         if (rp[i] == bcol) {
15787e67e3f9SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
15797e67e3f9SSatish Balay           goto finished;
15807e67e3f9SSatish Balay         }
15817e67e3f9SSatish Balay       }
15827e67e3f9SSatish Balay       *v++ = zero;
15837e67e3f9SSatish Balay       finished:;
15847e67e3f9SSatish Balay     }
15857e67e3f9SSatish Balay   }
15867e67e3f9SSatish Balay   return 0;
15877e67e3f9SSatish Balay }
15887e67e3f9SSatish Balay 
15895a838052SSatish Balay static int MatGetBlockSize_SeqBAIJ(Mat mat, int *bs)
15905a838052SSatish Balay {
15915a838052SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *) mat->data;
15925a838052SSatish Balay   *bs = baij->bs;
15935a838052SSatish Balay   return 0;
15945a838052SSatish Balay }
15955a838052SSatish Balay 
1596d9b7c43dSSatish Balay /* idx should be of length atlease bs */
1597d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ_Check_Block(int *idx, int bs, PetscTruth *flg)
1598d9b7c43dSSatish Balay {
1599d9b7c43dSSatish Balay   int i,row;
1600d9b7c43dSSatish Balay   row = idx[0];
1601d9b7c43dSSatish Balay   if (row%bs!=0) { *flg = PETSC_FALSE; return 0; }
1602d9b7c43dSSatish Balay 
1603d9b7c43dSSatish Balay   for ( i=1; i<bs; i++ ) {
1604d9b7c43dSSatish Balay     if (row+i != idx[i]) { *flg = PETSC_FALSE; return 0; }
1605d9b7c43dSSatish Balay   }
1606d9b7c43dSSatish Balay   *flg = PETSC_TRUE;
1607d9b7c43dSSatish Balay   return 0;
1608d9b7c43dSSatish Balay }
1609d9b7c43dSSatish Balay 
1610d9b7c43dSSatish Balay static int MatZeroRows_SeqBAIJ(Mat A,IS is, Scalar *diag)
1611d9b7c43dSSatish Balay {
1612d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
1613d9b7c43dSSatish Balay   IS          is_local;
1614d9b7c43dSSatish Balay   int         ierr,i,j,count,m=baij->m,is_n,*is_idx,*rows,bs=baij->bs,bs2=baij->bs2;
1615d9b7c43dSSatish Balay   PetscTruth  flg;
1616d9b7c43dSSatish Balay   Scalar      zero = 0.0,*aa;
1617d9b7c43dSSatish Balay 
1618d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
1619d9b7c43dSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
1620d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
1621537820f0SBarry Smith   ierr = ISCreateGeneral(A->comm,is_n,is_idx,&is_local); CHKERRQ(ierr);
1622d9b7c43dSSatish Balay   ierr = ISSort(is_local); CHKERRQ(ierr);
1623d9b7c43dSSatish Balay   ierr = ISGetIndices(is_local,&rows); CHKERRQ(ierr);
1624d9b7c43dSSatish Balay 
1625d9b7c43dSSatish Balay   i = 0;
1626d9b7c43dSSatish Balay   while (i < is_n) {
1627d9b7c43dSSatish Balay     if (rows[i]<0 || rows[i]>m) SETERRQ(1,"MatZeroRows_SeqBAIJ:row out of range");
1628d9b7c43dSSatish Balay     flg = PETSC_FALSE;
1629d9b7c43dSSatish Balay     if (i+bs <= is_n) {ierr = MatZeroRows_SeqBAIJ_Check_Block(rows+i,bs,&flg); CHKERRQ(ierr); }
1630d9b7c43dSSatish Balay     if (flg) { /* There exists a block of rows to be Zerowed */
1631d9b7c43dSSatish Balay       baij->ilen[rows[i]/bs] = 0;
1632d9b7c43dSSatish Balay       i += bs;
1633d9b7c43dSSatish Balay     } else { /* Zero out only the requested row */
1634d9b7c43dSSatish Balay       count = (baij->i[rows[i]/bs +1] - baij->i[rows[i]/bs])*bs;
1635d9b7c43dSSatish Balay       aa    = baij->a + baij->i[rows[i]/bs]*bs2 + (rows[i]%bs);
1636d9b7c43dSSatish Balay       for ( j=0; j<count; j++ ) {
1637d9b7c43dSSatish Balay         aa[0] = zero;
1638d9b7c43dSSatish Balay         aa+=bs;
1639d9b7c43dSSatish Balay       }
1640d9b7c43dSSatish Balay       i++;
1641d9b7c43dSSatish Balay     }
1642d9b7c43dSSatish Balay   }
1643d9b7c43dSSatish Balay   if (diag) {
1644d9b7c43dSSatish Balay     for ( j=0; j<is_n; j++ ) {
1645d9b7c43dSSatish Balay       ierr = (*A->ops.setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1646d9b7c43dSSatish Balay     }
1647d9b7c43dSSatish Balay   }
1648d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is,&is_idx); CHKERRQ(ierr);
1649d9b7c43dSSatish Balay   ierr = ISRestoreIndices(is_local,&rows); CHKERRQ(ierr);
1650d9b7c43dSSatish Balay   ierr = ISDestroy(is_local); CHKERRQ(ierr);
1651d9b7c43dSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1652d9b7c43dSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1653d9b7c43dSSatish Balay 
1654d9b7c43dSSatish Balay   return 0;
1655d9b7c43dSSatish Balay }
1656ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1657ba4ca20aSSatish Balay {
1658ba4ca20aSSatish Balay   static int called = 0;
1659ba4ca20aSSatish Balay   MPI_Comm   comm = A->comm;
1660ba4ca20aSSatish Balay 
1661ba4ca20aSSatish Balay   if (called) return 0; else called = 1;
1662ba4ca20aSSatish Balay   PetscPrintf(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");
1663ba4ca20aSSatish Balay   PetscPrintf(comm,"  -mat_block_size <block_size>\n");
1664ba4ca20aSSatish Balay   return 0;
1665ba4ca20aSSatish Balay }
1666d9b7c43dSSatish Balay 
16672593348eSBarry Smith /* -------------------------------------------------------------------*/
1668cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
16699867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
1670f44d3a6dSSatish Balay        MatMult_SeqBAIJ_N,MatMultAdd_SeqBAIJ_N,
1671faf6580fSSatish Balay        MatMultTrans_SeqBAIJ,MatMultTransAdd_SeqBAIJ,
16721a6a6d98SBarry Smith        MatSolve_SeqBAIJ_N,0,
1673b6490206SBarry Smith        0,0,
1674de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
1675b6490206SBarry Smith        0,
1676f2501298SSatish Balay        MatTranspose_SeqBAIJ,
167717e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
16789867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
1679584200bdSSatish Balay        0,MatAssemblyEnd_SeqBAIJ,
1680b6490206SBarry Smith        0,
1681d9b7c43dSSatish Balay        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,MatZeroRows_SeqBAIJ,
16827fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
1683b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
1684de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
1685b6490206SBarry Smith        0,0,/*  MatConvert_SeqBAIJ  */ 0,
1686b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
1687b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
1688af015674SSatish Balay        MatGetSubMatrices_SeqBAIJ,MatIncreaseOverlap_SeqBAIJ,
16897e67e3f9SSatish Balay        MatGetValues_SeqBAIJ,0,
1690ba4ca20aSSatish Balay        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
16913b2fbd54SBarry Smith        0,0,0,MatGetBlockSize_SeqBAIJ,
16923b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
16933b2fbd54SBarry Smith        MatRestoreRowIJ_SeqBAIJ};
16942593348eSBarry Smith 
16952593348eSBarry Smith /*@C
169644cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
169744cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
169844cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
16992bd5e0b2SLois Curfman McInnes    (or the array nzz).  By setting these parameters accurately, performance
17002bd5e0b2SLois Curfman McInnes    during matrix assembly can be increased by more than a factor of 50.
17012593348eSBarry Smith 
17022593348eSBarry Smith    Input Parameters:
17032593348eSBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
1704b6490206SBarry Smith .  bs - size of block
17052593348eSBarry Smith .  m - number of rows
17062593348eSBarry Smith .  n - number of columns
1707b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
17082bd5e0b2SLois Curfman McInnes .  nzz - array containing the number of block nonzeros in the various block rows
17092bd5e0b2SLois Curfman McInnes          (possibly different for each block row) or PETSC_NULL
17102593348eSBarry Smith 
17112593348eSBarry Smith    Output Parameter:
17122593348eSBarry Smith .  A - the matrix
17132593348eSBarry Smith 
17142593348eSBarry Smith    Notes:
171544cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
17162593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
171744cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
17182593348eSBarry Smith 
17192593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
17202593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
17212593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
17222593348eSBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
17232593348eSBarry Smith 
172444cd7ae7SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
17252593348eSBarry Smith @*/
1726b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
17272593348eSBarry Smith {
17282593348eSBarry Smith   Mat         B;
1729b6490206SBarry Smith   Mat_SeqBAIJ *b;
17303b2fbd54SBarry Smith   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs,size;
17313b2fbd54SBarry Smith 
17323b2fbd54SBarry Smith   MPI_Comm_size(comm,&size);
17333b2fbd54SBarry Smith   if (size > 1) SETERRQ(1,"MatCreateSeqBAIJ:Comm must be of size 1");
1734b6490206SBarry Smith 
1735f2501298SSatish Balay   if (mbs*bs!=m || nbs*bs!=n)
1736f2501298SSatish Balay     SETERRQ(1,"MatCreateSeqBAIJ:Number rows, cols must be divisible by blocksize");
17372593348eSBarry Smith 
17382593348eSBarry Smith   *A = 0;
1739b6490206SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
17402593348eSBarry Smith   PLogObjectCreate(B);
1741b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
174244cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
17432593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
17441a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
17451a6a6d98SBarry Smith   if (!flg) {
17467fc0212eSBarry Smith     switch (bs) {
17477fc0212eSBarry Smith       case 1:
17487fc0212eSBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
17491a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_1;
175039b95ed1SBarry Smith         B->ops.mult            = MatMult_SeqBAIJ_1;
1751f44d3a6dSSatish Balay         B->ops.multadd         = MatMultAdd_SeqBAIJ_1;
17527fc0212eSBarry Smith        break;
17534eeb42bcSBarry Smith       case 2:
17544eeb42bcSBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
17551a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_2;
175639b95ed1SBarry Smith         B->ops.mult            = MatMult_SeqBAIJ_2;
1757f44d3a6dSSatish Balay         B->ops.multadd         = MatMultAdd_SeqBAIJ_2;
17586d84be18SBarry Smith         break;
1759f327dd97SBarry Smith       case 3:
1760f327dd97SBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
17611a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_3;
176239b95ed1SBarry Smith         B->ops.mult            = MatMult_SeqBAIJ_3;
1763f44d3a6dSSatish Balay         B->ops.multadd         = MatMultAdd_SeqBAIJ_3;
17644eeb42bcSBarry Smith         break;
17656d84be18SBarry Smith       case 4:
17666d84be18SBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
17671a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_4;
176839b95ed1SBarry Smith         B->ops.mult            = MatMult_SeqBAIJ_4;
1769f44d3a6dSSatish Balay         B->ops.multadd         = MatMultAdd_SeqBAIJ_4;
17706d84be18SBarry Smith         break;
17716d84be18SBarry Smith       case 5:
17726d84be18SBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
17731a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_5;
177439b95ed1SBarry Smith         B->ops.mult            = MatMult_SeqBAIJ_5;
1775f44d3a6dSSatish Balay         B->ops.multadd         = MatMultAdd_SeqBAIJ_5;
17766d84be18SBarry Smith         break;
17777fc0212eSBarry Smith     }
17781a6a6d98SBarry Smith   }
1779b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
1780b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
17812593348eSBarry Smith   B->factor           = 0;
17822593348eSBarry Smith   B->lupivotthreshold = 1.0;
17832593348eSBarry Smith   b->row              = 0;
17842593348eSBarry Smith   b->col              = 0;
17852593348eSBarry Smith   b->reallocs         = 0;
17862593348eSBarry Smith 
178744cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
178844cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1789b6490206SBarry Smith   b->mbs     = mbs;
1790f2501298SSatish Balay   b->nbs     = nbs;
1791b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
17922593348eSBarry Smith   if (nnz == PETSC_NULL) {
1793119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
17942593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1795b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1796b6490206SBarry Smith     nz = nz*mbs;
17972593348eSBarry Smith   }
17982593348eSBarry Smith   else {
17992593348eSBarry Smith     nz = 0;
1800b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
18012593348eSBarry Smith   }
18022593348eSBarry Smith 
18032593348eSBarry Smith   /* allocate the matrix space */
18047e67e3f9SSatish Balay   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
18052593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
18067e67e3f9SSatish Balay   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
18077e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
18082593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
18092593348eSBarry Smith   b->i  = b->j + nz;
18102593348eSBarry Smith   b->singlemalloc = 1;
18112593348eSBarry Smith 
1812de6a44a3SBarry Smith   b->i[0] = 0;
1813b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
18142593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
18152593348eSBarry Smith   }
18162593348eSBarry Smith 
1817b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1818b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1819b6490206SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
1820b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
18212593348eSBarry Smith 
1822b6490206SBarry Smith   b->bs               = bs;
18239df24120SSatish Balay   b->bs2              = bs2;
1824b6490206SBarry Smith   b->mbs              = mbs;
18252593348eSBarry Smith   b->nz               = 0;
182618e7b8e4SLois Curfman McInnes   b->maxnz            = nz*bs2;
18272593348eSBarry Smith   b->sorted           = 0;
18282593348eSBarry Smith   b->roworiented      = 1;
18292593348eSBarry Smith   b->nonew            = 0;
18302593348eSBarry Smith   b->diag             = 0;
18312593348eSBarry Smith   b->solve_work       = 0;
1832de6a44a3SBarry Smith   b->mult_work        = 0;
18332593348eSBarry Smith   b->spptr            = 0;
18344e220ebcSLois Curfman McInnes   B->info.nz_unneeded = (double)b->maxnz;
18354e220ebcSLois Curfman McInnes 
18362593348eSBarry Smith   *A = B;
18372593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
18382593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
18392593348eSBarry Smith   return 0;
18402593348eSBarry Smith }
18412593348eSBarry Smith 
1842b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
18432593348eSBarry Smith {
18442593348eSBarry Smith   Mat         C;
1845b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
18469df24120SSatish Balay   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
1847de6a44a3SBarry Smith 
1848de6a44a3SBarry Smith   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
18492593348eSBarry Smith 
18502593348eSBarry Smith   *B = 0;
1851b6490206SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
18522593348eSBarry Smith   PLogObjectCreate(C);
1853b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
18542593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1855b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
1856b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
18572593348eSBarry Smith   C->factor     = A->factor;
18582593348eSBarry Smith   c->row        = 0;
18592593348eSBarry Smith   c->col        = 0;
18602593348eSBarry Smith   C->assembled  = PETSC_TRUE;
18612593348eSBarry Smith 
186244cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
186344cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
186444cd7ae7SLois Curfman McInnes   C->M          = a->m;
186544cd7ae7SLois Curfman McInnes   C->N          = a->n;
186644cd7ae7SLois Curfman McInnes 
1867b6490206SBarry Smith   c->bs         = a->bs;
18689df24120SSatish Balay   c->bs2        = a->bs2;
1869b6490206SBarry Smith   c->mbs        = a->mbs;
1870b6490206SBarry Smith   c->nbs        = a->nbs;
18712593348eSBarry Smith 
1872b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1873b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1874b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
18752593348eSBarry Smith     c->imax[i] = a->imax[i];
18762593348eSBarry Smith     c->ilen[i] = a->ilen[i];
18772593348eSBarry Smith   }
18782593348eSBarry Smith 
18792593348eSBarry Smith   /* allocate the matrix space */
18802593348eSBarry Smith   c->singlemalloc = 1;
18817e67e3f9SSatish Balay   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
18822593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
18837e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1884de6a44a3SBarry Smith   c->i  = c->j + nz;
1885b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1886b6490206SBarry Smith   if (mbs > 0) {
1887de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
18882593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
18897e67e3f9SSatish Balay       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
18902593348eSBarry Smith     }
18912593348eSBarry Smith   }
18922593348eSBarry Smith 
1893b6490206SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
18942593348eSBarry Smith   c->sorted      = a->sorted;
18952593348eSBarry Smith   c->roworiented = a->roworiented;
18962593348eSBarry Smith   c->nonew       = a->nonew;
18972593348eSBarry Smith 
18982593348eSBarry Smith   if (a->diag) {
1899b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1900b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1901b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
19022593348eSBarry Smith       c->diag[i] = a->diag[i];
19032593348eSBarry Smith     }
19042593348eSBarry Smith   }
19052593348eSBarry Smith   else c->diag          = 0;
19062593348eSBarry Smith   c->nz                 = a->nz;
19072593348eSBarry Smith   c->maxnz              = a->maxnz;
19082593348eSBarry Smith   c->solve_work         = 0;
19092593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
19107fc0212eSBarry Smith   c->mult_work          = 0;
19112593348eSBarry Smith   *B = C;
19122593348eSBarry Smith   return 0;
19132593348eSBarry Smith }
19142593348eSBarry Smith 
191519bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
19162593348eSBarry Smith {
1917b6490206SBarry Smith   Mat_SeqBAIJ  *a;
19182593348eSBarry Smith   Mat          B;
1919de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1920b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
192135aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1922a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1923b6490206SBarry Smith   Scalar       *aa;
192419bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
19252593348eSBarry Smith 
19261a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1927de6a44a3SBarry Smith   bs2  = bs*bs;
1928b6490206SBarry Smith 
19292593348eSBarry Smith   MPI_Comm_size(comm,&size);
1930b6490206SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
193190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
193277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1933de6a44a3SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
19342593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
19352593348eSBarry Smith 
1936b6490206SBarry Smith   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
193735aab85fSBarry Smith 
193835aab85fSBarry Smith   /*
193935aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
194035aab85fSBarry Smith     divisible by the blocksize
194135aab85fSBarry Smith   */
1942b6490206SBarry Smith   mbs        = M/bs;
194335aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
194435aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
194535aab85fSBarry Smith   else                  mbs++;
194635aab85fSBarry Smith   if (extra_rows) {
1947537820f0SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
194835aab85fSBarry Smith   }
1949b6490206SBarry Smith 
19502593348eSBarry Smith   /* read in row lengths */
195135aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
195277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
195335aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
19542593348eSBarry Smith 
1955b6490206SBarry Smith   /* read in column indices */
195635aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
195777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
195835aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1959b6490206SBarry Smith 
1960b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1961b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1962b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
196335aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
196435aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
196535aab85fSBarry Smith   masked = mask + mbs;
1966b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1967b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
196835aab85fSBarry Smith     nmask = 0;
1969b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1970b6490206SBarry Smith       kmax = rowlengths[rowcount];
1971b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
197235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
197335aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1974b6490206SBarry Smith       }
1975b6490206SBarry Smith       rowcount++;
1976b6490206SBarry Smith     }
197735aab85fSBarry Smith     browlengths[i] += nmask;
197835aab85fSBarry Smith     /* zero out the mask elements we set */
197935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1980b6490206SBarry Smith   }
1981b6490206SBarry Smith 
19822593348eSBarry Smith   /* create our matrix */
198335aab85fSBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
198435aab85fSBarry Smith          CHKERRQ(ierr);
19852593348eSBarry Smith   B = *A;
1986b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
19872593348eSBarry Smith 
19882593348eSBarry Smith   /* set matrix "i" values */
1989de6a44a3SBarry Smith   a->i[0] = 0;
1990b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1991b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1992b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
19932593348eSBarry Smith   }
1994b6490206SBarry Smith   a->nz         = 0;
1995b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
19962593348eSBarry Smith 
1997b6490206SBarry Smith   /* read in nonzero values */
199835aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
199977c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
200035aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
2001b6490206SBarry Smith 
2002b6490206SBarry Smith   /* set "a" and "j" values into matrix */
2003b6490206SBarry Smith   nzcount = 0; jcount = 0;
2004b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
2005b6490206SBarry Smith     nzcountb = nzcount;
200635aab85fSBarry Smith     nmask    = 0;
2007b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
2008b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
2009b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
201035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
201135aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
2012b6490206SBarry Smith       }
2013b6490206SBarry Smith       rowcount++;
2014b6490206SBarry Smith     }
2015de6a44a3SBarry Smith     /* sort the masked values */
201677c4ece6SBarry Smith     PetscSortInt(nmask,masked);
2017de6a44a3SBarry Smith 
2018b6490206SBarry Smith     /* set "j" values into matrix */
2019b6490206SBarry Smith     maskcount = 1;
202035aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
202135aab85fSBarry Smith       a->j[jcount++]  = masked[j];
2022de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
2023b6490206SBarry Smith     }
2024b6490206SBarry Smith     /* set "a" values into matrix */
2025de6a44a3SBarry Smith     ishift = bs2*a->i[i];
2026b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
2027b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
2028b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
2029de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
2030de6a44a3SBarry Smith         block  = mask[tmp] - 1;
2031de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
2032de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
2033b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
2034b6490206SBarry Smith       }
2035b6490206SBarry Smith     }
203635aab85fSBarry Smith     /* zero out the mask elements we set */
203735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
2038b6490206SBarry Smith   }
203935aab85fSBarry Smith   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
2040b6490206SBarry Smith 
2041b6490206SBarry Smith   PetscFree(rowlengths);
2042b6490206SBarry Smith   PetscFree(browlengths);
2043b6490206SBarry Smith   PetscFree(aa);
2044b6490206SBarry Smith   PetscFree(jj);
2045b6490206SBarry Smith   PetscFree(mask);
2046b6490206SBarry Smith 
2047b6490206SBarry Smith   B->assembled = PETSC_TRUE;
2048de6a44a3SBarry Smith 
2049de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
2050de6a44a3SBarry Smith   if (flg) {
205119bcc07fSBarry Smith     Viewer tviewer;
205219bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
2053639f9d9dSBarry Smith     ierr = ViewerSetFormat(tviewer,VIEWER_FORMAT_ASCII_INFO,0);CHKERRQ(ierr);
205419bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
205519bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
2056de6a44a3SBarry Smith   }
2057de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
2058de6a44a3SBarry Smith   if (flg) {
205919bcc07fSBarry Smith     Viewer tviewer;
206019bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
2061639f9d9dSBarry Smith     ierr = ViewerSetFormat(tviewer,VIEWER_FORMAT_ASCII_INFO_LONG,0);CHKERRQ(ierr);
206219bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
206319bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
2064de6a44a3SBarry Smith   }
2065de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
2066de6a44a3SBarry Smith   if (flg) {
206719bcc07fSBarry Smith     Viewer tviewer;
206819bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
206919bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
207019bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
2071de6a44a3SBarry Smith   }
2072de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
2073de6a44a3SBarry Smith   if (flg) {
207419bcc07fSBarry Smith     Viewer tviewer;
207519bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
2076639f9d9dSBarry Smith     ierr = ViewerSetFormat(tviewer,VIEWER_FORMAT_ASCII_MATLAB,"M");CHKERRQ(ierr);
207719bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
207819bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
2079de6a44a3SBarry Smith   }
2080de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
2081de6a44a3SBarry Smith   if (flg) {
208219bcc07fSBarry Smith     Viewer tviewer;
208319bcc07fSBarry Smith     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
208419bcc07fSBarry Smith     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
208519bcc07fSBarry Smith     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
208619bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
2087de6a44a3SBarry Smith   }
20882593348eSBarry Smith   return 0;
20892593348eSBarry Smith }
20902593348eSBarry Smith 
20912593348eSBarry Smith 
20922593348eSBarry Smith 
2093