xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 9867e2074ee79691714c5c21789531d1a36e00f9)
1b6490206SBarry Smith 
22593348eSBarry Smith #ifndef lint
3*9867e207SSatish Balay static char vcid[] = "$Id: baij.c,v 1.19 1996/03/26 18:56:10 bsmith 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*9867e207SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
225*9867e207SSatish Balay {
226*9867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
227*9867e207SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i;
228*9867e207SSatish Balay   Scalar      *aa,*v_i,*aa_i;
229*9867e207SSatish Balay 
230*9867e207SSatish Balay   bs = a->bs;
231*9867e207SSatish Balay   ai = a->i;
232*9867e207SSatish Balay   aj = a->j;
233*9867e207SSatish Balay   aa = a->a;
234*9867e207SSatish Balay 
235*9867e207SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range");
236*9867e207SSatish Balay 
237*9867e207SSatish Balay   bn  = row/bs;   /* Block number */
238*9867e207SSatish Balay   bp  = row % bs; /* Block Position */
239*9867e207SSatish Balay   M   = ai[bn+1] - ai[bn];
240*9867e207SSatish Balay   *nz = bs*M;
241*9867e207SSatish Balay 
242*9867e207SSatish Balay   if (v) {
243*9867e207SSatish Balay     *v = 0;
244*9867e207SSatish Balay     if (*nz) {
245*9867e207SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
246*9867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
247*9867e207SSatish Balay         v_i  = *v + i*bs;
248*9867e207SSatish Balay         aa_i = aa + bs*bs*(ai[bn] + i);
249*9867e207SSatish Balay         for ( j=bp,k=0; j<bs*bs; j+=bs,k++ ) {v_i[k] = aa_i[j];}
250*9867e207SSatish Balay       }
251*9867e207SSatish Balay     }
252*9867e207SSatish Balay   }
253*9867e207SSatish Balay 
254*9867e207SSatish Balay   if (idx) {
255*9867e207SSatish Balay     *idx = 0;
256*9867e207SSatish Balay     if (*nz) {
257*9867e207SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
258*9867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
259*9867e207SSatish Balay         idx_i = *idx + i*bs;
260*9867e207SSatish Balay         itmp  = bs*aj[ai[bn] + i];
261*9867e207SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
262*9867e207SSatish Balay       }
263*9867e207SSatish Balay     }
264*9867e207SSatish Balay   }
265*9867e207SSatish Balay   return 0;
266*9867e207SSatish Balay }
267*9867e207SSatish Balay 
268*9867e207SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
269*9867e207SSatish Balay {
270*9867e207SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
271*9867e207SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
272*9867e207SSatish Balay   return 0;
273*9867e207SSatish Balay }
274b6490206SBarry Smith 
275b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A)
2762593348eSBarry Smith {
277b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
278de6a44a3SBarry Smith   PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar));
2792593348eSBarry Smith   return 0;
2802593348eSBarry Smith }
2812593348eSBarry Smith 
282b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj)
2832593348eSBarry Smith {
2842593348eSBarry Smith   Mat        A  = (Mat) obj;
285b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2862593348eSBarry Smith 
2872593348eSBarry Smith #if defined(PETSC_LOG)
2882593348eSBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
2892593348eSBarry Smith #endif
2902593348eSBarry Smith   PetscFree(a->a);
2912593348eSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
2922593348eSBarry Smith   if (a->diag) PetscFree(a->diag);
2932593348eSBarry Smith   if (a->ilen) PetscFree(a->ilen);
2942593348eSBarry Smith   if (a->imax) PetscFree(a->imax);
2952593348eSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
296de6a44a3SBarry Smith   if (a->mult_work) PetscFree(a->mult_work);
2972593348eSBarry Smith   PetscFree(a);
2982593348eSBarry Smith   return 0;
2992593348eSBarry Smith }
3002593348eSBarry Smith 
301b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
3022593348eSBarry Smith {
303b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3042593348eSBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
3052593348eSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
3062593348eSBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
3072593348eSBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
3082593348eSBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
3092593348eSBarry Smith   else if (op == ROWS_SORTED ||
3102593348eSBarry Smith            op == SYMMETRIC_MATRIX ||
3112593348eSBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
3122593348eSBarry Smith            op == YES_NEW_DIAGONALS)
31394a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
3142593348eSBarry Smith   else if (op == NO_NEW_DIAGONALS)
315b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
3162593348eSBarry Smith   else
317b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
3182593348eSBarry Smith   return 0;
3192593348eSBarry Smith }
3202593348eSBarry Smith 
3212593348eSBarry Smith 
3222593348eSBarry Smith /* -------------------------------------------------------*/
3232593348eSBarry Smith /* Should check that shapes of vectors and matrices match */
3242593348eSBarry Smith /* -------------------------------------------------------*/
325b6490206SBarry Smith #include "pinclude/plapack.h"
326b6490206SBarry Smith 
327b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy)
3282593348eSBarry Smith {
329b6490206SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
330b6490206SBarry Smith   Scalar          *xg,*yg;
331b6490206SBarry Smith   register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5;
332b6490206SBarry Smith   register Scalar x1,x2,x3,x4,x5;
333b6490206SBarry Smith   int             mbs = a->mbs, m = a->m, i, *idx,*ii;
334b6490206SBarry Smith   int             bs = a->bs,j,n,bs2 = bs*bs;
3352593348eSBarry Smith 
336b6490206SBarry Smith   VecGetArray(xx,&xg); x = xg;  VecGetArray(yy,&yg); y = yg;
337b6490206SBarry Smith   PetscMemzero(y,m*sizeof(Scalar));
338b6490206SBarry Smith   x     = x;
3392593348eSBarry Smith   idx   = a->j;
3402593348eSBarry Smith   v     = a->a;
3412593348eSBarry Smith   ii    = a->i;
342b6490206SBarry Smith 
343b6490206SBarry Smith   switch (bs) {
344b6490206SBarry Smith     case 1:
3452593348eSBarry Smith      for ( i=0; i<m; i++ ) {
3462593348eSBarry Smith        n    = ii[1] - ii[0]; ii++;
3472593348eSBarry Smith        sum  = 0.0;
3482593348eSBarry Smith        while (n--) sum += *v++ * x[*idx++];
3492593348eSBarry Smith        y[i] = sum;
3502593348eSBarry Smith       }
3512593348eSBarry Smith       break;
352b6490206SBarry Smith     case 2:
353b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
354de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
355b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0;
356b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
357b6490206SBarry Smith           xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
358b6490206SBarry Smith           sum1 += v[0]*x1 + v[2]*x2;
359b6490206SBarry Smith           sum2 += v[1]*x1 + v[3]*x2;
360b6490206SBarry Smith           v += 4;
361b6490206SBarry Smith         }
362b6490206SBarry Smith         y[0] += sum1; y[1] += sum2;
363b6490206SBarry Smith         y += 2;
364b6490206SBarry Smith       }
365b6490206SBarry Smith       break;
366b6490206SBarry Smith     case 3:
367b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
368de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
369b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
370b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
371b6490206SBarry Smith           xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
372b6490206SBarry Smith           sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
373b6490206SBarry Smith           sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
374b6490206SBarry Smith           sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
375b6490206SBarry Smith           v += 9;
376b6490206SBarry Smith         }
377b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3;
378b6490206SBarry Smith         y += 3;
379b6490206SBarry Smith       }
380b6490206SBarry Smith       break;
381b6490206SBarry Smith     case 4:
382b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
383de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
384b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
385b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
386b6490206SBarry Smith           xb = x + 4*(*idx++);
387b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
388b6490206SBarry Smith           sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
389b6490206SBarry Smith           sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
390b6490206SBarry Smith           sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
391b6490206SBarry Smith           sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
392b6490206SBarry Smith           v += 16;
393b6490206SBarry Smith         }
394b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4;
395b6490206SBarry Smith         y += 4;
396b6490206SBarry Smith       }
397b6490206SBarry Smith       break;
398b6490206SBarry Smith     case 5:
399b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
400de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
401b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
402b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
403b6490206SBarry Smith           xb = x + 5*(*idx++);
404b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
405b6490206SBarry Smith           sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
406b6490206SBarry Smith           sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
407b6490206SBarry Smith           sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
408b6490206SBarry Smith           sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
409b6490206SBarry Smith           sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
410b6490206SBarry Smith           v += 25;
411b6490206SBarry Smith         }
412b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5;
413b6490206SBarry Smith         y += 5;
414b6490206SBarry Smith       }
415b6490206SBarry Smith       break;
416b6490206SBarry Smith       /* block sizes larger then 5 by 5 are handled by BLAS */
417b6490206SBarry Smith     default: {
418de6a44a3SBarry Smith       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
419de6a44a3SBarry Smith       if (!a->mult_work) {
420de6a44a3SBarry Smith         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
421de6a44a3SBarry Smith         CHKPTRQ(a->mult_work);
422de6a44a3SBarry Smith       }
423de6a44a3SBarry Smith       work = a->mult_work;
424b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
425de6a44a3SBarry Smith         n     = ii[1] - ii[0]; ii++;
426de6a44a3SBarry Smith         ncols = n*bs;
427de6a44a3SBarry Smith         workt = work;
428b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
429b6490206SBarry Smith           xb = x + bs*(*idx++);
430de6a44a3SBarry Smith           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
431de6a44a3SBarry Smith           workt += bs;
432b6490206SBarry Smith         }
433de6a44a3SBarry Smith         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One);
434de6a44a3SBarry Smith         v += n*bs2;
435b6490206SBarry Smith         y += bs;
4362593348eSBarry Smith       }
4372593348eSBarry Smith     }
4382593348eSBarry Smith   }
439b6490206SBarry Smith   PLogFlops(2*bs2*a->nz - m);
4402593348eSBarry Smith   return 0;
4412593348eSBarry Smith }
4422593348eSBarry Smith 
443de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
4442593348eSBarry Smith {
445b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
446bcd2baecSBarry Smith   if (nz)  *nz  = a->bs*a->bs*a->nz;
447bcd2baecSBarry Smith   if (nza) *nza = a->maxnz;
448bcd2baecSBarry Smith   if (mem) *mem = (int)A->mem;
4492593348eSBarry Smith   return 0;
4502593348eSBarry Smith }
4512593348eSBarry Smith 
45291d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
45391d316f6SSatish Balay {
45491d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
45591d316f6SSatish Balay 
45691d316f6SSatish Balay   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
45791d316f6SSatish Balay 
45891d316f6SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
45991d316f6SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
46091d316f6SSatish Balay       (a->nz != b->nz)||(a->indexshift != b->indexshift)) {
46191d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
46291d316f6SSatish Balay   }
46391d316f6SSatish Balay 
46491d316f6SSatish Balay   /* if the a->i are the same */
46591d316f6SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
46691d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
46791d316f6SSatish Balay   }
46891d316f6SSatish Balay 
46991d316f6SSatish Balay   /* if a->j are the same */
47091d316f6SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
47191d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
47291d316f6SSatish Balay   }
47391d316f6SSatish Balay 
47491d316f6SSatish Balay   /* if a->a are the same */
47591d316f6SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
47691d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
47791d316f6SSatish Balay   }
47891d316f6SSatish Balay   *flg = PETSC_TRUE;
47991d316f6SSatish Balay   return 0;
48091d316f6SSatish Balay 
48191d316f6SSatish Balay }
48291d316f6SSatish Balay 
48391d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
48491d316f6SSatish Balay {
48591d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
48617e48fc4SSatish Balay   int         i,j,k,n,row,bs,*ai,*aj,ambs;
48717e48fc4SSatish Balay   Scalar      *x, zero = 0.0,*aa,*aa_j;
48817e48fc4SSatish Balay 
48917e48fc4SSatish Balay   bs  = a->bs;
49017e48fc4SSatish Balay   aa   = a->a;
49117e48fc4SSatish Balay   ai   = a->i;
49217e48fc4SSatish Balay   aj   = a->j;
49317e48fc4SSatish Balay   ambs = a->mbs;
49491d316f6SSatish Balay 
49591d316f6SSatish Balay   VecSet(&zero,v);
49691d316f6SSatish Balay   VecGetArray(v,&x); VecGetLocalSize(v,&n);
497*9867e207SSatish Balay   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector");
49817e48fc4SSatish Balay   for ( i=0; i<ambs; i++ ) {
49917e48fc4SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
50017e48fc4SSatish Balay       if (aj[j] == i) {
50117e48fc4SSatish Balay         row  = i*bs;
50217e48fc4SSatish Balay         aa_j = aa+j*bs*bs;
50317e48fc4SSatish Balay         for (k=0; k<bs*bs; k+=(bs+1),row++) x[row] = aa_j[k];
50491d316f6SSatish Balay         break;
50591d316f6SSatish Balay       }
50691d316f6SSatish Balay     }
50791d316f6SSatish Balay   }
50891d316f6SSatish Balay   return 0;
50991d316f6SSatish Balay }
51091d316f6SSatish Balay 
511*9867e207SSatish Balay static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
512*9867e207SSatish Balay {
513*9867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
514*9867e207SSatish Balay   Scalar      *l,*r,x,*v,*aa,*li,*ri;
515*9867e207SSatish Balay   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs;
516*9867e207SSatish Balay 
517*9867e207SSatish Balay   ai  = a->i;
518*9867e207SSatish Balay   aj  = a->j;
519*9867e207SSatish Balay   aa  = a->a;
520*9867e207SSatish Balay   m   = a->m;
521*9867e207SSatish Balay   n   = a->n;
522*9867e207SSatish Balay   bs  = a->bs;
523*9867e207SSatish Balay   mbs = a->mbs;
524*9867e207SSatish Balay 
525*9867e207SSatish Balay   if (ll) {
526*9867e207SSatish Balay     VecGetArray(ll,&l); VecGetSize(ll,&lm);
527*9867e207SSatish Balay     if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length");
528*9867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
529*9867e207SSatish Balay       M  = ai[i+1] - ai[i];
530*9867e207SSatish Balay       li = l + i*bs;
531*9867e207SSatish Balay       v  = aa + bs*bs*ai[i];
532*9867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
533*9867e207SSatish Balay         for ( k=0; k<bs*bs; k++ ) {
534*9867e207SSatish Balay           (*v++) *= li[k%bs];
535*9867e207SSatish Balay         }
536*9867e207SSatish Balay       }
537*9867e207SSatish Balay     }
538*9867e207SSatish Balay   }
539*9867e207SSatish Balay 
540*9867e207SSatish Balay   if (rr) {
541*9867e207SSatish Balay     VecGetArray(rr,&r); VecGetSize(rr,&rn);
542*9867e207SSatish Balay     if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length");
543*9867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
544*9867e207SSatish Balay       M  = ai[i+1] - ai[i];
545*9867e207SSatish Balay       v  = aa + bs*bs*ai[i];
546*9867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
547*9867e207SSatish Balay         ri = r + bs*aj[ai[i]+j];
548*9867e207SSatish Balay         for ( k=0; k<bs; k++ ) {
549*9867e207SSatish Balay           x = ri[k];
550*9867e207SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
551*9867e207SSatish Balay         }
552*9867e207SSatish Balay       }
553*9867e207SSatish Balay     }
554*9867e207SSatish Balay   }
555*9867e207SSatish Balay   return 0;
556*9867e207SSatish Balay }
557*9867e207SSatish Balay 
558*9867e207SSatish Balay 
559b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
560b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
561b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec);
5627fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
5637fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
5647fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
5657fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
5667fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
5677fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
5682593348eSBarry Smith 
569b6490206SBarry Smith static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
5702593348eSBarry Smith {
571b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5722593348eSBarry Smith   *m = a->m; *n = a->n;
5732593348eSBarry Smith   return 0;
5742593348eSBarry Smith }
5752593348eSBarry Smith 
576b6490206SBarry Smith static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
5772593348eSBarry Smith {
578b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5792593348eSBarry Smith   *m = 0; *n = a->m;
5802593348eSBarry Smith   return 0;
5812593348eSBarry Smith }
582b6490206SBarry Smith 
583b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
5842593348eSBarry Smith {
585b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5862593348eSBarry Smith   Scalar      *v = a->a;
5872593348eSBarry Smith   double      sum = 0.0;
588b6490206SBarry Smith   int         i;
5892593348eSBarry Smith 
5902593348eSBarry Smith   if (type == NORM_FROBENIUS) {
5912593348eSBarry Smith     for (i=0; i<a->nz; i++ ) {
5922593348eSBarry Smith #if defined(PETSC_COMPLEX)
5932593348eSBarry Smith       sum += real(conj(*v)*(*v)); v++;
5942593348eSBarry Smith #else
5952593348eSBarry Smith       sum += (*v)*(*v); v++;
5962593348eSBarry Smith #endif
5972593348eSBarry Smith     }
5982593348eSBarry Smith     *norm = sqrt(sum);
5992593348eSBarry Smith   }
6002593348eSBarry Smith   else {
601b6490206SBarry Smith     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
6022593348eSBarry Smith   }
6032593348eSBarry Smith   return 0;
6042593348eSBarry Smith }
6052593348eSBarry Smith 
6062593348eSBarry Smith /*
6072593348eSBarry Smith      note: This can only work for identity for row and col. It would
6082593348eSBarry Smith    be good to check this and otherwise generate an error.
6092593348eSBarry Smith */
610b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
6112593348eSBarry Smith {
612b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
6132593348eSBarry Smith   Mat         outA;
614de6a44a3SBarry Smith   int         ierr;
6152593348eSBarry Smith 
616b6490206SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
6172593348eSBarry Smith 
6182593348eSBarry Smith   outA          = inA;
6192593348eSBarry Smith   inA->factor   = FACTOR_LU;
6202593348eSBarry Smith   a->row        = row;
6212593348eSBarry Smith   a->col        = col;
6222593348eSBarry Smith 
623de6a44a3SBarry Smith   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
6242593348eSBarry Smith 
6252593348eSBarry Smith   if (!a->diag) {
626b6490206SBarry Smith     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
6272593348eSBarry Smith   }
6287fc0212eSBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
6292593348eSBarry Smith   return 0;
6302593348eSBarry Smith }
6312593348eSBarry Smith 
632b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
6332593348eSBarry Smith {
634b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
635b6490206SBarry Smith   int         one = 1, totalnz = a->bs*a->bs*a->nz;
636b6490206SBarry Smith   BLscal_( &totalnz, alpha, a->a, &one );
637b6490206SBarry Smith   PLogFlops(totalnz);
6382593348eSBarry Smith   return 0;
6392593348eSBarry Smith }
6402593348eSBarry Smith 
641b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A)
6422593348eSBarry Smith {
6432593348eSBarry Smith   static int called = 0;
6442593348eSBarry Smith 
6452593348eSBarry Smith   if (called) return 0; else called = 1;
6462593348eSBarry Smith   return 0;
6472593348eSBarry Smith }
6482593348eSBarry Smith /* -------------------------------------------------------------------*/
649b6490206SBarry Smith static struct _MatOps MatOps = {0,
650*9867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
651b6490206SBarry Smith        MatMult_SeqBAIJ,0,
652b6490206SBarry Smith        0,0,
653de6a44a3SBarry Smith        MatSolve_SeqBAIJ,0,
654b6490206SBarry Smith        0,0,
655de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
656b6490206SBarry Smith        0,
657b6490206SBarry Smith        0,
65817e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
659*9867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
660b6490206SBarry Smith        0,0,
661b6490206SBarry Smith        0,
662b6490206SBarry Smith        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
663b6490206SBarry Smith        MatGetReordering_SeqBAIJ,
6647fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
665b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
666de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
667b6490206SBarry Smith        0,0,/*  MatConvert_SeqBAIJ  */ 0,
668b6490206SBarry Smith        0,0,
669b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
670b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
671b6490206SBarry Smith        0,0,
672b6490206SBarry Smith        0,0,
673b6490206SBarry Smith        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
674b6490206SBarry Smith        0};
6752593348eSBarry Smith 
6762593348eSBarry Smith /*@C
677b6490206SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format
6782593348eSBarry Smith    (the default parallel PETSc format).  For good matrix assembly performance
6792593348eSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
6802593348eSBarry Smith    (or nzz).  By setting these parameters accurately, performance can be
6812593348eSBarry Smith    increased by more than a factor of 50.
6822593348eSBarry Smith 
6832593348eSBarry Smith    Input Parameters:
6842593348eSBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
685b6490206SBarry Smith .  bs - size of block
6862593348eSBarry Smith .  m - number of rows
6872593348eSBarry Smith .  n - number of columns
688b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
689b6490206SBarry Smith .  nzz - number of block nonzeros per block row or PETSC_NULL
690b6490206SBarry Smith          (possibly different for each row)
6912593348eSBarry Smith 
6922593348eSBarry Smith    Output Parameter:
6932593348eSBarry Smith .  A - the matrix
6942593348eSBarry Smith 
6952593348eSBarry Smith    Notes:
696b6490206SBarry Smith    The BAIJ format, is fully compatible with standard Fortran 77
6972593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
6982593348eSBarry Smith    either one (as in Fortran) or zero.  See the users manual for details.
6992593348eSBarry Smith 
7002593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
7012593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
7022593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
7032593348eSBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
7042593348eSBarry Smith 
7052593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
7062593348eSBarry Smith @*/
707b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
7082593348eSBarry Smith {
7092593348eSBarry Smith   Mat         B;
710b6490206SBarry Smith   Mat_SeqBAIJ *b;
711b6490206SBarry Smith   int         i,len,ierr, flg,mbs = m/bs;
712b6490206SBarry Smith 
713b6490206SBarry Smith   if (mbs*bs != m)
714b6490206SBarry Smith     SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize");
7152593348eSBarry Smith 
7162593348eSBarry Smith   *A      = 0;
717b6490206SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
7182593348eSBarry Smith   PLogObjectCreate(B);
719b6490206SBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
7202593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
7217fc0212eSBarry Smith   switch (bs) {
7227fc0212eSBarry Smith     case 1:
7237fc0212eSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
7247fc0212eSBarry Smith       break;
7254eeb42bcSBarry Smith     case 2:
7264eeb42bcSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
7276d84be18SBarry Smith       break;
728f327dd97SBarry Smith     case 3:
729f327dd97SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
7304eeb42bcSBarry Smith       break;
7316d84be18SBarry Smith     case 4:
7326d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
7336d84be18SBarry Smith       break;
7346d84be18SBarry Smith     case 5:
7356d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
7366d84be18SBarry Smith       break;
7377fc0212eSBarry Smith   }
738b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
739b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
7402593348eSBarry Smith   B->factor           = 0;
7412593348eSBarry Smith   B->lupivotthreshold = 1.0;
7422593348eSBarry Smith   b->row              = 0;
7432593348eSBarry Smith   b->col              = 0;
7442593348eSBarry Smith   b->reallocs         = 0;
7452593348eSBarry Smith 
7462593348eSBarry Smith   b->m       = m;
747b6490206SBarry Smith   b->mbs     = mbs;
7482593348eSBarry Smith   b->n       = n;
749b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
7502593348eSBarry Smith   if (nnz == PETSC_NULL) {
751de6a44a3SBarry Smith     if (nz == PETSC_DEFAULT) nz = 5;
7522593348eSBarry Smith     else if (nz <= 0)        nz = 1;
753b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
754b6490206SBarry Smith     nz = nz*mbs;
7552593348eSBarry Smith   }
7562593348eSBarry Smith   else {
7572593348eSBarry Smith     nz = 0;
758b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
7592593348eSBarry Smith   }
7602593348eSBarry Smith 
7612593348eSBarry Smith   /* allocate the matrix space */
762b6490206SBarry Smith   len   = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int);
7632593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
764b6490206SBarry Smith   PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar));
765b6490206SBarry Smith   b->j  = (int *) (b->a + nz*bs*bs);
7662593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
7672593348eSBarry Smith   b->i  = b->j + nz;
7682593348eSBarry Smith   b->singlemalloc = 1;
7692593348eSBarry Smith 
770de6a44a3SBarry Smith   b->i[0] = 0;
771b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
7722593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
7732593348eSBarry Smith   }
7742593348eSBarry Smith 
775b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
776b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
777b6490206SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
778b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
7792593348eSBarry Smith 
780b6490206SBarry Smith   b->bs               = bs;
781b6490206SBarry Smith   b->mbs              = mbs;
7822593348eSBarry Smith   b->nz               = 0;
7832593348eSBarry Smith   b->maxnz            = nz;
7842593348eSBarry Smith   b->sorted           = 0;
7852593348eSBarry Smith   b->roworiented      = 1;
7862593348eSBarry Smith   b->nonew            = 0;
7872593348eSBarry Smith   b->diag             = 0;
7882593348eSBarry Smith   b->solve_work       = 0;
789de6a44a3SBarry Smith   b->mult_work        = 0;
7902593348eSBarry Smith   b->spptr            = 0;
7912593348eSBarry Smith 
7922593348eSBarry Smith   *A = B;
7932593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
7942593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
7952593348eSBarry Smith   return 0;
7962593348eSBarry Smith }
7972593348eSBarry Smith 
798b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
7992593348eSBarry Smith {
8002593348eSBarry Smith   Mat         C;
801b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
802de6a44a3SBarry Smith   int         i,len, mbs = a->mbs, bs = a->bs,nz = a->nz;
803de6a44a3SBarry Smith 
804de6a44a3SBarry Smith   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
8052593348eSBarry Smith 
8062593348eSBarry Smith   *B = 0;
807b6490206SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
8082593348eSBarry Smith   PLogObjectCreate(C);
809b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
8102593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
811b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
812b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
8132593348eSBarry Smith   C->factor     = A->factor;
8142593348eSBarry Smith   c->row        = 0;
8152593348eSBarry Smith   c->col        = 0;
8162593348eSBarry Smith   C->assembled  = PETSC_TRUE;
8172593348eSBarry Smith 
8182593348eSBarry Smith   c->m          = a->m;
8192593348eSBarry Smith   c->n          = a->n;
820b6490206SBarry Smith   c->bs         = a->bs;
821b6490206SBarry Smith   c->mbs        = a->mbs;
822b6490206SBarry Smith   c->nbs        = a->nbs;
8232593348eSBarry Smith 
824b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
825b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
826b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
8272593348eSBarry Smith     c->imax[i] = a->imax[i];
8282593348eSBarry Smith     c->ilen[i] = a->ilen[i];
8292593348eSBarry Smith   }
8302593348eSBarry Smith 
8312593348eSBarry Smith   /* allocate the matrix space */
8322593348eSBarry Smith   c->singlemalloc = 1;
833de6a44a3SBarry Smith   len   = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int));
8342593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
835de6a44a3SBarry Smith   c->j  = (int *) (c->a + nz*bs*bs);
836de6a44a3SBarry Smith   c->i  = c->j + nz;
837b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
838b6490206SBarry Smith   if (mbs > 0) {
839de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
8402593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
841de6a44a3SBarry Smith       PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar));
8422593348eSBarry Smith     }
8432593348eSBarry Smith   }
8442593348eSBarry Smith 
845b6490206SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
8462593348eSBarry Smith   c->sorted      = a->sorted;
8472593348eSBarry Smith   c->roworiented = a->roworiented;
8482593348eSBarry Smith   c->nonew       = a->nonew;
8492593348eSBarry Smith 
8502593348eSBarry Smith   if (a->diag) {
851b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
852b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
853b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
8542593348eSBarry Smith       c->diag[i] = a->diag[i];
8552593348eSBarry Smith     }
8562593348eSBarry Smith   }
8572593348eSBarry Smith   else c->diag          = 0;
8582593348eSBarry Smith   c->nz                 = a->nz;
8592593348eSBarry Smith   c->maxnz              = a->maxnz;
8602593348eSBarry Smith   c->solve_work         = 0;
8612593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
8627fc0212eSBarry Smith   c->mult_work          = 0;
8632593348eSBarry Smith   *B = C;
8642593348eSBarry Smith   return 0;
8652593348eSBarry Smith }
8662593348eSBarry Smith 
86719bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
8682593348eSBarry Smith {
869b6490206SBarry Smith   Mat_SeqBAIJ  *a;
8702593348eSBarry Smith   Mat          B;
871de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
872b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
87335aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
874de6a44a3SBarry Smith   int          *masked, nmask,tmp,ishift,bs2;
875b6490206SBarry Smith   Scalar       *aa;
87619bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
8772593348eSBarry Smith 
87835aab85fSBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
879de6a44a3SBarry Smith   bs2  = bs*bs;
880b6490206SBarry Smith 
8812593348eSBarry Smith   MPI_Comm_size(comm,&size);
882b6490206SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
88390ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
88477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
885de6a44a3SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
8862593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
8872593348eSBarry Smith 
888b6490206SBarry Smith   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
88935aab85fSBarry Smith 
89035aab85fSBarry Smith   /*
89135aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
89235aab85fSBarry Smith     divisible by the blocksize
89335aab85fSBarry Smith   */
894b6490206SBarry Smith   mbs        = M/bs;
89535aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
89635aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
89735aab85fSBarry Smith   else                  mbs++;
89835aab85fSBarry Smith   if (extra_rows) {
89935aab85fSBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize");
90035aab85fSBarry Smith   }
901b6490206SBarry Smith 
9022593348eSBarry Smith   /* read in row lengths */
90335aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
90477c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
90535aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
9062593348eSBarry Smith 
907b6490206SBarry Smith   /* read in column indices */
90835aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
90977c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
91035aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
911b6490206SBarry Smith 
912b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
913b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
914b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
91535aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
91635aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
91735aab85fSBarry Smith   masked = mask + mbs;
918b6490206SBarry Smith   rowcount = 0; nzcount = 0;
919b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
92035aab85fSBarry Smith     nmask = 0;
921b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
922b6490206SBarry Smith       kmax = rowlengths[rowcount];
923b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
92435aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
92535aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
926b6490206SBarry Smith       }
927b6490206SBarry Smith       rowcount++;
928b6490206SBarry Smith     }
92935aab85fSBarry Smith     browlengths[i] += nmask;
93035aab85fSBarry Smith     /* zero out the mask elements we set */
93135aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
932b6490206SBarry Smith   }
933b6490206SBarry Smith 
9342593348eSBarry Smith   /* create our matrix */
93535aab85fSBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
93635aab85fSBarry Smith          CHKERRQ(ierr);
9372593348eSBarry Smith   B = *A;
938b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
9392593348eSBarry Smith 
9402593348eSBarry Smith   /* set matrix "i" values */
941de6a44a3SBarry Smith   a->i[0] = 0;
942b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
943b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
944b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
9452593348eSBarry Smith   }
946*9867e207SSatish Balay   a->indexshift = 0;
947b6490206SBarry Smith   a->nz         = 0;
948b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
9492593348eSBarry Smith 
950b6490206SBarry Smith   /* read in nonzero values */
95135aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
95277c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
95335aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
954b6490206SBarry Smith 
955b6490206SBarry Smith   /* set "a" and "j" values into matrix */
956b6490206SBarry Smith   nzcount = 0; jcount = 0;
957b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
958b6490206SBarry Smith     nzcountb = nzcount;
95935aab85fSBarry Smith     nmask    = 0;
960b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
961b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
962b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
96335aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
96435aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
965b6490206SBarry Smith       }
966b6490206SBarry Smith       rowcount++;
967b6490206SBarry Smith     }
968de6a44a3SBarry Smith     /* sort the masked values */
96977c4ece6SBarry Smith     PetscSortInt(nmask,masked);
970de6a44a3SBarry Smith 
971b6490206SBarry Smith     /* set "j" values into matrix */
972b6490206SBarry Smith     maskcount = 1;
97335aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
97435aab85fSBarry Smith       a->j[jcount++]  = masked[j];
975de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
976b6490206SBarry Smith     }
977b6490206SBarry Smith     /* set "a" values into matrix */
978de6a44a3SBarry Smith     ishift = bs2*a->i[i];
979b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
980b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
981b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
982de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
983de6a44a3SBarry Smith         block  = mask[tmp] - 1;
984de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
985de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
986b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
987b6490206SBarry Smith       }
988b6490206SBarry Smith     }
98935aab85fSBarry Smith     /* zero out the mask elements we set */
99035aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
991b6490206SBarry Smith   }
99235aab85fSBarry Smith   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
993b6490206SBarry Smith 
994b6490206SBarry Smith   PetscFree(rowlengths);
995b6490206SBarry Smith   PetscFree(browlengths);
996b6490206SBarry Smith   PetscFree(aa);
997b6490206SBarry Smith   PetscFree(jj);
998b6490206SBarry Smith   PetscFree(mask);
999b6490206SBarry Smith 
1000b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1001de6a44a3SBarry Smith 
1002de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
1003de6a44a3SBarry Smith   if (flg) {
100419bcc07fSBarry Smith     Viewer tviewer;
100519bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
100690ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr);
100719bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
100819bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1009de6a44a3SBarry Smith   }
1010de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
1011de6a44a3SBarry Smith   if (flg) {
101219bcc07fSBarry Smith     Viewer tviewer;
101319bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
101490ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
101519bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
101619bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1017de6a44a3SBarry Smith   }
1018de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
1019de6a44a3SBarry Smith   if (flg) {
102019bcc07fSBarry Smith     Viewer tviewer;
102119bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
102219bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
102319bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1024de6a44a3SBarry Smith   }
1025de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
1026de6a44a3SBarry Smith   if (flg) {
102719bcc07fSBarry Smith     Viewer tviewer;
102819bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
102990ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr);
103019bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
103119bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1032de6a44a3SBarry Smith   }
1033de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
1034de6a44a3SBarry Smith   if (flg) {
103519bcc07fSBarry Smith     Viewer tviewer;
103619bcc07fSBarry Smith     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
103719bcc07fSBarry Smith     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
103819bcc07fSBarry Smith     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
103919bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1040de6a44a3SBarry Smith   }
10412593348eSBarry Smith   return 0;
10422593348eSBarry Smith }
10432593348eSBarry Smith 
10442593348eSBarry Smith 
10452593348eSBarry Smith 
1046