xref: /petsc/src/mat/impls/baij/seq/baij.c (revision bb1453f09d95c8e1a2680f8e579ecf9418beaf69)
1b6490206SBarry Smith 
22593348eSBarry Smith #ifndef lint
3*bb1453f0SSatish Balay static char vcid[] = "$Id: baij.c,v 1.43 1996/04/30 23:04:41 balay Exp balay $";
42593348eSBarry Smith #endif
52593348eSBarry Smith 
62593348eSBarry Smith /*
7b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
82593348eSBarry Smith   matrix storage format.
92593348eSBarry Smith */
10b6490206SBarry Smith #include "baij.h"
111a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
121a6a6d98SBarry Smith #include "src/inline/spops.h"
132593348eSBarry Smith #include "petsc.h"
143b547af2SSatish Balay 
15bcd2baecSBarry Smith extern   int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
162593348eSBarry Smith 
17b6490206SBarry Smith static int MatGetReordering_SeqBAIJ(Mat A,MatOrdering type,IS *rperm,IS *cperm)
182593348eSBarry Smith {
19b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
20bcd2baecSBarry Smith   int         ierr, *ia, *ja,n = a->mbs,*idx,i,ishift,oshift;
212593348eSBarry Smith 
222593348eSBarry Smith   /*
232593348eSBarry Smith      this is tacky: In the future when we have written special factorization
242593348eSBarry Smith      and solve routines for the identity permutation we should use a
252593348eSBarry Smith      stride index set instead of the general one.
262593348eSBarry Smith   */
272593348eSBarry Smith   if (type  == ORDER_NATURAL) {
282593348eSBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
292593348eSBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
302593348eSBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
312593348eSBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
322593348eSBarry Smith     PetscFree(idx);
332593348eSBarry Smith     ISSetPermutation(*rperm);
342593348eSBarry Smith     ISSetPermutation(*cperm);
352593348eSBarry Smith     ISSetIdentity(*rperm);
362593348eSBarry Smith     ISSetIdentity(*cperm);
372593348eSBarry Smith     return 0;
382593348eSBarry Smith   }
392593348eSBarry Smith 
40bcd2baecSBarry Smith   MatReorderingRegisterAll();
41a7c10996SSatish Balay   ishift = 0;
42bcd2baecSBarry Smith   oshift = -MatReorderingIndexShift[(int)type];
43bcd2baecSBarry Smith   if (MatReorderingRequiresSymmetric[(int)type]) {
441a6a6d98SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,ishift,oshift,&ia,&ja);CHKERRQ(ierr);
451a6a6d98SBarry Smith     ierr = MatGetReordering_IJ(n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
462593348eSBarry Smith     PetscFree(ia); PetscFree(ja);
47bcd2baecSBarry Smith   } else {
48bcd2baecSBarry Smith     if (ishift == oshift) {
491a6a6d98SBarry Smith       ierr = MatGetReordering_IJ(n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
50bcd2baecSBarry Smith     }
51bcd2baecSBarry Smith     else if (ishift == -1) {
52bcd2baecSBarry Smith       /* temporarily subtract 1 from i and j indices */
531a6a6d98SBarry Smith       int nz = a->i[n] - 1;
54bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
551a6a6d98SBarry Smith       for ( i=0; i<n+1; i++ ) a->i[i]--;
561a6a6d98SBarry Smith       ierr = MatGetReordering_IJ(n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
57bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
581a6a6d98SBarry Smith       for ( i=0; i<n+1; i++ ) a->i[i]++;
59bcd2baecSBarry Smith     } else {
60bcd2baecSBarry Smith       /* temporarily add 1 to i and j indices */
611a6a6d98SBarry Smith       int nz = a->i[n] - 1;
62bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
631a6a6d98SBarry Smith       for ( i=0; i<n+1; i++ ) a->i[i]++;
641a6a6d98SBarry Smith       ierr = MatGetReordering_IJ(n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
65bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
661a6a6d98SBarry Smith       for ( i=0; i<n+1; i++ ) a->i[i]--;
67bcd2baecSBarry Smith     }
68bcd2baecSBarry Smith   }
692593348eSBarry Smith   return 0;
702593348eSBarry Smith }
712593348eSBarry Smith 
72de6a44a3SBarry Smith /*
73de6a44a3SBarry Smith      Adds diagonal pointers to sparse matrix structure.
74de6a44a3SBarry Smith */
75de6a44a3SBarry Smith 
76de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A)
77de6a44a3SBarry Smith {
78de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
797fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
80de6a44a3SBarry Smith 
81de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
82de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
837fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
84de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
85de6a44a3SBarry Smith       if (a->j[j] == i) {
86de6a44a3SBarry Smith         diag[i] = j;
87de6a44a3SBarry Smith         break;
88de6a44a3SBarry Smith       }
89de6a44a3SBarry Smith     }
90de6a44a3SBarry Smith   }
91de6a44a3SBarry Smith   a->diag = diag;
92de6a44a3SBarry Smith   return 0;
93de6a44a3SBarry Smith }
942593348eSBarry Smith 
952593348eSBarry Smith #include "draw.h"
962593348eSBarry Smith #include "pinclude/pviewer.h"
9777c4ece6SBarry Smith #include "sys.h"
982593348eSBarry Smith 
99b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
1002593348eSBarry Smith {
101b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1029df24120SSatish Balay   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
103b6490206SBarry Smith   Scalar      *aa;
1042593348eSBarry Smith 
10590ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1062593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
1072593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
1082593348eSBarry Smith   col_lens[1] = a->m;
1092593348eSBarry Smith   col_lens[2] = a->n;
1107e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
1112593348eSBarry Smith 
1122593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
113b6490206SBarry Smith   count = 0;
114b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
115b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
116b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
117b6490206SBarry Smith     }
1182593348eSBarry Smith   }
11977c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
1202593348eSBarry Smith   PetscFree(col_lens);
1212593348eSBarry Smith 
1222593348eSBarry Smith   /* store column indices (zero start index) */
1237e67e3f9SSatish Balay   jj = (int *) PetscMalloc( a->nz*bs2*sizeof(int) ); CHKPTRQ(jj);
124b6490206SBarry Smith   count = 0;
125b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
126b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
127b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
128b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
129b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
1302593348eSBarry Smith         }
1312593348eSBarry Smith       }
132b6490206SBarry Smith     }
133b6490206SBarry Smith   }
1347e67e3f9SSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,BINARY_INT,0); CHKERRQ(ierr);
135b6490206SBarry Smith   PetscFree(jj);
1362593348eSBarry Smith 
1372593348eSBarry Smith   /* store nonzero values */
1387e67e3f9SSatish Balay   aa = (Scalar *) PetscMalloc(a->nz*bs2*sizeof(Scalar)); CHKPTRQ(aa);
139b6490206SBarry Smith   count = 0;
140b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
141b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
142b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
143b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
1447e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
145b6490206SBarry Smith         }
146b6490206SBarry Smith       }
147b6490206SBarry Smith     }
148b6490206SBarry Smith   }
1497e67e3f9SSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
150b6490206SBarry Smith   PetscFree(aa);
1512593348eSBarry Smith   return 0;
1522593348eSBarry Smith }
1532593348eSBarry Smith 
154b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
1552593348eSBarry Smith {
156b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1579df24120SSatish Balay   int         ierr, i,j,format,bs = a->bs,k,l,bs2=a->bs2;
1582593348eSBarry Smith   FILE        *fd;
1592593348eSBarry Smith   char        *outputname;
1602593348eSBarry Smith 
16190ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
1622593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
16390ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
1647ddc982cSLois Curfman McInnes   if (format == ASCII_FORMAT_INFO || format == ASCII_FORMAT_INFO_DETAILED) {
1657ddc982cSLois Curfman McInnes     fprintf(fd,"  block size is %d\n",bs);
1662593348eSBarry Smith   }
16790ace30eSBarry Smith   else if (format == ASCII_FORMAT_MATLAB) {
168b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported");
1692593348eSBarry Smith   }
17044cd7ae7SLois Curfman McInnes   else if (format == ASCII_FORMAT_COMMON) {
17144cd7ae7SLois Curfman McInnes     for ( i=0; i<a->mbs; i++ ) {
17244cd7ae7SLois Curfman McInnes       for ( j=0; j<bs; j++ ) {
17344cd7ae7SLois Curfman McInnes         fprintf(fd,"row %d:",i*bs+j);
17444cd7ae7SLois Curfman McInnes         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
17544cd7ae7SLois Curfman McInnes           for ( l=0; l<bs; l++ ) {
17644cd7ae7SLois Curfman McInnes #if defined(PETSC_COMPLEX)
17744cd7ae7SLois Curfman McInnes           if (imag(a->a[bs2*k + l*bs + j]) != 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
17844cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
17944cd7ae7SLois Curfman McInnes               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
18044cd7ae7SLois Curfman McInnes           else if (real(a->a[bs2*k + l*bs + j]) != 0.0)
18144cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
18244cd7ae7SLois Curfman McInnes #else
18344cd7ae7SLois Curfman McInnes           if (a->a[bs2*k + l*bs + j] != 0.0)
18444cd7ae7SLois Curfman McInnes             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
18544cd7ae7SLois Curfman McInnes #endif
18644cd7ae7SLois Curfman McInnes           }
18744cd7ae7SLois Curfman McInnes         }
18844cd7ae7SLois Curfman McInnes         fprintf(fd,"\n");
18944cd7ae7SLois Curfman McInnes       }
19044cd7ae7SLois Curfman McInnes     }
19144cd7ae7SLois Curfman McInnes   }
1922593348eSBarry Smith   else {
193b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
194b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
195b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
196b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
197b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
19888685aaeSLois Curfman McInnes #if defined(PETSC_COMPLEX)
1997e67e3f9SSatish Balay           if (imag(a->a[bs2*k + l*bs + j]) != 0.0) {
20088685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
2017e67e3f9SSatish Balay               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
20288685aaeSLois Curfman McInnes           }
20388685aaeSLois Curfman McInnes           else {
2047e67e3f9SSatish Balay             fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
20588685aaeSLois Curfman McInnes           }
20688685aaeSLois Curfman McInnes #else
2077e67e3f9SSatish Balay             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
20888685aaeSLois Curfman McInnes #endif
2092593348eSBarry Smith           }
2102593348eSBarry Smith         }
2112593348eSBarry Smith         fprintf(fd,"\n");
2122593348eSBarry Smith       }
2132593348eSBarry Smith     }
214b6490206SBarry Smith   }
2152593348eSBarry Smith   fflush(fd);
2162593348eSBarry Smith   return 0;
2172593348eSBarry Smith }
2182593348eSBarry Smith 
219b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
2202593348eSBarry Smith {
2212593348eSBarry Smith   Mat         A = (Mat) obj;
22219bcc07fSBarry Smith   ViewerType  vtype;
22319bcc07fSBarry Smith   int         ierr;
2242593348eSBarry Smith 
2252593348eSBarry Smith   if (!viewer) {
22619bcc07fSBarry Smith     viewer = STDOUT_VIEWER_SELF;
2272593348eSBarry Smith   }
22819bcc07fSBarry Smith 
22919bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
23019bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
231b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported");
2322593348eSBarry Smith   }
23319bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
234b6490206SBarry Smith     return MatView_SeqBAIJ_ASCII(A,viewer);
2352593348eSBarry Smith   }
23619bcc07fSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
237b6490206SBarry Smith     return MatView_SeqBAIJ_Binary(A,viewer);
2382593348eSBarry Smith   }
23919bcc07fSBarry Smith   else if (vtype == DRAW_VIEWER) {
240b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported");
2412593348eSBarry Smith   }
2422593348eSBarry Smith   return 0;
2432593348eSBarry Smith }
244b6490206SBarry Smith 
245119a934aSSatish Balay #define CHUNKSIZE  10
246cd0e1443SSatish Balay 
247cd0e1443SSatish Balay /* This version has row oriented v  */
248cd0e1443SSatish Balay static int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
249cd0e1443SSatish Balay {
250cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
251cd0e1443SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
252cd0e1443SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
253a7c10996SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
2549df24120SSatish Balay   int          ridx,cidx,bs2=a->bs2;
255cd0e1443SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
256cd0e1443SSatish Balay 
257cd0e1443SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
258cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
259cd0e1443SSatish Balay     if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row");
260cd0e1443SSatish Balay     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large");
261a7c10996SSatish Balay     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
262cd0e1443SSatish Balay     rmax = imax[brow]; nrow = ailen[brow];
263cd0e1443SSatish Balay     low = 0;
264cd0e1443SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
265cd0e1443SSatish Balay       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column");
266cd0e1443SSatish Balay       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large");
267a7c10996SSatish Balay       col = in[l]; bcol = col/bs;
268cd0e1443SSatish Balay       ridx = row % bs; cidx = col % bs;
269cd0e1443SSatish Balay       if (roworiented) {
270cd0e1443SSatish Balay         value = *v++;
271cd0e1443SSatish Balay       }
272cd0e1443SSatish Balay       else {
273cd0e1443SSatish Balay         value = v[k + l*m];
274cd0e1443SSatish Balay       }
275cd0e1443SSatish Balay       if (!sorted) low = 0; high = nrow;
276cd0e1443SSatish Balay       while (high-low > 5) {
277cd0e1443SSatish Balay         t = (low+high)/2;
278cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
279cd0e1443SSatish Balay         else             low  = t;
280cd0e1443SSatish Balay       }
281cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
282cd0e1443SSatish Balay         if (rp[i] > bcol) break;
283cd0e1443SSatish Balay         if (rp[i] == bcol) {
2847e67e3f9SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
285cd0e1443SSatish Balay           if (is == ADD_VALUES) *bap += value;
286cd0e1443SSatish Balay           else                  *bap  = value;
287cd0e1443SSatish Balay           goto noinsert;
288cd0e1443SSatish Balay         }
289cd0e1443SSatish Balay       }
290cd0e1443SSatish Balay       if (nonew) goto noinsert;
291cd0e1443SSatish Balay       if (nrow >= rmax) {
292cd0e1443SSatish Balay         /* there is no extra room in row, therefore enlarge */
293cd0e1443SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
294cd0e1443SSatish Balay         Scalar *new_a;
295cd0e1443SSatish Balay 
296cd0e1443SSatish Balay         /* malloc new storage space */
2977e67e3f9SSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
298cd0e1443SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
2997e67e3f9SSatish Balay         new_j   = (int *) (new_a + bs2*new_nz);
300cd0e1443SSatish Balay         new_i   = new_j + new_nz;
301cd0e1443SSatish Balay 
302cd0e1443SSatish Balay         /* copy over old data into new slots */
303cd0e1443SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
3047e67e3f9SSatish Balay         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
305a7c10996SSatish Balay         PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));
306a7c10996SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow);
307a7c10996SSatish Balay         PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,
308cd0e1443SSatish Balay                                                            len*sizeof(int));
309a7c10996SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(Scalar));
310a7c10996SSatish Balay         PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
311a7c10996SSatish Balay         PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),
312a7c10996SSatish Balay                     aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(Scalar));
313cd0e1443SSatish Balay         /* free up old matrix storage */
314cd0e1443SSatish Balay         PetscFree(a->a);
315cd0e1443SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
316cd0e1443SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
317cd0e1443SSatish Balay         a->singlemalloc = 1;
318cd0e1443SSatish Balay 
319a7c10996SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
320cd0e1443SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
3217e67e3f9SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
322119a934aSSatish Balay         a->maxnz += CHUNKSIZE;
323cd0e1443SSatish Balay         a->reallocs++;
324119a934aSSatish Balay         a->nz++;
325cd0e1443SSatish Balay       }
3267e67e3f9SSatish Balay       N = nrow++ - 1;
327cd0e1443SSatish Balay       /* shift up all the later entries in this row */
328cd0e1443SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
329cd0e1443SSatish Balay         rp[ii+1] = rp[ii];
3307e67e3f9SSatish Balay          PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
331cd0e1443SSatish Balay       }
3327e67e3f9SSatish Balay       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
333cd0e1443SSatish Balay       rp[i] = bcol;
3347e67e3f9SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
335cd0e1443SSatish Balay       noinsert:;
336cd0e1443SSatish Balay       low = i;
337cd0e1443SSatish Balay     }
338cd0e1443SSatish Balay     ailen[brow] = nrow;
339cd0e1443SSatish Balay   }
340cd0e1443SSatish Balay   return 0;
341cd0e1443SSatish Balay }
342cd0e1443SSatish Balay 
3430b824a48SSatish Balay static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
3440b824a48SSatish Balay {
3450b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3460b824a48SSatish Balay   *m = a->m; *n = a->n;
3470b824a48SSatish Balay   return 0;
3480b824a48SSatish Balay }
3490b824a48SSatish Balay 
3500b824a48SSatish Balay static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
3510b824a48SSatish Balay {
3520b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3530b824a48SSatish Balay   *m = 0; *n = a->m;
3540b824a48SSatish Balay   return 0;
3550b824a48SSatish Balay }
3560b824a48SSatish Balay 
3579867e207SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
3589867e207SSatish Balay {
3599867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3607e67e3f9SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
3619867e207SSatish Balay   Scalar      *aa,*v_i,*aa_i;
3629867e207SSatish Balay 
3639867e207SSatish Balay   bs  = a->bs;
3649867e207SSatish Balay   ai  = a->i;
3659867e207SSatish Balay   aj  = a->j;
3669867e207SSatish Balay   aa  = a->a;
3679df24120SSatish Balay   bs2 = a->bs2;
3689867e207SSatish Balay 
3699867e207SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range");
3709867e207SSatish Balay 
3719867e207SSatish Balay   bn  = row/bs;   /* Block number */
3729867e207SSatish Balay   bp  = row % bs; /* Block Position */
3739867e207SSatish Balay   M   = ai[bn+1] - ai[bn];
3749867e207SSatish Balay   *nz = bs*M;
3759867e207SSatish Balay 
3769867e207SSatish Balay   if (v) {
3779867e207SSatish Balay     *v = 0;
3789867e207SSatish Balay     if (*nz) {
3799867e207SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
3809867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
3819867e207SSatish Balay         v_i  = *v + i*bs;
3827e67e3f9SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
3837e67e3f9SSatish Balay         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
3849867e207SSatish Balay       }
3859867e207SSatish Balay     }
3869867e207SSatish Balay   }
3879867e207SSatish Balay 
3889867e207SSatish Balay   if (idx) {
3899867e207SSatish Balay     *idx = 0;
3909867e207SSatish Balay     if (*nz) {
3919867e207SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
3929867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
3939867e207SSatish Balay         idx_i = *idx + i*bs;
3949867e207SSatish Balay         itmp  = bs*aj[ai[bn] + i];
3959867e207SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
3969867e207SSatish Balay       }
3979867e207SSatish Balay     }
3989867e207SSatish Balay   }
3999867e207SSatish Balay   return 0;
4009867e207SSatish Balay }
4019867e207SSatish Balay 
4029867e207SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
4039867e207SSatish Balay {
4049867e207SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
4059867e207SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
4069867e207SSatish Balay   return 0;
4079867e207SSatish Balay }
408b6490206SBarry Smith 
409f2501298SSatish Balay static int MatTranspose_SeqBAIJ(Mat A,Mat *B)
410f2501298SSatish Balay {
411f2501298SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
412f2501298SSatish Balay   Mat         C;
413f2501298SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
4149df24120SSatish Balay   int         *rows,*cols,bs2=a->bs2;
415f2501298SSatish Balay   Scalar      *array=a->a;
416f2501298SSatish Balay 
417f2501298SSatish Balay   if (B==PETSC_NULL && mbs!=nbs)
418f2501298SSatish Balay     SETERRQ(1,"MatTranspose_SeqBAIJ:Square matrix only for in-place");
419f2501298SSatish Balay   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
420f2501298SSatish Balay   PetscMemzero(col,(1+nbs)*sizeof(int));
421a7c10996SSatish Balay 
422a7c10996SSatish Balay   for ( i=0; i<ai[mbs]; i++ ) col[aj[i]] += 1;
423f2501298SSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
424f2501298SSatish Balay   PetscFree(col);
425f2501298SSatish Balay   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
426f2501298SSatish Balay   cols = rows + bs;
427f2501298SSatish Balay   for ( i=0; i<mbs; i++ ) {
428f2501298SSatish Balay     cols[0] = i*bs;
429f2501298SSatish Balay     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
430f2501298SSatish Balay     len = ai[i+1] - ai[i];
431f2501298SSatish Balay     for ( j=0; j<len; j++ ) {
432f2501298SSatish Balay       rows[0] = (*aj++)*bs;
433f2501298SSatish Balay       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
434f2501298SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
435f2501298SSatish Balay       array += bs2;
436f2501298SSatish Balay     }
437f2501298SSatish Balay   }
4381073c447SSatish Balay   PetscFree(rows);
439f2501298SSatish Balay 
440f2501298SSatish Balay   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
441f2501298SSatish Balay   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
442f2501298SSatish Balay 
443f2501298SSatish Balay   if (B != PETSC_NULL) {
444f2501298SSatish Balay     *B = C;
445f2501298SSatish Balay   } else {
446f2501298SSatish Balay     /* This isn't really an in-place transpose */
447f2501298SSatish Balay     PetscFree(a->a);
448f2501298SSatish Balay     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
449f2501298SSatish Balay     if (a->diag) PetscFree(a->diag);
450f2501298SSatish Balay     if (a->ilen) PetscFree(a->ilen);
451f2501298SSatish Balay     if (a->imax) PetscFree(a->imax);
452f2501298SSatish Balay     if (a->solve_work) PetscFree(a->solve_work);
453f2501298SSatish Balay     PetscFree(a);
454f2501298SSatish Balay     PetscMemcpy(A,C,sizeof(struct _Mat));
455f2501298SSatish Balay     PetscHeaderDestroy(C);
456f2501298SSatish Balay   }
457f2501298SSatish Balay   return 0;
458f2501298SSatish Balay }
459f2501298SSatish Balay 
460f2501298SSatish Balay 
461584200bdSSatish Balay static int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
462584200bdSSatish Balay {
463584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
464584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
465a7c10996SSatish Balay   int        m = a->m,*ip, N, *ailen = a->ilen;
4669df24120SSatish Balay   int        mbs = a->mbs, bs2 = a->bs2;
467584200bdSSatish Balay   Scalar     *aa = a->a, *ap;
468584200bdSSatish Balay 
469584200bdSSatish Balay   if (mode == FLUSH_ASSEMBLY) return 0;
470584200bdSSatish Balay 
471584200bdSSatish Balay   for ( i=1; i<mbs; i++ ) {
472584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
473584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
474584200bdSSatish Balay     if (fshift) {
475a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
476584200bdSSatish Balay       N = ailen[i];
477584200bdSSatish Balay       for ( j=0; j<N; j++ ) {
478584200bdSSatish Balay         ip[j-fshift] = ip[j];
4797e67e3f9SSatish Balay         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
480584200bdSSatish Balay       }
481584200bdSSatish Balay     }
482584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
483584200bdSSatish Balay   }
484584200bdSSatish Balay   if (mbs) {
485584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
486584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
487584200bdSSatish Balay   }
488584200bdSSatish Balay   /* reset ilen and imax for each row */
489584200bdSSatish Balay   for ( i=0; i<mbs; i++ ) {
490584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
491584200bdSSatish Balay   }
492a7c10996SSatish Balay   a->nz = ai[mbs];
493584200bdSSatish Balay 
494584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
495584200bdSSatish Balay   if (fshift && a->diag) {
496584200bdSSatish Balay     PetscFree(a->diag);
497584200bdSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
498584200bdSSatish Balay     a->diag = 0;
499584200bdSSatish Balay   }
50067790700SSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ: Unneed storage space(blocks) %d used %d, rows %d, block size %d\n", fshift*bs2,a->nz*bs2,m,a->bs);
501584200bdSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ: Number of mallocs during MatSetValues %d\n",
502584200bdSSatish Balay            a->reallocs);
503584200bdSSatish Balay   return 0;
504584200bdSSatish Balay }
505584200bdSSatish Balay 
506b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A)
5072593348eSBarry Smith {
508b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5099df24120SSatish Balay   PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(Scalar));
5102593348eSBarry Smith   return 0;
5112593348eSBarry Smith }
5122593348eSBarry Smith 
513b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj)
5142593348eSBarry Smith {
5152593348eSBarry Smith   Mat         A  = (Mat) obj;
516b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5172593348eSBarry Smith 
5182593348eSBarry Smith #if defined(PETSC_LOG)
5192593348eSBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
5202593348eSBarry Smith #endif
5212593348eSBarry Smith   PetscFree(a->a);
5222593348eSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
5232593348eSBarry Smith   if (a->diag) PetscFree(a->diag);
5242593348eSBarry Smith   if (a->ilen) PetscFree(a->ilen);
5252593348eSBarry Smith   if (a->imax) PetscFree(a->imax);
5262593348eSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
527de6a44a3SBarry Smith   if (a->mult_work) PetscFree(a->mult_work);
5282593348eSBarry Smith   PetscFree(a);
529f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
530f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
5312593348eSBarry Smith   return 0;
5322593348eSBarry Smith }
5332593348eSBarry Smith 
534b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
5352593348eSBarry Smith {
536b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5372593348eSBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
5382593348eSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
5392593348eSBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
5402593348eSBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
5412593348eSBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
5422593348eSBarry Smith   else if (op == ROWS_SORTED ||
5432593348eSBarry Smith            op == SYMMETRIC_MATRIX ||
5442593348eSBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
5452593348eSBarry Smith            op == YES_NEW_DIAGONALS)
54694a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
5472593348eSBarry Smith   else if (op == NO_NEW_DIAGONALS)
548b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
5492593348eSBarry Smith   else
550b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
5512593348eSBarry Smith   return 0;
5522593348eSBarry Smith }
5532593348eSBarry Smith 
5542593348eSBarry Smith 
5552593348eSBarry Smith /* -------------------------------------------------------*/
5562593348eSBarry Smith /* Should check that shapes of vectors and matrices match */
5572593348eSBarry Smith /* -------------------------------------------------------*/
558b6490206SBarry Smith #include "pinclude/plapack.h"
559b6490206SBarry Smith 
5601a6a6d98SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec zz)
5612593348eSBarry Smith {
562b6490206SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5631a6a6d98SBarry Smith   Scalar          *xg,*zg;
564bb42667fSSatish Balay   register Scalar *x,*z,*v,sum,*xb,sum1,sum2,sum3,sum4,sum5;
565b6490206SBarry Smith   register Scalar x1,x2,x3,x4,x5;
5661a6a6d98SBarry Smith   int             mbs=a->mbs,i,*idx,*ii;
5679df24120SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,ierr;
5682593348eSBarry Smith 
569bb42667fSSatish Balay   ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg;
570bb42667fSSatish Balay   ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg;
571b6490206SBarry Smith 
572119a934aSSatish Balay   idx   = a->j;
573119a934aSSatish Balay   v     = a->a;
574119a934aSSatish Balay   ii    = a->i;
575119a934aSSatish Balay 
576119a934aSSatish Balay   switch (bs) {
577119a934aSSatish Balay   case 1:
578119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
579119a934aSSatish Balay       n    = ii[1] - ii[0]; ii++;
580119a934aSSatish Balay       sum  = 0.0;
581119a934aSSatish Balay       while (n--) sum += *v++ * x[*idx++];
5821a6a6d98SBarry Smith       z[i] = sum;
583119a934aSSatish Balay     }
584119a934aSSatish Balay     break;
585119a934aSSatish Balay   case 2:
586119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
587119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
588119a934aSSatish Balay       sum1 = 0.0; sum2 = 0.0;
589119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
590119a934aSSatish Balay         xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
591119a934aSSatish Balay         sum1 += v[0]*x1 + v[2]*x2;
592119a934aSSatish Balay         sum2 += v[1]*x1 + v[3]*x2;
593119a934aSSatish Balay         v += 4;
594119a934aSSatish Balay       }
5951a6a6d98SBarry Smith       z[0] = sum1; z[1] = sum2;
596119a934aSSatish Balay       z += 2;
597119a934aSSatish Balay     }
598119a934aSSatish Balay     break;
599119a934aSSatish Balay   case 3:
600119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
601119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
602119a934aSSatish Balay       sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
603119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
604119a934aSSatish Balay         xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
605119a934aSSatish Balay         sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
606119a934aSSatish Balay         sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
607119a934aSSatish Balay         sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
608119a934aSSatish Balay         v += 9;
609119a934aSSatish Balay       }
6101a6a6d98SBarry Smith       z[0] = sum1; z[1] = sum2; z[2] = sum3;
611119a934aSSatish Balay       z += 3;
612119a934aSSatish Balay     }
613119a934aSSatish Balay     break;
614119a934aSSatish Balay   case 4:
615119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
616119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
617119a934aSSatish Balay       sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
618119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
619119a934aSSatish Balay         xb = x + 4*(*idx++);
620119a934aSSatish Balay         x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
621119a934aSSatish Balay         sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
622119a934aSSatish Balay         sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
623119a934aSSatish Balay         sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
624119a934aSSatish Balay         sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
625119a934aSSatish Balay         v += 16;
626119a934aSSatish Balay       }
6271a6a6d98SBarry Smith       z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
628119a934aSSatish Balay       z += 4;
629119a934aSSatish Balay     }
630119a934aSSatish Balay     break;
631119a934aSSatish Balay   case 5:
632119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
633119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
634119a934aSSatish Balay       sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
635119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
636119a934aSSatish Balay         xb = x + 5*(*idx++);
637119a934aSSatish Balay         x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
638119a934aSSatish Balay         sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
639119a934aSSatish Balay         sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
640119a934aSSatish Balay         sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
641119a934aSSatish Balay         sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
642119a934aSSatish Balay         sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
643119a934aSSatish Balay         v += 25;
644119a934aSSatish Balay       }
6451a6a6d98SBarry Smith       z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
646119a934aSSatish Balay       z += 5;
647119a934aSSatish Balay     }
648119a934aSSatish Balay     break;
649119a934aSSatish Balay     /* block sizes larger then 5 by 5 are handled by BLAS */
650119a934aSSatish Balay   default: {
6511a6a6d98SBarry Smith       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt, _DZero = 0.0;
652119a934aSSatish Balay       if (!a->mult_work) {
6533b547af2SSatish Balay         k = PetscMax(a->m,a->n);
654bb42667fSSatish Balay         a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));
655119a934aSSatish Balay         CHKPTRQ(a->mult_work);
656119a934aSSatish Balay       }
657119a934aSSatish Balay       work = a->mult_work;
658119a934aSSatish Balay       for ( i=0; i<mbs; i++ ) {
659119a934aSSatish Balay         n     = ii[1] - ii[0]; ii++;
660119a934aSSatish Balay         ncols = n*bs;
661119a934aSSatish Balay         workt = work;
662119a934aSSatish Balay         for ( j=0; j<n; j++ ) {
663119a934aSSatish Balay           xb = x + bs*(*idx++);
664119a934aSSatish Balay           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
665119a934aSSatish Balay           workt += bs;
666119a934aSSatish Balay         }
6671a6a6d98SBarry Smith         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One);
668119a934aSSatish Balay         v += n*bs2;
669119a934aSSatish Balay         z += bs;
670119a934aSSatish Balay       }
671119a934aSSatish Balay     }
672119a934aSSatish Balay   }
6730419e6cdSSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
6740419e6cdSSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
6751a6a6d98SBarry Smith   PLogFlops(2*a->nz*bs2 - a->m);
676bb42667fSSatish Balay   return 0;
677bb42667fSSatish Balay }
678bb42667fSSatish Balay 
6791a6a6d98SBarry Smith static int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz)
680bb42667fSSatish Balay {
681bb42667fSSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6821a6a6d98SBarry Smith   Scalar          *xg,*zg,*zb;
683bb42667fSSatish Balay   register Scalar *x,*z,*v,*xb,x1,x2,x3,x4,x5;
684*bb1453f0SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
6859df24120SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
686119a934aSSatish Balay 
687119a934aSSatish Balay 
688bb42667fSSatish Balay   ierr = VecGetArray(xx,&xg); CHKERRQ(ierr); x = xg;
689bb42667fSSatish Balay   ierr = VecGetArray(zz,&zg); CHKERRQ(ierr); z = zg;
690*bb1453f0SSatish Balay   PetscMemzero(z,N*sizeof(Scalar));
691bb42667fSSatish Balay 
692119a934aSSatish Balay   idx   = a->j;
693119a934aSSatish Balay   v     = a->a;
694119a934aSSatish Balay   ii    = a->i;
695119a934aSSatish Balay 
696119a934aSSatish Balay   switch (bs) {
697119a934aSSatish Balay   case 1:
698119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
699119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
700119a934aSSatish Balay       xb = x + i; x1 = xb[0];
701119a934aSSatish Balay       ib = idx + ai[i];
702119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
703*bb1453f0SSatish Balay         rval    = ib[j];
704*bb1453f0SSatish Balay         z[rval] += *v++ * x1;
705119a934aSSatish Balay       }
706119a934aSSatish Balay     }
707119a934aSSatish Balay     break;
708119a934aSSatish Balay   case 2:
709119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
710119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
711119a934aSSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
712119a934aSSatish Balay       ib = idx + ai[i];
713119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
714119a934aSSatish Balay         rval      = ib[j]*2;
715*bb1453f0SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
716*bb1453f0SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
717119a934aSSatish Balay         v += 4;
718119a934aSSatish Balay       }
719119a934aSSatish Balay     }
720119a934aSSatish Balay     break;
721119a934aSSatish Balay   case 3:
722119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
723119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
724119a934aSSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
725119a934aSSatish Balay       ib = idx + ai[i];
726119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
727119a934aSSatish Balay         rval      = ib[j]*3;
728*bb1453f0SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
729*bb1453f0SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
730*bb1453f0SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
731119a934aSSatish Balay         v += 9;
732119a934aSSatish Balay       }
733119a934aSSatish Balay     }
734119a934aSSatish Balay     break;
735119a934aSSatish Balay   case 4:
736119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
737119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
738119a934aSSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
739119a934aSSatish Balay       ib = idx + ai[i];
740119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
741119a934aSSatish Balay         rval      = ib[j]*4;
742*bb1453f0SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
743*bb1453f0SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
744*bb1453f0SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
745*bb1453f0SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
746119a934aSSatish Balay         v += 16;
747119a934aSSatish Balay       }
748119a934aSSatish Balay     }
749119a934aSSatish Balay     break;
750119a934aSSatish Balay   case 5:
751119a934aSSatish Balay     for ( i=0; i<mbs; i++ ) {
752119a934aSSatish Balay       n  = ii[1] - ii[0]; ii++;
753119a934aSSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
754119a934aSSatish Balay       x4 = xb[3];   x5 = xb[4];
755119a934aSSatish Balay       ib = idx + ai[i];
756119a934aSSatish Balay       for ( j=0; j<n; j++ ) {
757119a934aSSatish Balay         rval      = ib[j]*5;
758*bb1453f0SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
759*bb1453f0SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
760*bb1453f0SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
761*bb1453f0SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
762*bb1453f0SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
763119a934aSSatish Balay         v += 25;
764119a934aSSatish Balay       }
765119a934aSSatish Balay     }
766119a934aSSatish Balay     break;
767119a934aSSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
768119a934aSSatish Balay     default: {
769119a934aSSatish Balay       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
770119a934aSSatish Balay       if (!a->mult_work) {
7713b547af2SSatish Balay         k = PetscMax(a->m,a->n);
772bb42667fSSatish Balay         a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));
773119a934aSSatish Balay         CHKPTRQ(a->mult_work);
774119a934aSSatish Balay       }
775119a934aSSatish Balay       work = a->mult_work;
776119a934aSSatish Balay       for ( i=0; i<mbs; i++ ) {
777119a934aSSatish Balay         n     = ii[1] - ii[0]; ii++;
778119a934aSSatish Balay         ncols = n*bs;
779119a934aSSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
780119a934aSSatish Balay         LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One);
781119a934aSSatish Balay         v += n*bs2;
782119a934aSSatish Balay         x += bs;
783119a934aSSatish Balay         workt = work;
784119a934aSSatish Balay         for ( j=0; j<n; j++ ) {
785119a934aSSatish Balay           zb = z + bs*(*idx++);
786*bb1453f0SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
787119a934aSSatish Balay           workt += bs;
788119a934aSSatish Balay         }
789119a934aSSatish Balay       }
790119a934aSSatish Balay     }
791119a934aSSatish Balay   }
7920419e6cdSSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
7930419e6cdSSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
7942593348eSBarry Smith   return 0;
7952593348eSBarry Smith }
7962593348eSBarry Smith 
797de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
7982593348eSBarry Smith {
799b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
8009df24120SSatish Balay   if (nz)  *nz  = a->bs2*a->nz;
801bcd2baecSBarry Smith   if (nza) *nza = a->maxnz;
802bcd2baecSBarry Smith   if (mem) *mem = (int)A->mem;
8032593348eSBarry Smith   return 0;
8042593348eSBarry Smith }
8052593348eSBarry Smith 
80691d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
80791d316f6SSatish Balay {
80891d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
80991d316f6SSatish Balay 
81091d316f6SSatish Balay   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
81191d316f6SSatish Balay 
81291d316f6SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
81391d316f6SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
814a7c10996SSatish Balay       (a->nz != b->nz)) {
81591d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
81691d316f6SSatish Balay   }
81791d316f6SSatish Balay 
81891d316f6SSatish Balay   /* if the a->i are the same */
81991d316f6SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
82091d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
82191d316f6SSatish Balay   }
82291d316f6SSatish Balay 
82391d316f6SSatish Balay   /* if a->j are the same */
82491d316f6SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
82591d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
82691d316f6SSatish Balay   }
82791d316f6SSatish Balay 
82891d316f6SSatish Balay   /* if a->a are the same */
82991d316f6SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
83091d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
83191d316f6SSatish Balay   }
83291d316f6SSatish Balay   *flg = PETSC_TRUE;
83391d316f6SSatish Balay   return 0;
83491d316f6SSatish Balay 
83591d316f6SSatish Balay }
83691d316f6SSatish Balay 
83791d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
83891d316f6SSatish Balay {
83991d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
8407e67e3f9SSatish Balay   int         i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
84117e48fc4SSatish Balay   Scalar      *x, zero = 0.0,*aa,*aa_j;
84217e48fc4SSatish Balay 
84317e48fc4SSatish Balay   bs   = a->bs;
84417e48fc4SSatish Balay   aa   = a->a;
84517e48fc4SSatish Balay   ai   = a->i;
84617e48fc4SSatish Balay   aj   = a->j;
84717e48fc4SSatish Balay   ambs = a->mbs;
8489df24120SSatish Balay   bs2  = a->bs2;
84991d316f6SSatish Balay 
85091d316f6SSatish Balay   VecSet(&zero,v);
85191d316f6SSatish Balay   VecGetArray(v,&x); VecGetLocalSize(v,&n);
8529867e207SSatish Balay   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector");
85317e48fc4SSatish Balay   for ( i=0; i<ambs; i++ ) {
85417e48fc4SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
85517e48fc4SSatish Balay       if (aj[j] == i) {
85617e48fc4SSatish Balay         row  = i*bs;
8577e67e3f9SSatish Balay         aa_j = aa+j*bs2;
8587e67e3f9SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
85991d316f6SSatish Balay         break;
86091d316f6SSatish Balay       }
86191d316f6SSatish Balay     }
86291d316f6SSatish Balay   }
86391d316f6SSatish Balay   return 0;
86491d316f6SSatish Balay }
86591d316f6SSatish Balay 
8669867e207SSatish Balay static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
8679867e207SSatish Balay {
8689867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
8699867e207SSatish Balay   Scalar      *l,*r,x,*v,*aa,*li,*ri;
8707e67e3f9SSatish Balay   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
8719867e207SSatish Balay 
8729867e207SSatish Balay   ai  = a->i;
8739867e207SSatish Balay   aj  = a->j;
8749867e207SSatish Balay   aa  = a->a;
8759867e207SSatish Balay   m   = a->m;
8769867e207SSatish Balay   n   = a->n;
8779867e207SSatish Balay   bs  = a->bs;
8789867e207SSatish Balay   mbs = a->mbs;
8799df24120SSatish Balay   bs2 = a->bs2;
8809867e207SSatish Balay   if (ll) {
8819867e207SSatish Balay     VecGetArray(ll,&l); VecGetSize(ll,&lm);
8829867e207SSatish Balay     if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length");
8839867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
8849867e207SSatish Balay       M  = ai[i+1] - ai[i];
8859867e207SSatish Balay       li = l + i*bs;
8867e67e3f9SSatish Balay       v  = aa + bs2*ai[i];
8879867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
8887e67e3f9SSatish Balay         for ( k=0; k<bs2; k++ ) {
8899867e207SSatish Balay           (*v++) *= li[k%bs];
8909867e207SSatish Balay         }
8919867e207SSatish Balay       }
8929867e207SSatish Balay     }
8939867e207SSatish Balay   }
8949867e207SSatish Balay 
8959867e207SSatish Balay   if (rr) {
8969867e207SSatish Balay     VecGetArray(rr,&r); VecGetSize(rr,&rn);
8979867e207SSatish Balay     if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length");
8989867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
8999867e207SSatish Balay       M  = ai[i+1] - ai[i];
9007e67e3f9SSatish Balay       v  = aa + bs2*ai[i];
9019867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
9029867e207SSatish Balay         ri = r + bs*aj[ai[i]+j];
9039867e207SSatish Balay         for ( k=0; k<bs; k++ ) {
9049867e207SSatish Balay           x = ri[k];
9059867e207SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
9069867e207SSatish Balay         }
9079867e207SSatish Balay       }
9089867e207SSatish Balay     }
9099867e207SSatish Balay   }
9109867e207SSatish Balay   return 0;
9119867e207SSatish Balay }
9129867e207SSatish Balay 
9139867e207SSatish Balay 
914b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
915b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
9162a38ea04SSatish Balay extern int MatIncreaseOverlap_SeqBAIJ(Mat,int,IS*,int);
917736121d4SSatish Balay extern int MatGetSubMatrix_SeqBAIJ(Mat,IS,IS,MatGetSubMatrixCall,Mat*);
918736121d4SSatish Balay extern int MatGetSubMatrices_SeqBAIJ(Mat,int,IS*,IS*,MatGetSubMatrixCall,Mat**);
9191a6a6d98SBarry Smith 
9201a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_N(Mat,Vec,Vec);
9211a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_1(Mat,Vec,Vec);
9221a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_2(Mat,Vec,Vec);
9231a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_3(Mat,Vec,Vec);
9241a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_4(Mat,Vec,Vec);
9251a6a6d98SBarry Smith extern int MatSolve_SeqBAIJ_5(Mat,Vec,Vec);
9261a6a6d98SBarry Smith 
9277fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
9287fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
9297fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
9307fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
9317fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
9327fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
9332593348eSBarry Smith 
934b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
9352593348eSBarry Smith {
936b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9372593348eSBarry Smith   Scalar      *v = a->a;
9382593348eSBarry Smith   double      sum = 0.0;
9399df24120SSatish Balay   int         i,nz=a->nz,bs2=a->bs2;
9402593348eSBarry Smith 
9412593348eSBarry Smith   if (type == NORM_FROBENIUS) {
9420419e6cdSSatish Balay     for (i=0; i< bs2*nz; i++ ) {
9432593348eSBarry Smith #if defined(PETSC_COMPLEX)
9442593348eSBarry Smith       sum += real(conj(*v)*(*v)); v++;
9452593348eSBarry Smith #else
9462593348eSBarry Smith       sum += (*v)*(*v); v++;
9472593348eSBarry Smith #endif
9482593348eSBarry Smith     }
9492593348eSBarry Smith     *norm = sqrt(sum);
9502593348eSBarry Smith   }
9512593348eSBarry Smith   else {
952b6490206SBarry Smith     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
9532593348eSBarry Smith   }
9542593348eSBarry Smith   return 0;
9552593348eSBarry Smith }
9562593348eSBarry Smith 
9572593348eSBarry Smith /*
9582593348eSBarry Smith      note: This can only work for identity for row and col. It would
9592593348eSBarry Smith    be good to check this and otherwise generate an error.
9602593348eSBarry Smith */
961b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
9622593348eSBarry Smith {
963b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
9642593348eSBarry Smith   Mat         outA;
965de6a44a3SBarry Smith   int         ierr;
9662593348eSBarry Smith 
967b6490206SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
9682593348eSBarry Smith 
9692593348eSBarry Smith   outA          = inA;
9702593348eSBarry Smith   inA->factor   = FACTOR_LU;
9712593348eSBarry Smith   a->row        = row;
9722593348eSBarry Smith   a->col        = col;
9732593348eSBarry Smith 
974de6a44a3SBarry Smith   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
9752593348eSBarry Smith 
9762593348eSBarry Smith   if (!a->diag) {
977b6490206SBarry Smith     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
9782593348eSBarry Smith   }
9797fc0212eSBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
9802593348eSBarry Smith   return 0;
9812593348eSBarry Smith }
9822593348eSBarry Smith 
983b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
9842593348eSBarry Smith {
985b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
9869df24120SSatish Balay   int         one = 1, totalnz = a->bs2*a->nz;
987b6490206SBarry Smith   BLscal_( &totalnz, alpha, a->a, &one );
988b6490206SBarry Smith   PLogFlops(totalnz);
9892593348eSBarry Smith   return 0;
9902593348eSBarry Smith }
9912593348eSBarry Smith 
9927e67e3f9SSatish Balay static int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
9937e67e3f9SSatish Balay {
9947e67e3f9SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
9957e67e3f9SSatish Balay   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
996a7c10996SSatish Balay   int        *ai = a->i, *ailen = a->ilen;
9979df24120SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
9987e67e3f9SSatish Balay   Scalar     *ap, *aa = a->a, zero = 0.0;
9997e67e3f9SSatish Balay 
10007e67e3f9SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over rows */
10017e67e3f9SSatish Balay     row  = im[k]; brow = row/bs;
10027e67e3f9SSatish Balay     if (row < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative row");
10037e67e3f9SSatish Balay     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqBAIJ:Row too large");
1004a7c10996SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
10057e67e3f9SSatish Balay     nrow = ailen[brow];
10067e67e3f9SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over columns */
10077e67e3f9SSatish Balay       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative column");
10087e67e3f9SSatish Balay       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqBAIJ:Column too large");
1009a7c10996SSatish Balay       col  = in[l] ;
10107e67e3f9SSatish Balay       bcol = col/bs;
10117e67e3f9SSatish Balay       cidx = col%bs;
10127e67e3f9SSatish Balay       ridx = row%bs;
10137e67e3f9SSatish Balay       high = nrow;
10147e67e3f9SSatish Balay       low  = 0; /* assume unsorted */
10157e67e3f9SSatish Balay       while (high-low > 5) {
10167e67e3f9SSatish Balay         t = (low+high)/2;
10177e67e3f9SSatish Balay         if (rp[t] > bcol) high = t;
10187e67e3f9SSatish Balay         else             low  = t;
10197e67e3f9SSatish Balay       }
10207e67e3f9SSatish Balay       for ( i=low; i<high; i++ ) {
10217e67e3f9SSatish Balay         if (rp[i] > bcol) break;
10227e67e3f9SSatish Balay         if (rp[i] == bcol) {
10237e67e3f9SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
10247e67e3f9SSatish Balay           goto finished;
10257e67e3f9SSatish Balay         }
10267e67e3f9SSatish Balay       }
10277e67e3f9SSatish Balay       *v++ = zero;
10287e67e3f9SSatish Balay       finished:;
10297e67e3f9SSatish Balay     }
10307e67e3f9SSatish Balay   }
10317e67e3f9SSatish Balay   return 0;
10327e67e3f9SSatish Balay }
10337e67e3f9SSatish Balay 
10342593348eSBarry Smith /* -------------------------------------------------------------------*/
1035cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
10369867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
10371a6a6d98SBarry Smith        MatMult_SeqBAIJ,0,
10381a6a6d98SBarry Smith        MatMultTrans_SeqBAIJ,0,
10391a6a6d98SBarry Smith        MatSolve_SeqBAIJ_N,0,
1040b6490206SBarry Smith        0,0,
1041de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
1042b6490206SBarry Smith        0,
1043f2501298SSatish Balay        MatTranspose_SeqBAIJ,
104417e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
10459867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
1046584200bdSSatish Balay        0,MatAssemblyEnd_SeqBAIJ,
1047b6490206SBarry Smith        0,
1048b6490206SBarry Smith        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
1049b6490206SBarry Smith        MatGetReordering_SeqBAIJ,
10507fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
1051b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
1052de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
1053b6490206SBarry Smith        0,0,/*  MatConvert_SeqBAIJ  */ 0,
1054736121d4SSatish Balay        MatGetSubMatrix_SeqBAIJ,0,
1055b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
1056b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
10577c7cbea0SSatish Balay        0,MatIncreaseOverlap_SeqBAIJ,
10587e67e3f9SSatish Balay        MatGetValues_SeqBAIJ,0,
10591a6a6d98SBarry Smith        0,MatScale_SeqBAIJ,
1060b6490206SBarry Smith        0};
10612593348eSBarry Smith 
10622593348eSBarry Smith /*@C
106344cd7ae7SLois Curfman McInnes    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
106444cd7ae7SLois Curfman McInnes    compressed row) format.  For good matrix assembly performance the
106544cd7ae7SLois Curfman McInnes    user should preallocate the matrix storage by setting the parameter nz
10662593348eSBarry Smith    (or nzz).  By setting these parameters accurately, performance can be
10672593348eSBarry Smith    increased by more than a factor of 50.
10682593348eSBarry Smith 
10692593348eSBarry Smith    Input Parameters:
10702593348eSBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
1071b6490206SBarry Smith .  bs - size of block
10722593348eSBarry Smith .  m - number of rows
10732593348eSBarry Smith .  n - number of columns
1074b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1075b6490206SBarry Smith .  nzz - number of block nonzeros per block row or PETSC_NULL
1076b6490206SBarry Smith          (possibly different for each row)
10772593348eSBarry Smith 
10782593348eSBarry Smith    Output Parameter:
10792593348eSBarry Smith .  A - the matrix
10802593348eSBarry Smith 
10812593348eSBarry Smith    Notes:
108244cd7ae7SLois Curfman McInnes    The block AIJ format is fully compatible with standard Fortran 77
10832593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
108444cd7ae7SLois Curfman McInnes    either one (as in Fortran) or zero.  See the users' manual for details.
10852593348eSBarry Smith 
10862593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
10872593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
10882593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
10892593348eSBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
10902593348eSBarry Smith 
109144cd7ae7SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
10922593348eSBarry Smith @*/
1093b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
10942593348eSBarry Smith {
10952593348eSBarry Smith   Mat         B;
1096b6490206SBarry Smith   Mat_SeqBAIJ *b;
1097f2501298SSatish Balay   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs;
1098b6490206SBarry Smith 
1099f2501298SSatish Balay   if (mbs*bs!=m || nbs*bs!=n)
1100f2501298SSatish Balay     SETERRQ(1,"MatCreateSeqBAIJ:Number rows, cols must be divisible by blocksize");
11012593348eSBarry Smith 
11022593348eSBarry Smith   *A = 0;
1103b6490206SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
11042593348eSBarry Smith   PLogObjectCreate(B);
1105b6490206SBarry Smith   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
110644cd7ae7SLois Curfman McInnes   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
11072593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
11081a6a6d98SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg); CHKERRQ(ierr);
11091a6a6d98SBarry Smith   if (!flg) {
11107fc0212eSBarry Smith     switch (bs) {
11117fc0212eSBarry Smith       case 1:
11127fc0212eSBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
11131a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_1;
11147fc0212eSBarry Smith         break;
11154eeb42bcSBarry Smith       case 2:
11164eeb42bcSBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
11171a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_2;
11186d84be18SBarry Smith         break;
1119f327dd97SBarry Smith       case 3:
1120f327dd97SBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
11211a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_3;
11224eeb42bcSBarry Smith         break;
11236d84be18SBarry Smith       case 4:
11246d84be18SBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
11251a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_4;
11266d84be18SBarry Smith         break;
11276d84be18SBarry Smith       case 5:
11286d84be18SBarry Smith         B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
11291a6a6d98SBarry Smith         B->ops.solve           = MatSolve_SeqBAIJ_5;
11306d84be18SBarry Smith         break;
11317fc0212eSBarry Smith     }
11321a6a6d98SBarry Smith   }
1133b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
1134b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
11352593348eSBarry Smith   B->factor           = 0;
11362593348eSBarry Smith   B->lupivotthreshold = 1.0;
11372593348eSBarry Smith   b->row              = 0;
11382593348eSBarry Smith   b->col              = 0;
11392593348eSBarry Smith   b->reallocs         = 0;
11402593348eSBarry Smith 
114144cd7ae7SLois Curfman McInnes   b->m       = m; B->m = m; B->M = m;
114244cd7ae7SLois Curfman McInnes   b->n       = n; B->n = n; B->N = n;
1143b6490206SBarry Smith   b->mbs     = mbs;
1144f2501298SSatish Balay   b->nbs     = nbs;
1145b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
11462593348eSBarry Smith   if (nnz == PETSC_NULL) {
1147119a934aSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
11482593348eSBarry Smith     else if (nz <= 0)        nz = 1;
1149b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1150b6490206SBarry Smith     nz = nz*mbs;
11512593348eSBarry Smith   }
11522593348eSBarry Smith   else {
11532593348eSBarry Smith     nz = 0;
1154b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
11552593348eSBarry Smith   }
11562593348eSBarry Smith 
11572593348eSBarry Smith   /* allocate the matrix space */
11587e67e3f9SSatish Balay   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
11592593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
11607e67e3f9SSatish Balay   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
11617e67e3f9SSatish Balay   b->j  = (int *) (b->a + nz*bs2);
11622593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
11632593348eSBarry Smith   b->i  = b->j + nz;
11642593348eSBarry Smith   b->singlemalloc = 1;
11652593348eSBarry Smith 
1166de6a44a3SBarry Smith   b->i[0] = 0;
1167b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
11682593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
11692593348eSBarry Smith   }
11702593348eSBarry Smith 
1171b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
1172b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1173b6490206SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
1174b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
11752593348eSBarry Smith 
1176b6490206SBarry Smith   b->bs               = bs;
11779df24120SSatish Balay   b->bs2              = bs2;
1178b6490206SBarry Smith   b->mbs              = mbs;
11792593348eSBarry Smith   b->nz               = 0;
11802593348eSBarry Smith   b->maxnz            = nz;
11812593348eSBarry Smith   b->sorted           = 0;
11822593348eSBarry Smith   b->roworiented      = 1;
11832593348eSBarry Smith   b->nonew            = 0;
11842593348eSBarry Smith   b->diag             = 0;
11852593348eSBarry Smith   b->solve_work       = 0;
1186de6a44a3SBarry Smith   b->mult_work        = 0;
11872593348eSBarry Smith   b->spptr            = 0;
11882593348eSBarry Smith   *A = B;
11892593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
11902593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
11912593348eSBarry Smith   return 0;
11922593348eSBarry Smith }
11932593348eSBarry Smith 
1194b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
11952593348eSBarry Smith {
11962593348eSBarry Smith   Mat         C;
1197b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
11989df24120SSatish Balay   int         i,len, mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
1199de6a44a3SBarry Smith 
1200de6a44a3SBarry Smith   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
12012593348eSBarry Smith 
12022593348eSBarry Smith   *B = 0;
1203b6490206SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
12042593348eSBarry Smith   PLogObjectCreate(C);
1205b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
12062593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1207b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
1208b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
12092593348eSBarry Smith   C->factor     = A->factor;
12102593348eSBarry Smith   c->row        = 0;
12112593348eSBarry Smith   c->col        = 0;
12122593348eSBarry Smith   C->assembled  = PETSC_TRUE;
12132593348eSBarry Smith 
121444cd7ae7SLois Curfman McInnes   c->m = C->m   = a->m;
121544cd7ae7SLois Curfman McInnes   c->n = C->n   = a->n;
121644cd7ae7SLois Curfman McInnes   C->M          = a->m;
121744cd7ae7SLois Curfman McInnes   C->N          = a->n;
121844cd7ae7SLois Curfman McInnes 
1219b6490206SBarry Smith   c->bs         = a->bs;
12209df24120SSatish Balay   c->bs2        = a->bs2;
1221b6490206SBarry Smith   c->mbs        = a->mbs;
1222b6490206SBarry Smith   c->nbs        = a->nbs;
12232593348eSBarry Smith 
1224b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1225b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1226b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
12272593348eSBarry Smith     c->imax[i] = a->imax[i];
12282593348eSBarry Smith     c->ilen[i] = a->ilen[i];
12292593348eSBarry Smith   }
12302593348eSBarry Smith 
12312593348eSBarry Smith   /* allocate the matrix space */
12322593348eSBarry Smith   c->singlemalloc = 1;
12337e67e3f9SSatish Balay   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
12342593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
12357e67e3f9SSatish Balay   c->j  = (int *) (c->a + nz*bs2);
1236de6a44a3SBarry Smith   c->i  = c->j + nz;
1237b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1238b6490206SBarry Smith   if (mbs > 0) {
1239de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
12402593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
12417e67e3f9SSatish Balay       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
12422593348eSBarry Smith     }
12432593348eSBarry Smith   }
12442593348eSBarry Smith 
1245b6490206SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
12462593348eSBarry Smith   c->sorted      = a->sorted;
12472593348eSBarry Smith   c->roworiented = a->roworiented;
12482593348eSBarry Smith   c->nonew       = a->nonew;
12492593348eSBarry Smith 
12502593348eSBarry Smith   if (a->diag) {
1251b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1252b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1253b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
12542593348eSBarry Smith       c->diag[i] = a->diag[i];
12552593348eSBarry Smith     }
12562593348eSBarry Smith   }
12572593348eSBarry Smith   else c->diag          = 0;
12582593348eSBarry Smith   c->nz                 = a->nz;
12592593348eSBarry Smith   c->maxnz              = a->maxnz;
12602593348eSBarry Smith   c->solve_work         = 0;
12612593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
12627fc0212eSBarry Smith   c->mult_work          = 0;
12632593348eSBarry Smith   *B = C;
12642593348eSBarry Smith   return 0;
12652593348eSBarry Smith }
12662593348eSBarry Smith 
126719bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
12682593348eSBarry Smith {
1269b6490206SBarry Smith   Mat_SeqBAIJ  *a;
12702593348eSBarry Smith   Mat          B;
1271de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1272b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
127335aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1274a7c10996SSatish Balay   int          *masked, nmask,tmp,bs2,ishift;
1275b6490206SBarry Smith   Scalar       *aa;
127619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
12772593348eSBarry Smith 
12781a6a6d98SBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,&flg);CHKERRQ(ierr);
1279de6a44a3SBarry Smith   bs2  = bs*bs;
1280b6490206SBarry Smith 
12812593348eSBarry Smith   MPI_Comm_size(comm,&size);
1282b6490206SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
128390ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
128477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1285de6a44a3SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
12862593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
12872593348eSBarry Smith 
1288b6490206SBarry Smith   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
128935aab85fSBarry Smith 
129035aab85fSBarry Smith   /*
129135aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
129235aab85fSBarry Smith     divisible by the blocksize
129335aab85fSBarry Smith   */
1294b6490206SBarry Smith   mbs        = M/bs;
129535aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
129635aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
129735aab85fSBarry Smith   else                  mbs++;
129835aab85fSBarry Smith   if (extra_rows) {
12991a6a6d98SBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize");
130035aab85fSBarry Smith   }
1301b6490206SBarry Smith 
13022593348eSBarry Smith   /* read in row lengths */
130335aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
130477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
130535aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
13062593348eSBarry Smith 
1307b6490206SBarry Smith   /* read in column indices */
130835aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
130977c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
131035aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1311b6490206SBarry Smith 
1312b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1313b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1314b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
131535aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
131635aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
131735aab85fSBarry Smith   masked = mask + mbs;
1318b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1319b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
132035aab85fSBarry Smith     nmask = 0;
1321b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1322b6490206SBarry Smith       kmax = rowlengths[rowcount];
1323b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
132435aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
132535aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1326b6490206SBarry Smith       }
1327b6490206SBarry Smith       rowcount++;
1328b6490206SBarry Smith     }
132935aab85fSBarry Smith     browlengths[i] += nmask;
133035aab85fSBarry Smith     /* zero out the mask elements we set */
133135aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1332b6490206SBarry Smith   }
1333b6490206SBarry Smith 
13342593348eSBarry Smith   /* create our matrix */
133535aab85fSBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
133635aab85fSBarry Smith          CHKERRQ(ierr);
13372593348eSBarry Smith   B = *A;
1338b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
13392593348eSBarry Smith 
13402593348eSBarry Smith   /* set matrix "i" values */
1341de6a44a3SBarry Smith   a->i[0] = 0;
1342b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1343b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1344b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
13452593348eSBarry Smith   }
1346b6490206SBarry Smith   a->nz         = 0;
1347b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
13482593348eSBarry Smith 
1349b6490206SBarry Smith   /* read in nonzero values */
135035aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
135177c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
135235aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1353b6490206SBarry Smith 
1354b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1355b6490206SBarry Smith   nzcount = 0; jcount = 0;
1356b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1357b6490206SBarry Smith     nzcountb = nzcount;
135835aab85fSBarry Smith     nmask    = 0;
1359b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1360b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1361b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
136235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
136335aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1364b6490206SBarry Smith       }
1365b6490206SBarry Smith       rowcount++;
1366b6490206SBarry Smith     }
1367de6a44a3SBarry Smith     /* sort the masked values */
136877c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1369de6a44a3SBarry Smith 
1370b6490206SBarry Smith     /* set "j" values into matrix */
1371b6490206SBarry Smith     maskcount = 1;
137235aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
137335aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1374de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1375b6490206SBarry Smith     }
1376b6490206SBarry Smith     /* set "a" values into matrix */
1377de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1378b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1379b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1380b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1381de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
1382de6a44a3SBarry Smith         block  = mask[tmp] - 1;
1383de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
1384de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1385b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1386b6490206SBarry Smith       }
1387b6490206SBarry Smith     }
138835aab85fSBarry Smith     /* zero out the mask elements we set */
138935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1390b6490206SBarry Smith   }
139135aab85fSBarry Smith   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
1392b6490206SBarry Smith 
1393b6490206SBarry Smith   PetscFree(rowlengths);
1394b6490206SBarry Smith   PetscFree(browlengths);
1395b6490206SBarry Smith   PetscFree(aa);
1396b6490206SBarry Smith   PetscFree(jj);
1397b6490206SBarry Smith   PetscFree(mask);
1398b6490206SBarry Smith 
1399b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1400de6a44a3SBarry Smith 
1401de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
1402de6a44a3SBarry Smith   if (flg) {
140319bcc07fSBarry Smith     Viewer tviewer;
140419bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
140590ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr);
140619bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
140719bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1408de6a44a3SBarry Smith   }
1409de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
1410de6a44a3SBarry Smith   if (flg) {
141119bcc07fSBarry Smith     Viewer tviewer;
141219bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
141390ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
141419bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
141519bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1416de6a44a3SBarry Smith   }
1417de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
1418de6a44a3SBarry Smith   if (flg) {
141919bcc07fSBarry Smith     Viewer tviewer;
142019bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
142119bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
142219bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1423de6a44a3SBarry Smith   }
1424de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
1425de6a44a3SBarry Smith   if (flg) {
142619bcc07fSBarry Smith     Viewer tviewer;
142719bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
142890ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr);
142919bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
143019bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1431de6a44a3SBarry Smith   }
1432de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
1433de6a44a3SBarry Smith   if (flg) {
143419bcc07fSBarry Smith     Viewer tviewer;
143519bcc07fSBarry Smith     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
143619bcc07fSBarry Smith     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
143719bcc07fSBarry Smith     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
143819bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1439de6a44a3SBarry Smith   }
14402593348eSBarry Smith   return 0;
14412593348eSBarry Smith }
14422593348eSBarry Smith 
14432593348eSBarry Smith 
14442593348eSBarry Smith 
1445