xref: /petsc/src/mat/impls/baij/seq/baij.c (revision cd0e1443fd6b0c0b4585f16770ffe0b64fb28823)
1b6490206SBarry Smith 
22593348eSBarry Smith #ifndef lint
3*cd0e1443SSatish Balay static char vcid[] = "$Id: baij.c,v 1.24 1996/03/31 19:59:04 curfman 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"
112593348eSBarry Smith #include "vec/vecimpl.h"
122593348eSBarry Smith #include "inline/spops.h"
132593348eSBarry Smith #include "petsc.h"
142593348eSBarry Smith 
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();
41bcd2baecSBarry Smith   ishift = a->indexshift;
42bcd2baecSBarry Smith   oshift = -MatReorderingIndexShift[(int)type];
43bcd2baecSBarry Smith   if (MatReorderingRequiresSymmetric[(int)type]) {
44bcd2baecSBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja);
45bcd2baecSBarry Smith     CHKERRQ(ierr);
46bcd2baecSBarry Smith     ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
472593348eSBarry Smith     PetscFree(ia); PetscFree(ja);
48bcd2baecSBarry Smith   } else {
49bcd2baecSBarry Smith     if (ishift == oshift) {
50bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
51bcd2baecSBarry Smith     }
52bcd2baecSBarry Smith     else if (ishift == -1) {
53bcd2baecSBarry Smith       /* temporarily subtract 1 from i and j indices */
54bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
55bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
56bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
57bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
58bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
59bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
60bcd2baecSBarry Smith     } else {
61bcd2baecSBarry Smith       /* temporarily add 1 to i and j indices */
62bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
63bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
64bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
65bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
66bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
67bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
68bcd2baecSBarry Smith     }
69bcd2baecSBarry Smith   }
702593348eSBarry Smith   return 0;
712593348eSBarry Smith }
722593348eSBarry Smith 
73de6a44a3SBarry Smith /*
74de6a44a3SBarry Smith      Adds diagonal pointers to sparse matrix structure.
75de6a44a3SBarry Smith */
76de6a44a3SBarry Smith 
77de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A)
78de6a44a3SBarry Smith {
79de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
807fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
81de6a44a3SBarry Smith 
82de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
83de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
847fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
85de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
86de6a44a3SBarry Smith       if (a->j[j] == i) {
87de6a44a3SBarry Smith         diag[i] = j;
88de6a44a3SBarry Smith         break;
89de6a44a3SBarry Smith       }
90de6a44a3SBarry Smith     }
91de6a44a3SBarry Smith   }
92de6a44a3SBarry Smith   a->diag = diag;
93de6a44a3SBarry Smith   return 0;
94de6a44a3SBarry Smith }
952593348eSBarry Smith 
962593348eSBarry Smith #include "draw.h"
972593348eSBarry Smith #include "pinclude/pviewer.h"
9877c4ece6SBarry Smith #include "sys.h"
992593348eSBarry Smith 
100b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
1012593348eSBarry Smith {
102b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
103b6490206SBarry Smith   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l;
104b6490206SBarry Smith   Scalar      *aa;
1052593348eSBarry Smith 
10690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1072593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
1082593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
1092593348eSBarry Smith   col_lens[1] = a->m;
1102593348eSBarry Smith   col_lens[2] = a->n;
111b6490206SBarry Smith   col_lens[3] = a->nz*bs*bs;
1122593348eSBarry Smith 
1132593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
114b6490206SBarry Smith   count = 0;
115b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
116b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
117b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
118b6490206SBarry Smith     }
1192593348eSBarry Smith   }
12077c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
1212593348eSBarry Smith   PetscFree(col_lens);
1222593348eSBarry Smith 
1232593348eSBarry Smith   /* store column indices (zero start index) */
124b6490206SBarry Smith   jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj);
125b6490206SBarry Smith   count = 0;
126b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
127b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
128b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
129b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
130b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
1312593348eSBarry Smith         }
1322593348eSBarry Smith       }
133b6490206SBarry Smith     }
134b6490206SBarry Smith   }
13577c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs*bs*a->nz,BINARY_INT,0); CHKERRQ(ierr);
136b6490206SBarry Smith   PetscFree(jj);
1372593348eSBarry Smith 
1382593348eSBarry Smith   /* store nonzero values */
139b6490206SBarry Smith   aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa);
140b6490206SBarry Smith   count = 0;
141b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
142b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
143b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
144b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
145b6490206SBarry Smith           aa[count++] = a->a[bs*bs*k + l*bs + j];
146b6490206SBarry Smith         }
147b6490206SBarry Smith       }
148b6490206SBarry Smith     }
149b6490206SBarry Smith   }
15077c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs*bs*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
151b6490206SBarry Smith   PetscFree(aa);
1522593348eSBarry Smith   return 0;
1532593348eSBarry Smith }
1542593348eSBarry Smith 
155b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
1562593348eSBarry Smith {
157b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
158de6a44a3SBarry Smith   int         ierr, i,j,format,bs = a->bs,k,l;
1592593348eSBarry Smith   FILE        *fd;
1602593348eSBarry Smith   char        *outputname;
1612593348eSBarry Smith 
16290ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
1632593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
16490ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
16590ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO) {
1662593348eSBarry Smith     /* no need to print additional information */ ;
1672593348eSBarry Smith   }
16890ace30eSBarry Smith   else if (format == ASCII_FORMAT_MATLAB) {
169b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported");
1702593348eSBarry Smith   }
1712593348eSBarry Smith   else {
172b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
173b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
174b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
175b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
176b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
17788685aaeSLois Curfman McInnes #if defined(PETSC_COMPLEX)
17888685aaeSLois Curfman McInnes           if (imag(a->a[bs*bs*k + l*bs + j]) != 0.0) {
17988685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
18088685aaeSLois Curfman McInnes               real(a->a[bs*bs*k + l*bs + j]),imag(a->a[bs*bs*k + l*bs + j]));
18188685aaeSLois Curfman McInnes           }
18288685aaeSLois Curfman McInnes           else {
18388685aaeSLois Curfman McInnes             fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs*bs*k + l*bs + j]));
18488685aaeSLois Curfman McInnes           }
18588685aaeSLois Curfman McInnes #else
186de6a44a3SBarry Smith             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs*bs*k + l*bs + j]);
18788685aaeSLois Curfman McInnes #endif
1882593348eSBarry Smith           }
1892593348eSBarry Smith         }
1902593348eSBarry Smith         fprintf(fd,"\n");
1912593348eSBarry Smith       }
1922593348eSBarry Smith     }
193b6490206SBarry Smith   }
1942593348eSBarry Smith   fflush(fd);
1952593348eSBarry Smith   return 0;
1962593348eSBarry Smith }
1972593348eSBarry Smith 
198b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
1992593348eSBarry Smith {
2002593348eSBarry Smith   Mat         A = (Mat) obj;
20119bcc07fSBarry Smith   ViewerType  vtype;
20219bcc07fSBarry Smith   int         ierr;
2032593348eSBarry Smith 
2042593348eSBarry Smith   if (!viewer) {
20519bcc07fSBarry Smith     viewer = STDOUT_VIEWER_SELF;
2062593348eSBarry Smith   }
20719bcc07fSBarry Smith 
20819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
20919bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
210b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported");
2112593348eSBarry Smith   }
21219bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
213b6490206SBarry Smith     return MatView_SeqBAIJ_ASCII(A,viewer);
2142593348eSBarry Smith   }
21519bcc07fSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
216b6490206SBarry Smith     return MatView_SeqBAIJ_Binary(A,viewer);
2172593348eSBarry Smith   }
21819bcc07fSBarry Smith   else if (vtype == DRAW_VIEWER) {
219b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported");
2202593348eSBarry Smith   }
2212593348eSBarry Smith   return 0;
2222593348eSBarry Smith }
223b6490206SBarry Smith 
224*cd0e1443SSatish Balay #define CHUNKSIZE   10
225*cd0e1443SSatish Balay 
226*cd0e1443SSatish Balay /* This version has row oriented v  */
227*cd0e1443SSatish Balay static int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
228*cd0e1443SSatish Balay {
229*cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
230*cd0e1443SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
231*cd0e1443SSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
232*cd0e1443SSatish Balay   int         *aj=a->j,nonew=a->nonew,shift=a->indexshift,bs=a->bs,brow,bcol;
233*cd0e1443SSatish Balay   int          ridx,cidx;
234*cd0e1443SSatish Balay   Scalar      *ap,value,*aa=a->a,*bap;
235*cd0e1443SSatish Balay 
236*cd0e1443SSatish Balay   for ( k=0; k<m; k++ ) { /* loop over added rows */
237*cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
238*cd0e1443SSatish Balay     if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row");
239*cd0e1443SSatish Balay     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large");
240*cd0e1443SSatish Balay     rp   = aj + ai[brow] + shift; ap = aa + bs*bs*ai[brow] + shift;
241*cd0e1443SSatish Balay     rmax = imax[brow]; nrow = ailen[brow];
242*cd0e1443SSatish Balay     low = 0;
243*cd0e1443SSatish Balay     for ( l=0; l<n; l++ ) { /* loop over added columns */
244*cd0e1443SSatish Balay       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column");
245*cd0e1443SSatish Balay       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large");
246*cd0e1443SSatish Balay       col = in[l] - shift; bcol = col/bs;
247*cd0e1443SSatish Balay       ridx = row % bs; cidx = col % bs;
248*cd0e1443SSatish Balay       if (roworiented) {
249*cd0e1443SSatish Balay         value = *v++;
250*cd0e1443SSatish Balay       }
251*cd0e1443SSatish Balay       else {
252*cd0e1443SSatish Balay         value = v[k + l*m];
253*cd0e1443SSatish Balay       }
254*cd0e1443SSatish Balay       if (!sorted) low = 0; high = nrow;
255*cd0e1443SSatish Balay       while (high-low > 5) {
256*cd0e1443SSatish Balay         t = (low+high)/2;
257*cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
258*cd0e1443SSatish Balay         else             low  = t;
259*cd0e1443SSatish Balay       }
260*cd0e1443SSatish Balay       for ( i=low; i<high; i++ ) {
261*cd0e1443SSatish Balay         if (rp[i] > bcol) break;
262*cd0e1443SSatish Balay         if (rp[i] == bcol) {
263*cd0e1443SSatish Balay           bap  = ap +  bs*bs*i + bs*ridx + cidx;
264*cd0e1443SSatish Balay           if (is == ADD_VALUES) *bap += value;
265*cd0e1443SSatish Balay           else                  *bap  = value;
266*cd0e1443SSatish Balay           goto noinsert;
267*cd0e1443SSatish Balay         }
268*cd0e1443SSatish Balay       }
269*cd0e1443SSatish Balay       if (nonew) goto noinsert;
270*cd0e1443SSatish Balay       if (nrow >= rmax) {
271*cd0e1443SSatish Balay         /* there is no extra room in row, therefore enlarge */
272*cd0e1443SSatish Balay         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
273*cd0e1443SSatish Balay         Scalar *new_a;
274*cd0e1443SSatish Balay 
275*cd0e1443SSatish Balay         /* malloc new storage space */
276*cd0e1443SSatish Balay         len     = new_nz*(sizeof(int)+bs*bs*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
277*cd0e1443SSatish Balay         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
278*cd0e1443SSatish Balay         new_j   = (int *) (new_a + bs*bs*new_nz);
279*cd0e1443SSatish Balay         new_i   = new_j + new_nz;
280*cd0e1443SSatish Balay 
281*cd0e1443SSatish Balay         /* copy over old data into new slots */
282*cd0e1443SSatish Balay         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
283*cd0e1443SSatish Balay         for ( ii=row+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
284*cd0e1443SSatish Balay         PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int));
285*cd0e1443SSatish Balay         len = (new_nz - CHUNKSIZE - ai[brow] - nrow - shift);
286*cd0e1443SSatish Balay         PetscMemcpy(new_j+ai[brow]+shift+nrow+CHUNKSIZE,aj+ai[brow]+shift+nrow,
287*cd0e1443SSatish Balay                                                            len*sizeof(int));
288*cd0e1443SSatish Balay         PetscMemcpy(new_a,aa,(ai[brow]+nrow+shift)*bs*bs*sizeof(Scalar));
289*cd0e1443SSatish Balay         PetscMemzero(new_a+bs*bs*(ai[brow]+shift+nrow),bs*bs*CHUNKSIZE);
290*cd0e1443SSatish Balay         PetscMemcpy(new_a+bs*bs*(ai[brow]+shift+nrow+CHUNKSIZE),
291*cd0e1443SSatish Balay                     aa+bs*bs*(ai[brow]+shift+nrow),bs*bs*len*sizeof(Scalar));
292*cd0e1443SSatish Balay         /* free up old matrix storage */
293*cd0e1443SSatish Balay         PetscFree(a->a);
294*cd0e1443SSatish Balay         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
295*cd0e1443SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
296*cd0e1443SSatish Balay         a->singlemalloc = 1;
297*cd0e1443SSatish Balay 
298*cd0e1443SSatish Balay         rp   = aj + ai[brow] + shift; ap = aa + ai[brow] + shift;
299*cd0e1443SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
300*cd0e1443SSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs*bs*sizeof(Scalar)));
301*cd0e1443SSatish Balay         a->maxnz += bs*bs*CHUNKSIZE;
302*cd0e1443SSatish Balay         a->reallocs++;
303*cd0e1443SSatish Balay       }
304*cd0e1443SSatish Balay       N = nrow++ - 1; a->nz++;
305*cd0e1443SSatish Balay       /* shift up all the later entries in this row */
306*cd0e1443SSatish Balay       for ( ii=N; ii>=i; ii-- ) {
307*cd0e1443SSatish Balay         rp[ii+1] = rp[ii];
308*cd0e1443SSatish Balay         PetscMemcpy(ap+bs*bs*(ii+1),ap+bs*bs*(ii),bs*bs);
309*cd0e1443SSatish Balay       }
310*cd0e1443SSatish Balay       rp[i] = bcol;
311*cd0e1443SSatish Balay       ap[bs*bs*i + bs*ridx + cidx] = value;
312*cd0e1443SSatish Balay       noinsert:;
313*cd0e1443SSatish Balay       low = i;
314*cd0e1443SSatish Balay     }
315*cd0e1443SSatish Balay     ailen[brow] = nrow;
316*cd0e1443SSatish Balay   }
317*cd0e1443SSatish Balay   return 0;
318*cd0e1443SSatish Balay }
319*cd0e1443SSatish Balay 
3200b824a48SSatish Balay static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
3210b824a48SSatish Balay {
3220b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3230b824a48SSatish Balay   *m = a->m; *n = a->n;
3240b824a48SSatish Balay   return 0;
3250b824a48SSatish Balay }
3260b824a48SSatish Balay 
3270b824a48SSatish Balay static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
3280b824a48SSatish Balay {
3290b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3300b824a48SSatish Balay   *m = 0; *n = a->m;
3310b824a48SSatish Balay   return 0;
3320b824a48SSatish Balay }
3330b824a48SSatish Balay 
3349867e207SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
3359867e207SSatish Balay {
3369867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3379867e207SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i;
3389867e207SSatish Balay   Scalar      *aa,*v_i,*aa_i;
3399867e207SSatish Balay 
3409867e207SSatish Balay   bs = a->bs;
3419867e207SSatish Balay   ai = a->i;
3429867e207SSatish Balay   aj = a->j;
3439867e207SSatish Balay   aa = a->a;
3449867e207SSatish Balay 
3459867e207SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range");
3469867e207SSatish Balay 
3479867e207SSatish Balay   bn  = row/bs;   /* Block number */
3489867e207SSatish Balay   bp  = row % bs; /* Block Position */
3499867e207SSatish Balay   M   = ai[bn+1] - ai[bn];
3509867e207SSatish Balay   *nz = bs*M;
3519867e207SSatish Balay 
3529867e207SSatish Balay   if (v) {
3539867e207SSatish Balay     *v = 0;
3549867e207SSatish Balay     if (*nz) {
3559867e207SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
3569867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
3579867e207SSatish Balay         v_i  = *v + i*bs;
3589867e207SSatish Balay         aa_i = aa + bs*bs*(ai[bn] + i);
3599867e207SSatish Balay         for ( j=bp,k=0; j<bs*bs; j+=bs,k++ ) {v_i[k] = aa_i[j];}
3609867e207SSatish Balay       }
3619867e207SSatish Balay     }
3629867e207SSatish Balay   }
3639867e207SSatish Balay 
3649867e207SSatish Balay   if (idx) {
3659867e207SSatish Balay     *idx = 0;
3669867e207SSatish Balay     if (*nz) {
3679867e207SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
3689867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
3699867e207SSatish Balay         idx_i = *idx + i*bs;
3709867e207SSatish Balay         itmp  = bs*aj[ai[bn] + i];
3719867e207SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
3729867e207SSatish Balay       }
3739867e207SSatish Balay     }
3749867e207SSatish Balay   }
3759867e207SSatish Balay   return 0;
3769867e207SSatish Balay }
3779867e207SSatish Balay 
3789867e207SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
3799867e207SSatish Balay {
3809867e207SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
3819867e207SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
3829867e207SSatish Balay   return 0;
3839867e207SSatish Balay }
384b6490206SBarry Smith 
385b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A)
3862593348eSBarry Smith {
387b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
388de6a44a3SBarry Smith   PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar));
3892593348eSBarry Smith   return 0;
3902593348eSBarry Smith }
3912593348eSBarry Smith 
392b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj)
3932593348eSBarry Smith {
3942593348eSBarry Smith   Mat         A  = (Mat) obj;
395b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3962593348eSBarry Smith 
3972593348eSBarry Smith #if defined(PETSC_LOG)
3982593348eSBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
3992593348eSBarry Smith #endif
4002593348eSBarry Smith   PetscFree(a->a);
4012593348eSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
4022593348eSBarry Smith   if (a->diag) PetscFree(a->diag);
4032593348eSBarry Smith   if (a->ilen) PetscFree(a->ilen);
4042593348eSBarry Smith   if (a->imax) PetscFree(a->imax);
4052593348eSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
406de6a44a3SBarry Smith   if (a->mult_work) PetscFree(a->mult_work);
4072593348eSBarry Smith   PetscFree(a);
408f2655603SLois Curfman McInnes   PLogObjectDestroy(A);
409f2655603SLois Curfman McInnes   PetscHeaderDestroy(A);
4102593348eSBarry Smith   return 0;
4112593348eSBarry Smith }
4122593348eSBarry Smith 
413b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
4142593348eSBarry Smith {
415b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4162593348eSBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
4172593348eSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
4182593348eSBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
4192593348eSBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
4202593348eSBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
4212593348eSBarry Smith   else if (op == ROWS_SORTED ||
4222593348eSBarry Smith            op == SYMMETRIC_MATRIX ||
4232593348eSBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
4242593348eSBarry Smith            op == YES_NEW_DIAGONALS)
42594a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
4262593348eSBarry Smith   else if (op == NO_NEW_DIAGONALS)
427b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
4282593348eSBarry Smith   else
429b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
4302593348eSBarry Smith   return 0;
4312593348eSBarry Smith }
4322593348eSBarry Smith 
4332593348eSBarry Smith 
4342593348eSBarry Smith /* -------------------------------------------------------*/
4352593348eSBarry Smith /* Should check that shapes of vectors and matrices match */
4362593348eSBarry Smith /* -------------------------------------------------------*/
437b6490206SBarry Smith #include "pinclude/plapack.h"
438b6490206SBarry Smith 
439b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy)
4402593348eSBarry Smith {
441b6490206SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
442b6490206SBarry Smith   Scalar          *xg,*yg;
443b6490206SBarry Smith   register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5;
444b6490206SBarry Smith   register Scalar x1,x2,x3,x4,x5;
445b6490206SBarry Smith   int             mbs = a->mbs, m = a->m, i, *idx,*ii;
446b6490206SBarry Smith   int             bs = a->bs,j,n,bs2 = bs*bs;
4472593348eSBarry Smith 
448b6490206SBarry Smith   VecGetArray(xx,&xg); x = xg;  VecGetArray(yy,&yg); y = yg;
449b6490206SBarry Smith   PetscMemzero(y,m*sizeof(Scalar));
450b6490206SBarry Smith   x     = x;
4512593348eSBarry Smith   idx   = a->j;
4522593348eSBarry Smith   v     = a->a;
4532593348eSBarry Smith   ii    = a->i;
454b6490206SBarry Smith 
455b6490206SBarry Smith   switch (bs) {
456b6490206SBarry Smith     case 1:
4572593348eSBarry Smith      for ( i=0; i<m; i++ ) {
4582593348eSBarry Smith        n    = ii[1] - ii[0]; ii++;
4592593348eSBarry Smith        sum  = 0.0;
4602593348eSBarry Smith        while (n--) sum += *v++ * x[*idx++];
4612593348eSBarry Smith        y[i] = sum;
4622593348eSBarry Smith       }
4632593348eSBarry Smith       break;
464b6490206SBarry Smith     case 2:
465b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
466de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
467b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0;
468b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
469b6490206SBarry Smith           xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
470b6490206SBarry Smith           sum1 += v[0]*x1 + v[2]*x2;
471b6490206SBarry Smith           sum2 += v[1]*x1 + v[3]*x2;
472b6490206SBarry Smith           v += 4;
473b6490206SBarry Smith         }
474b6490206SBarry Smith         y[0] += sum1; y[1] += sum2;
475b6490206SBarry Smith         y += 2;
476b6490206SBarry Smith       }
477b6490206SBarry Smith       break;
478b6490206SBarry Smith     case 3:
479b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
480de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
481b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
482b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
483b6490206SBarry Smith           xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
484b6490206SBarry Smith           sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
485b6490206SBarry Smith           sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
486b6490206SBarry Smith           sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
487b6490206SBarry Smith           v += 9;
488b6490206SBarry Smith         }
489b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3;
490b6490206SBarry Smith         y += 3;
491b6490206SBarry Smith       }
492b6490206SBarry Smith       break;
493b6490206SBarry Smith     case 4:
494b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
495de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
496b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
497b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
498b6490206SBarry Smith           xb = x + 4*(*idx++);
499b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
500b6490206SBarry Smith           sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
501b6490206SBarry Smith           sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
502b6490206SBarry Smith           sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
503b6490206SBarry Smith           sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
504b6490206SBarry Smith           v += 16;
505b6490206SBarry Smith         }
506b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4;
507b6490206SBarry Smith         y += 4;
508b6490206SBarry Smith       }
509b6490206SBarry Smith       break;
510b6490206SBarry Smith     case 5:
511b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
512de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
513b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
514b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
515b6490206SBarry Smith           xb = x + 5*(*idx++);
516b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
517b6490206SBarry Smith           sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
518b6490206SBarry Smith           sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
519b6490206SBarry Smith           sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
520b6490206SBarry Smith           sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
521b6490206SBarry Smith           sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
522b6490206SBarry Smith           v += 25;
523b6490206SBarry Smith         }
524b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5;
525b6490206SBarry Smith         y += 5;
526b6490206SBarry Smith       }
527b6490206SBarry Smith       break;
528b6490206SBarry Smith       /* block sizes larger then 5 by 5 are handled by BLAS */
529b6490206SBarry Smith     default: {
530de6a44a3SBarry Smith       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
531de6a44a3SBarry Smith       if (!a->mult_work) {
532de6a44a3SBarry Smith         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
533de6a44a3SBarry Smith         CHKPTRQ(a->mult_work);
534de6a44a3SBarry Smith       }
535de6a44a3SBarry Smith       work = a->mult_work;
536b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
537de6a44a3SBarry Smith         n     = ii[1] - ii[0]; ii++;
538de6a44a3SBarry Smith         ncols = n*bs;
539de6a44a3SBarry Smith         workt = work;
540b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
541b6490206SBarry Smith           xb = x + bs*(*idx++);
542de6a44a3SBarry Smith           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
543de6a44a3SBarry Smith           workt += bs;
544b6490206SBarry Smith         }
545de6a44a3SBarry Smith         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One);
546de6a44a3SBarry Smith         v += n*bs2;
547b6490206SBarry Smith         y += bs;
5482593348eSBarry Smith       }
5492593348eSBarry Smith     }
5502593348eSBarry Smith   }
551b6490206SBarry Smith   PLogFlops(2*bs2*a->nz - m);
5522593348eSBarry Smith   return 0;
5532593348eSBarry Smith }
5542593348eSBarry Smith 
555de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
5562593348eSBarry Smith {
557b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
558bcd2baecSBarry Smith   if (nz)  *nz  = a->bs*a->bs*a->nz;
559bcd2baecSBarry Smith   if (nza) *nza = a->maxnz;
560bcd2baecSBarry Smith   if (mem) *mem = (int)A->mem;
5612593348eSBarry Smith   return 0;
5622593348eSBarry Smith }
5632593348eSBarry Smith 
56491d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
56591d316f6SSatish Balay {
56691d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
56791d316f6SSatish Balay 
56891d316f6SSatish Balay   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
56991d316f6SSatish Balay 
57091d316f6SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
57191d316f6SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
57291d316f6SSatish Balay       (a->nz != b->nz)||(a->indexshift != b->indexshift)) {
57391d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
57491d316f6SSatish Balay   }
57591d316f6SSatish Balay 
57691d316f6SSatish Balay   /* if the a->i are the same */
57791d316f6SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
57891d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
57991d316f6SSatish Balay   }
58091d316f6SSatish Balay 
58191d316f6SSatish Balay   /* if a->j are the same */
58291d316f6SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
58391d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
58491d316f6SSatish Balay   }
58591d316f6SSatish Balay 
58691d316f6SSatish Balay   /* if a->a are the same */
58791d316f6SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
58891d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
58991d316f6SSatish Balay   }
59091d316f6SSatish Balay   *flg = PETSC_TRUE;
59191d316f6SSatish Balay   return 0;
59291d316f6SSatish Balay 
59391d316f6SSatish Balay }
59491d316f6SSatish Balay 
59591d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
59691d316f6SSatish Balay {
59791d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
59817e48fc4SSatish Balay   int         i,j,k,n,row,bs,*ai,*aj,ambs;
59917e48fc4SSatish Balay   Scalar      *x, zero = 0.0,*aa,*aa_j;
60017e48fc4SSatish Balay 
60117e48fc4SSatish Balay   bs  = a->bs;
60217e48fc4SSatish Balay   aa   = a->a;
60317e48fc4SSatish Balay   ai   = a->i;
60417e48fc4SSatish Balay   aj   = a->j;
60517e48fc4SSatish Balay   ambs = a->mbs;
60691d316f6SSatish Balay 
60791d316f6SSatish Balay   VecSet(&zero,v);
60891d316f6SSatish Balay   VecGetArray(v,&x); VecGetLocalSize(v,&n);
6099867e207SSatish Balay   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector");
61017e48fc4SSatish Balay   for ( i=0; i<ambs; i++ ) {
61117e48fc4SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
61217e48fc4SSatish Balay       if (aj[j] == i) {
61317e48fc4SSatish Balay         row  = i*bs;
61417e48fc4SSatish Balay         aa_j = aa+j*bs*bs;
61517e48fc4SSatish Balay         for (k=0; k<bs*bs; k+=(bs+1),row++) x[row] = aa_j[k];
61691d316f6SSatish Balay         break;
61791d316f6SSatish Balay       }
61891d316f6SSatish Balay     }
61991d316f6SSatish Balay   }
62091d316f6SSatish Balay   return 0;
62191d316f6SSatish Balay }
62291d316f6SSatish Balay 
6239867e207SSatish Balay static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
6249867e207SSatish Balay {
6259867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6269867e207SSatish Balay   Scalar      *l,*r,x,*v,*aa,*li,*ri;
6279867e207SSatish Balay   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs;
6289867e207SSatish Balay 
6299867e207SSatish Balay   ai  = a->i;
6309867e207SSatish Balay   aj  = a->j;
6319867e207SSatish Balay   aa  = a->a;
6329867e207SSatish Balay   m   = a->m;
6339867e207SSatish Balay   n   = a->n;
6349867e207SSatish Balay   bs  = a->bs;
6359867e207SSatish Balay   mbs = a->mbs;
6369867e207SSatish Balay 
6379867e207SSatish Balay   if (ll) {
6389867e207SSatish Balay     VecGetArray(ll,&l); VecGetSize(ll,&lm);
6399867e207SSatish Balay     if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length");
6409867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
6419867e207SSatish Balay       M  = ai[i+1] - ai[i];
6429867e207SSatish Balay       li = l + i*bs;
6439867e207SSatish Balay       v  = aa + bs*bs*ai[i];
6449867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
6459867e207SSatish Balay         for ( k=0; k<bs*bs; k++ ) {
6469867e207SSatish Balay           (*v++) *= li[k%bs];
6479867e207SSatish Balay         }
6489867e207SSatish Balay       }
6499867e207SSatish Balay     }
6509867e207SSatish Balay   }
6519867e207SSatish Balay 
6529867e207SSatish Balay   if (rr) {
6539867e207SSatish Balay     VecGetArray(rr,&r); VecGetSize(rr,&rn);
6549867e207SSatish Balay     if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length");
6559867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
6569867e207SSatish Balay       M  = ai[i+1] - ai[i];
6579867e207SSatish Balay       v  = aa + bs*bs*ai[i];
6589867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
6599867e207SSatish Balay         ri = r + bs*aj[ai[i]+j];
6609867e207SSatish Balay         for ( k=0; k<bs; k++ ) {
6619867e207SSatish Balay           x = ri[k];
6629867e207SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
6639867e207SSatish Balay         }
6649867e207SSatish Balay       }
6659867e207SSatish Balay     }
6669867e207SSatish Balay   }
6679867e207SSatish Balay   return 0;
6689867e207SSatish Balay }
6699867e207SSatish Balay 
6709867e207SSatish Balay 
671b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
672b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
673b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec);
6747fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
6757fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
6767fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
6777fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
6787fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
6797fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
6802593348eSBarry Smith 
681b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
6822593348eSBarry Smith {
683b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6842593348eSBarry Smith   Scalar      *v = a->a;
6852593348eSBarry Smith   double      sum = 0.0;
686b6490206SBarry Smith   int         i;
6872593348eSBarry Smith 
6882593348eSBarry Smith   if (type == NORM_FROBENIUS) {
6892593348eSBarry Smith     for (i=0; i<a->nz; i++ ) {
6902593348eSBarry Smith #if defined(PETSC_COMPLEX)
6912593348eSBarry Smith       sum += real(conj(*v)*(*v)); v++;
6922593348eSBarry Smith #else
6932593348eSBarry Smith       sum += (*v)*(*v); v++;
6942593348eSBarry Smith #endif
6952593348eSBarry Smith     }
6962593348eSBarry Smith     *norm = sqrt(sum);
6972593348eSBarry Smith   }
6982593348eSBarry Smith   else {
699b6490206SBarry Smith     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
7002593348eSBarry Smith   }
7012593348eSBarry Smith   return 0;
7022593348eSBarry Smith }
7032593348eSBarry Smith 
7042593348eSBarry Smith /*
7052593348eSBarry Smith      note: This can only work for identity for row and col. It would
7062593348eSBarry Smith    be good to check this and otherwise generate an error.
7072593348eSBarry Smith */
708b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
7092593348eSBarry Smith {
710b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
7112593348eSBarry Smith   Mat         outA;
712de6a44a3SBarry Smith   int         ierr;
7132593348eSBarry Smith 
714b6490206SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
7152593348eSBarry Smith 
7162593348eSBarry Smith   outA          = inA;
7172593348eSBarry Smith   inA->factor   = FACTOR_LU;
7182593348eSBarry Smith   a->row        = row;
7192593348eSBarry Smith   a->col        = col;
7202593348eSBarry Smith 
721de6a44a3SBarry Smith   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
7222593348eSBarry Smith 
7232593348eSBarry Smith   if (!a->diag) {
724b6490206SBarry Smith     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
7252593348eSBarry Smith   }
7267fc0212eSBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
7272593348eSBarry Smith   return 0;
7282593348eSBarry Smith }
7292593348eSBarry Smith 
730b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
7312593348eSBarry Smith {
732b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
733b6490206SBarry Smith   int         one = 1, totalnz = a->bs*a->bs*a->nz;
734b6490206SBarry Smith   BLscal_( &totalnz, alpha, a->a, &one );
735b6490206SBarry Smith   PLogFlops(totalnz);
7362593348eSBarry Smith   return 0;
7372593348eSBarry Smith }
7382593348eSBarry Smith 
739b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A)
7402593348eSBarry Smith {
7412593348eSBarry Smith   static int called = 0;
7422593348eSBarry Smith 
7432593348eSBarry Smith   if (called) return 0; else called = 1;
7442593348eSBarry Smith   return 0;
7452593348eSBarry Smith }
7462593348eSBarry Smith /* -------------------------------------------------------------------*/
747*cd0e1443SSatish Balay static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
7489867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
749b6490206SBarry Smith        MatMult_SeqBAIJ,0,
750b6490206SBarry Smith        0,0,
751de6a44a3SBarry Smith        MatSolve_SeqBAIJ,0,
752b6490206SBarry Smith        0,0,
753de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
754b6490206SBarry Smith        0,
755b6490206SBarry Smith        0,
75617e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
7579867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
758b6490206SBarry Smith        0,0,
759b6490206SBarry Smith        0,
760b6490206SBarry Smith        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
761b6490206SBarry Smith        MatGetReordering_SeqBAIJ,
7627fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
763b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
764de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
765b6490206SBarry Smith        0,0,/*  MatConvert_SeqBAIJ  */ 0,
766b6490206SBarry Smith        0,0,
767b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
768b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
769b6490206SBarry Smith        0,0,
770b6490206SBarry Smith        0,0,
771b6490206SBarry Smith        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
772b6490206SBarry Smith        0};
7732593348eSBarry Smith 
7742593348eSBarry Smith /*@C
775b6490206SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format
7762593348eSBarry Smith    (the default parallel PETSc format).  For good matrix assembly performance
7772593348eSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
7782593348eSBarry Smith    (or nzz).  By setting these parameters accurately, performance can be
7792593348eSBarry Smith    increased by more than a factor of 50.
7802593348eSBarry Smith 
7812593348eSBarry Smith    Input Parameters:
7822593348eSBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
783b6490206SBarry Smith .  bs - size of block
7842593348eSBarry Smith .  m - number of rows
7852593348eSBarry Smith .  n - number of columns
786b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
787b6490206SBarry Smith .  nzz - number of block nonzeros per block row or PETSC_NULL
788b6490206SBarry Smith          (possibly different for each row)
7892593348eSBarry Smith 
7902593348eSBarry Smith    Output Parameter:
7912593348eSBarry Smith .  A - the matrix
7922593348eSBarry Smith 
7932593348eSBarry Smith    Notes:
794b6490206SBarry Smith    The BAIJ format, is fully compatible with standard Fortran 77
7952593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
7962593348eSBarry Smith    either one (as in Fortran) or zero.  See the users manual for details.
7972593348eSBarry Smith 
7982593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
7992593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
8002593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
8012593348eSBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
8022593348eSBarry Smith 
8032593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
8042593348eSBarry Smith @*/
805b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
8062593348eSBarry Smith {
8072593348eSBarry Smith   Mat         B;
808b6490206SBarry Smith   Mat_SeqBAIJ *b;
809b6490206SBarry Smith   int         i,len,ierr, flg,mbs = m/bs;
810b6490206SBarry Smith 
811b6490206SBarry Smith   if (mbs*bs != m)
812b6490206SBarry Smith     SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize");
8132593348eSBarry Smith 
8142593348eSBarry Smith   *A      = 0;
815b6490206SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
8162593348eSBarry Smith   PLogObjectCreate(B);
817b6490206SBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
8182593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
8197fc0212eSBarry Smith   switch (bs) {
8207fc0212eSBarry Smith     case 1:
8217fc0212eSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
8227fc0212eSBarry Smith       break;
8234eeb42bcSBarry Smith     case 2:
8244eeb42bcSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
8256d84be18SBarry Smith       break;
826f327dd97SBarry Smith     case 3:
827f327dd97SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
8284eeb42bcSBarry Smith       break;
8296d84be18SBarry Smith     case 4:
8306d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
8316d84be18SBarry Smith       break;
8326d84be18SBarry Smith     case 5:
8336d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
8346d84be18SBarry Smith       break;
8357fc0212eSBarry Smith   }
836b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
837b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
8382593348eSBarry Smith   B->factor           = 0;
8392593348eSBarry Smith   B->lupivotthreshold = 1.0;
8402593348eSBarry Smith   b->row              = 0;
8412593348eSBarry Smith   b->col              = 0;
8422593348eSBarry Smith   b->reallocs         = 0;
8432593348eSBarry Smith 
8442593348eSBarry Smith   b->m       = m;
845b6490206SBarry Smith   b->mbs     = mbs;
8462593348eSBarry Smith   b->n       = n;
847b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
8482593348eSBarry Smith   if (nnz == PETSC_NULL) {
849de6a44a3SBarry Smith     if (nz == PETSC_DEFAULT) nz = 5;
8502593348eSBarry Smith     else if (nz <= 0)        nz = 1;
851b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
852b6490206SBarry Smith     nz = nz*mbs;
8532593348eSBarry Smith   }
8542593348eSBarry Smith   else {
8552593348eSBarry Smith     nz = 0;
856b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
8572593348eSBarry Smith   }
8582593348eSBarry Smith 
8592593348eSBarry Smith   /* allocate the matrix space */
860b6490206SBarry Smith   len   = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int);
8612593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
862b6490206SBarry Smith   PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar));
863b6490206SBarry Smith   b->j  = (int *) (b->a + nz*bs*bs);
8642593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
8652593348eSBarry Smith   b->i  = b->j + nz;
8662593348eSBarry Smith   b->singlemalloc = 1;
8672593348eSBarry Smith 
868de6a44a3SBarry Smith   b->i[0] = 0;
869b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
8702593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
8712593348eSBarry Smith   }
8722593348eSBarry Smith 
873b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
874b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
875b6490206SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
876b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
8772593348eSBarry Smith 
878b6490206SBarry Smith   b->bs               = bs;
879b6490206SBarry Smith   b->mbs              = mbs;
8802593348eSBarry Smith   b->nz               = 0;
8812593348eSBarry Smith   b->maxnz            = nz;
8822593348eSBarry Smith   b->sorted           = 0;
8832593348eSBarry Smith   b->roworiented      = 1;
8842593348eSBarry Smith   b->nonew            = 0;
8852593348eSBarry Smith   b->diag             = 0;
8862593348eSBarry Smith   b->solve_work       = 0;
887de6a44a3SBarry Smith   b->mult_work        = 0;
8882593348eSBarry Smith   b->spptr            = 0;
8892593348eSBarry Smith 
8902593348eSBarry Smith   *A = B;
8912593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
8922593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
8932593348eSBarry Smith   return 0;
8942593348eSBarry Smith }
8952593348eSBarry Smith 
896b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
8972593348eSBarry Smith {
8982593348eSBarry Smith   Mat         C;
899b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
900de6a44a3SBarry Smith   int         i,len, mbs = a->mbs, bs = a->bs,nz = a->nz;
901de6a44a3SBarry Smith 
902de6a44a3SBarry Smith   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
9032593348eSBarry Smith 
9042593348eSBarry Smith   *B = 0;
905b6490206SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
9062593348eSBarry Smith   PLogObjectCreate(C);
907b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
9082593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
909b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
910b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
9112593348eSBarry Smith   C->factor     = A->factor;
9122593348eSBarry Smith   c->row        = 0;
9132593348eSBarry Smith   c->col        = 0;
9142593348eSBarry Smith   C->assembled  = PETSC_TRUE;
9152593348eSBarry Smith 
9162593348eSBarry Smith   c->m          = a->m;
9172593348eSBarry Smith   c->n          = a->n;
918b6490206SBarry Smith   c->bs         = a->bs;
919b6490206SBarry Smith   c->mbs        = a->mbs;
920b6490206SBarry Smith   c->nbs        = a->nbs;
9212593348eSBarry Smith 
922b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
923b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
924b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
9252593348eSBarry Smith     c->imax[i] = a->imax[i];
9262593348eSBarry Smith     c->ilen[i] = a->ilen[i];
9272593348eSBarry Smith   }
9282593348eSBarry Smith 
9292593348eSBarry Smith   /* allocate the matrix space */
9302593348eSBarry Smith   c->singlemalloc = 1;
931de6a44a3SBarry Smith   len   = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int));
9322593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
933de6a44a3SBarry Smith   c->j  = (int *) (c->a + nz*bs*bs);
934de6a44a3SBarry Smith   c->i  = c->j + nz;
935b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
936b6490206SBarry Smith   if (mbs > 0) {
937de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
9382593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
939de6a44a3SBarry Smith       PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar));
9402593348eSBarry Smith     }
9412593348eSBarry Smith   }
9422593348eSBarry Smith 
943b6490206SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
9442593348eSBarry Smith   c->sorted      = a->sorted;
9452593348eSBarry Smith   c->roworiented = a->roworiented;
9462593348eSBarry Smith   c->nonew       = a->nonew;
9472593348eSBarry Smith 
9482593348eSBarry Smith   if (a->diag) {
949b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
950b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
951b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
9522593348eSBarry Smith       c->diag[i] = a->diag[i];
9532593348eSBarry Smith     }
9542593348eSBarry Smith   }
9552593348eSBarry Smith   else c->diag          = 0;
9562593348eSBarry Smith   c->nz                 = a->nz;
9572593348eSBarry Smith   c->maxnz              = a->maxnz;
9582593348eSBarry Smith   c->solve_work         = 0;
9592593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
9607fc0212eSBarry Smith   c->mult_work          = 0;
9612593348eSBarry Smith   *B = C;
9622593348eSBarry Smith   return 0;
9632593348eSBarry Smith }
9642593348eSBarry Smith 
96519bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
9662593348eSBarry Smith {
967b6490206SBarry Smith   Mat_SeqBAIJ  *a;
9682593348eSBarry Smith   Mat          B;
969de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
970b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
97135aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
972de6a44a3SBarry Smith   int          *masked, nmask,tmp,ishift,bs2;
973b6490206SBarry Smith   Scalar       *aa;
97419bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
9752593348eSBarry Smith 
97635aab85fSBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
977de6a44a3SBarry Smith   bs2  = bs*bs;
978b6490206SBarry Smith 
9792593348eSBarry Smith   MPI_Comm_size(comm,&size);
980b6490206SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
98190ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
98277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
983de6a44a3SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
9842593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
9852593348eSBarry Smith 
986b6490206SBarry Smith   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
98735aab85fSBarry Smith 
98835aab85fSBarry Smith   /*
98935aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
99035aab85fSBarry Smith     divisible by the blocksize
99135aab85fSBarry Smith   */
992b6490206SBarry Smith   mbs        = M/bs;
99335aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
99435aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
99535aab85fSBarry Smith   else                  mbs++;
99635aab85fSBarry Smith   if (extra_rows) {
99735aab85fSBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize");
99835aab85fSBarry Smith   }
999b6490206SBarry Smith 
10002593348eSBarry Smith   /* read in row lengths */
100135aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
100277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
100335aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
10042593348eSBarry Smith 
1005b6490206SBarry Smith   /* read in column indices */
100635aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
100777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
100835aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1009b6490206SBarry Smith 
1010b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1011b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1012b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
101335aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
101435aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
101535aab85fSBarry Smith   masked = mask + mbs;
1016b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1017b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
101835aab85fSBarry Smith     nmask = 0;
1019b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1020b6490206SBarry Smith       kmax = rowlengths[rowcount];
1021b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
102235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
102335aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1024b6490206SBarry Smith       }
1025b6490206SBarry Smith       rowcount++;
1026b6490206SBarry Smith     }
102735aab85fSBarry Smith     browlengths[i] += nmask;
102835aab85fSBarry Smith     /* zero out the mask elements we set */
102935aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1030b6490206SBarry Smith   }
1031b6490206SBarry Smith 
10322593348eSBarry Smith   /* create our matrix */
103335aab85fSBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
103435aab85fSBarry Smith          CHKERRQ(ierr);
10352593348eSBarry Smith   B = *A;
1036b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
10372593348eSBarry Smith 
10382593348eSBarry Smith   /* set matrix "i" values */
1039de6a44a3SBarry Smith   a->i[0] = 0;
1040b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
1041b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1042b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
10432593348eSBarry Smith   }
10449867e207SSatish Balay   a->indexshift = 0;
1045b6490206SBarry Smith   a->nz         = 0;
1046b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
10472593348eSBarry Smith 
1048b6490206SBarry Smith   /* read in nonzero values */
104935aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
105077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
105135aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1052b6490206SBarry Smith 
1053b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1054b6490206SBarry Smith   nzcount = 0; jcount = 0;
1055b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
1056b6490206SBarry Smith     nzcountb = nzcount;
105735aab85fSBarry Smith     nmask    = 0;
1058b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1059b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1060b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
106135aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
106235aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1063b6490206SBarry Smith       }
1064b6490206SBarry Smith       rowcount++;
1065b6490206SBarry Smith     }
1066de6a44a3SBarry Smith     /* sort the masked values */
106777c4ece6SBarry Smith     PetscSortInt(nmask,masked);
1068de6a44a3SBarry Smith 
1069b6490206SBarry Smith     /* set "j" values into matrix */
1070b6490206SBarry Smith     maskcount = 1;
107135aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
107235aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1073de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1074b6490206SBarry Smith     }
1075b6490206SBarry Smith     /* set "a" values into matrix */
1076de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1077b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
1078b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1079b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
1080de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
1081de6a44a3SBarry Smith         block  = mask[tmp] - 1;
1082de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
1083de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1084b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1085b6490206SBarry Smith       }
1086b6490206SBarry Smith     }
108735aab85fSBarry Smith     /* zero out the mask elements we set */
108835aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1089b6490206SBarry Smith   }
109035aab85fSBarry Smith   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
1091b6490206SBarry Smith 
1092b6490206SBarry Smith   PetscFree(rowlengths);
1093b6490206SBarry Smith   PetscFree(browlengths);
1094b6490206SBarry Smith   PetscFree(aa);
1095b6490206SBarry Smith   PetscFree(jj);
1096b6490206SBarry Smith   PetscFree(mask);
1097b6490206SBarry Smith 
1098b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1099de6a44a3SBarry Smith 
1100de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
1101de6a44a3SBarry Smith   if (flg) {
110219bcc07fSBarry Smith     Viewer tviewer;
110319bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
110490ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr);
110519bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
110619bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1107de6a44a3SBarry Smith   }
1108de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
1109de6a44a3SBarry Smith   if (flg) {
111019bcc07fSBarry Smith     Viewer tviewer;
111119bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
111290ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
111319bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
111419bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1115de6a44a3SBarry Smith   }
1116de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
1117de6a44a3SBarry Smith   if (flg) {
111819bcc07fSBarry Smith     Viewer tviewer;
111919bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
112019bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
112119bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1122de6a44a3SBarry Smith   }
1123de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
1124de6a44a3SBarry Smith   if (flg) {
112519bcc07fSBarry Smith     Viewer tviewer;
112619bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
112790ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr);
112819bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
112919bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1130de6a44a3SBarry Smith   }
1131de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
1132de6a44a3SBarry Smith   if (flg) {
113319bcc07fSBarry Smith     Viewer tviewer;
113419bcc07fSBarry Smith     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
113519bcc07fSBarry Smith     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
113619bcc07fSBarry Smith     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
113719bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1138de6a44a3SBarry Smith   }
11392593348eSBarry Smith   return 0;
11402593348eSBarry Smith }
11412593348eSBarry Smith 
11422593348eSBarry Smith 
11432593348eSBarry Smith 
1144