xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 0b824a480752e09606a2525fc849a9d89d7ee6f2)
1b6490206SBarry Smith 
22593348eSBarry Smith #ifndef lint
3*0b824a48SSatish Balay static char vcid[] = "$Id: baij.c,v 1.20 1996/03/27 00:06:28 balay Exp balay $";
42593348eSBarry Smith #endif
52593348eSBarry Smith 
62593348eSBarry Smith /*
7b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
82593348eSBarry Smith   matrix storage format.
92593348eSBarry Smith */
10b6490206SBarry Smith #include "baij.h"
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*0b824a48SSatish Balay static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
225*0b824a48SSatish Balay {
226*0b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
227*0b824a48SSatish Balay   *m = a->m; *n = a->n;
228*0b824a48SSatish Balay   return 0;
229*0b824a48SSatish Balay }
230*0b824a48SSatish Balay 
231*0b824a48SSatish Balay static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
232*0b824a48SSatish Balay {
233*0b824a48SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
234*0b824a48SSatish Balay   *m = 0; *n = a->m;
235*0b824a48SSatish Balay   return 0;
236*0b824a48SSatish Balay }
237*0b824a48SSatish Balay 
2389867e207SSatish Balay int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2399867e207SSatish Balay {
2409867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2419867e207SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i;
2429867e207SSatish Balay   Scalar      *aa,*v_i,*aa_i;
2439867e207SSatish Balay 
2449867e207SSatish Balay   bs = a->bs;
2459867e207SSatish Balay   ai = a->i;
2469867e207SSatish Balay   aj = a->j;
2479867e207SSatish Balay   aa = a->a;
2489867e207SSatish Balay 
2499867e207SSatish Balay   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range");
2509867e207SSatish Balay 
2519867e207SSatish Balay   bn  = row/bs;   /* Block number */
2529867e207SSatish Balay   bp  = row % bs; /* Block Position */
2539867e207SSatish Balay   M   = ai[bn+1] - ai[bn];
2549867e207SSatish Balay   *nz = bs*M;
2559867e207SSatish Balay 
2569867e207SSatish Balay   if (v) {
2579867e207SSatish Balay     *v = 0;
2589867e207SSatish Balay     if (*nz) {
2599867e207SSatish Balay       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
2609867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2619867e207SSatish Balay         v_i  = *v + i*bs;
2629867e207SSatish Balay         aa_i = aa + bs*bs*(ai[bn] + i);
2639867e207SSatish Balay         for ( j=bp,k=0; j<bs*bs; j+=bs,k++ ) {v_i[k] = aa_i[j];}
2649867e207SSatish Balay       }
2659867e207SSatish Balay     }
2669867e207SSatish Balay   }
2679867e207SSatish Balay 
2689867e207SSatish Balay   if (idx) {
2699867e207SSatish Balay     *idx = 0;
2709867e207SSatish Balay     if (*nz) {
2719867e207SSatish Balay       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
2729867e207SSatish Balay       for ( i=0; i<M; i++ ) { /* for each block in the block row */
2739867e207SSatish Balay         idx_i = *idx + i*bs;
2749867e207SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2759867e207SSatish Balay         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
2769867e207SSatish Balay       }
2779867e207SSatish Balay     }
2789867e207SSatish Balay   }
2799867e207SSatish Balay   return 0;
2809867e207SSatish Balay }
2819867e207SSatish Balay 
2829867e207SSatish Balay int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
2839867e207SSatish Balay {
2849867e207SSatish Balay   if (idx) {if (*idx) PetscFree(*idx);}
2859867e207SSatish Balay   if (v)   {if (*v)   PetscFree(*v);}
2869867e207SSatish Balay   return 0;
2879867e207SSatish Balay }
288b6490206SBarry Smith 
289b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A)
2902593348eSBarry Smith {
291b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
292de6a44a3SBarry Smith   PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar));
2932593348eSBarry Smith   return 0;
2942593348eSBarry Smith }
2952593348eSBarry Smith 
296b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj)
2972593348eSBarry Smith {
2982593348eSBarry Smith   Mat        A  = (Mat) obj;
299b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3002593348eSBarry Smith 
3012593348eSBarry Smith #if defined(PETSC_LOG)
3022593348eSBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
3032593348eSBarry Smith #endif
3042593348eSBarry Smith   PetscFree(a->a);
3052593348eSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
3062593348eSBarry Smith   if (a->diag) PetscFree(a->diag);
3072593348eSBarry Smith   if (a->ilen) PetscFree(a->ilen);
3082593348eSBarry Smith   if (a->imax) PetscFree(a->imax);
3092593348eSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
310de6a44a3SBarry Smith   if (a->mult_work) PetscFree(a->mult_work);
3112593348eSBarry Smith   PetscFree(a);
3122593348eSBarry Smith   return 0;
3132593348eSBarry Smith }
3142593348eSBarry Smith 
315b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
3162593348eSBarry Smith {
317b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
3182593348eSBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
3192593348eSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
3202593348eSBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
3212593348eSBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
3222593348eSBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
3232593348eSBarry Smith   else if (op == ROWS_SORTED ||
3242593348eSBarry Smith            op == SYMMETRIC_MATRIX ||
3252593348eSBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
3262593348eSBarry Smith            op == YES_NEW_DIAGONALS)
32794a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
3282593348eSBarry Smith   else if (op == NO_NEW_DIAGONALS)
329b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
3302593348eSBarry Smith   else
331b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
3322593348eSBarry Smith   return 0;
3332593348eSBarry Smith }
3342593348eSBarry Smith 
3352593348eSBarry Smith 
3362593348eSBarry Smith /* -------------------------------------------------------*/
3372593348eSBarry Smith /* Should check that shapes of vectors and matrices match */
3382593348eSBarry Smith /* -------------------------------------------------------*/
339b6490206SBarry Smith #include "pinclude/plapack.h"
340b6490206SBarry Smith 
341b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy)
3422593348eSBarry Smith {
343b6490206SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
344b6490206SBarry Smith   Scalar          *xg,*yg;
345b6490206SBarry Smith   register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5;
346b6490206SBarry Smith   register Scalar x1,x2,x3,x4,x5;
347b6490206SBarry Smith   int             mbs = a->mbs, m = a->m, i, *idx,*ii;
348b6490206SBarry Smith   int             bs = a->bs,j,n,bs2 = bs*bs;
3492593348eSBarry Smith 
350b6490206SBarry Smith   VecGetArray(xx,&xg); x = xg;  VecGetArray(yy,&yg); y = yg;
351b6490206SBarry Smith   PetscMemzero(y,m*sizeof(Scalar));
352b6490206SBarry Smith   x     = x;
3532593348eSBarry Smith   idx   = a->j;
3542593348eSBarry Smith   v     = a->a;
3552593348eSBarry Smith   ii    = a->i;
356b6490206SBarry Smith 
357b6490206SBarry Smith   switch (bs) {
358b6490206SBarry Smith     case 1:
3592593348eSBarry Smith      for ( i=0; i<m; i++ ) {
3602593348eSBarry Smith        n    = ii[1] - ii[0]; ii++;
3612593348eSBarry Smith        sum  = 0.0;
3622593348eSBarry Smith        while (n--) sum += *v++ * x[*idx++];
3632593348eSBarry Smith        y[i] = sum;
3642593348eSBarry Smith       }
3652593348eSBarry Smith       break;
366b6490206SBarry Smith     case 2:
367b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
368de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
369b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0;
370b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
371b6490206SBarry Smith           xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
372b6490206SBarry Smith           sum1 += v[0]*x1 + v[2]*x2;
373b6490206SBarry Smith           sum2 += v[1]*x1 + v[3]*x2;
374b6490206SBarry Smith           v += 4;
375b6490206SBarry Smith         }
376b6490206SBarry Smith         y[0] += sum1; y[1] += sum2;
377b6490206SBarry Smith         y += 2;
378b6490206SBarry Smith       }
379b6490206SBarry Smith       break;
380b6490206SBarry Smith     case 3:
381b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
382de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
383b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
384b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
385b6490206SBarry Smith           xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
386b6490206SBarry Smith           sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
387b6490206SBarry Smith           sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
388b6490206SBarry Smith           sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
389b6490206SBarry Smith           v += 9;
390b6490206SBarry Smith         }
391b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3;
392b6490206SBarry Smith         y += 3;
393b6490206SBarry Smith       }
394b6490206SBarry Smith       break;
395b6490206SBarry Smith     case 4:
396b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
397de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
398b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
399b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
400b6490206SBarry Smith           xb = x + 4*(*idx++);
401b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
402b6490206SBarry Smith           sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
403b6490206SBarry Smith           sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
404b6490206SBarry Smith           sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
405b6490206SBarry Smith           sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
406b6490206SBarry Smith           v += 16;
407b6490206SBarry Smith         }
408b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4;
409b6490206SBarry Smith         y += 4;
410b6490206SBarry Smith       }
411b6490206SBarry Smith       break;
412b6490206SBarry Smith     case 5:
413b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
414de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
415b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
416b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
417b6490206SBarry Smith           xb = x + 5*(*idx++);
418b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
419b6490206SBarry Smith           sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
420b6490206SBarry Smith           sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
421b6490206SBarry Smith           sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
422b6490206SBarry Smith           sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
423b6490206SBarry Smith           sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
424b6490206SBarry Smith           v += 25;
425b6490206SBarry Smith         }
426b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5;
427b6490206SBarry Smith         y += 5;
428b6490206SBarry Smith       }
429b6490206SBarry Smith       break;
430b6490206SBarry Smith       /* block sizes larger then 5 by 5 are handled by BLAS */
431b6490206SBarry Smith     default: {
432de6a44a3SBarry Smith       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
433de6a44a3SBarry Smith       if (!a->mult_work) {
434de6a44a3SBarry Smith         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
435de6a44a3SBarry Smith         CHKPTRQ(a->mult_work);
436de6a44a3SBarry Smith       }
437de6a44a3SBarry Smith       work = a->mult_work;
438b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
439de6a44a3SBarry Smith         n     = ii[1] - ii[0]; ii++;
440de6a44a3SBarry Smith         ncols = n*bs;
441de6a44a3SBarry Smith         workt = work;
442b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
443b6490206SBarry Smith           xb = x + bs*(*idx++);
444de6a44a3SBarry Smith           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
445de6a44a3SBarry Smith           workt += bs;
446b6490206SBarry Smith         }
447de6a44a3SBarry Smith         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One);
448de6a44a3SBarry Smith         v += n*bs2;
449b6490206SBarry Smith         y += bs;
4502593348eSBarry Smith       }
4512593348eSBarry Smith     }
4522593348eSBarry Smith   }
453b6490206SBarry Smith   PLogFlops(2*bs2*a->nz - m);
4542593348eSBarry Smith   return 0;
4552593348eSBarry Smith }
4562593348eSBarry Smith 
457de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
4582593348eSBarry Smith {
459b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
460bcd2baecSBarry Smith   if (nz)  *nz  = a->bs*a->bs*a->nz;
461bcd2baecSBarry Smith   if (nza) *nza = a->maxnz;
462bcd2baecSBarry Smith   if (mem) *mem = (int)A->mem;
4632593348eSBarry Smith   return 0;
4642593348eSBarry Smith }
4652593348eSBarry Smith 
46691d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
46791d316f6SSatish Balay {
46891d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
46991d316f6SSatish Balay 
47091d316f6SSatish Balay   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
47191d316f6SSatish Balay 
47291d316f6SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
47391d316f6SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
47491d316f6SSatish Balay       (a->nz != b->nz)||(a->indexshift != b->indexshift)) {
47591d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
47691d316f6SSatish Balay   }
47791d316f6SSatish Balay 
47891d316f6SSatish Balay   /* if the a->i are the same */
47991d316f6SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
48091d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
48191d316f6SSatish Balay   }
48291d316f6SSatish Balay 
48391d316f6SSatish Balay   /* if a->j are the same */
48491d316f6SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
48591d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
48691d316f6SSatish Balay   }
48791d316f6SSatish Balay 
48891d316f6SSatish Balay   /* if a->a are the same */
48991d316f6SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
49091d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
49191d316f6SSatish Balay   }
49291d316f6SSatish Balay   *flg = PETSC_TRUE;
49391d316f6SSatish Balay   return 0;
49491d316f6SSatish Balay 
49591d316f6SSatish Balay }
49691d316f6SSatish Balay 
49791d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
49891d316f6SSatish Balay {
49991d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
50017e48fc4SSatish Balay   int         i,j,k,n,row,bs,*ai,*aj,ambs;
50117e48fc4SSatish Balay   Scalar      *x, zero = 0.0,*aa,*aa_j;
50217e48fc4SSatish Balay 
50317e48fc4SSatish Balay   bs  = a->bs;
50417e48fc4SSatish Balay   aa   = a->a;
50517e48fc4SSatish Balay   ai   = a->i;
50617e48fc4SSatish Balay   aj   = a->j;
50717e48fc4SSatish Balay   ambs = a->mbs;
50891d316f6SSatish Balay 
50991d316f6SSatish Balay   VecSet(&zero,v);
51091d316f6SSatish Balay   VecGetArray(v,&x); VecGetLocalSize(v,&n);
5119867e207SSatish Balay   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector");
51217e48fc4SSatish Balay   for ( i=0; i<ambs; i++ ) {
51317e48fc4SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
51417e48fc4SSatish Balay       if (aj[j] == i) {
51517e48fc4SSatish Balay         row  = i*bs;
51617e48fc4SSatish Balay         aa_j = aa+j*bs*bs;
51717e48fc4SSatish Balay         for (k=0; k<bs*bs; k+=(bs+1),row++) x[row] = aa_j[k];
51891d316f6SSatish Balay         break;
51991d316f6SSatish Balay       }
52091d316f6SSatish Balay     }
52191d316f6SSatish Balay   }
52291d316f6SSatish Balay   return 0;
52391d316f6SSatish Balay }
52491d316f6SSatish Balay 
5259867e207SSatish Balay static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
5269867e207SSatish Balay {
5279867e207SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5289867e207SSatish Balay   Scalar      *l,*r,x,*v,*aa,*li,*ri;
5299867e207SSatish Balay   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs;
5309867e207SSatish Balay 
5319867e207SSatish Balay   ai  = a->i;
5329867e207SSatish Balay   aj  = a->j;
5339867e207SSatish Balay   aa  = a->a;
5349867e207SSatish Balay   m   = a->m;
5359867e207SSatish Balay   n   = a->n;
5369867e207SSatish Balay   bs  = a->bs;
5379867e207SSatish Balay   mbs = a->mbs;
5389867e207SSatish Balay 
5399867e207SSatish Balay   if (ll) {
5409867e207SSatish Balay     VecGetArray(ll,&l); VecGetSize(ll,&lm);
5419867e207SSatish Balay     if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length");
5429867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
5439867e207SSatish Balay       M  = ai[i+1] - ai[i];
5449867e207SSatish Balay       li = l + i*bs;
5459867e207SSatish Balay       v  = aa + bs*bs*ai[i];
5469867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
5479867e207SSatish Balay         for ( k=0; k<bs*bs; k++ ) {
5489867e207SSatish Balay           (*v++) *= li[k%bs];
5499867e207SSatish Balay         }
5509867e207SSatish Balay       }
5519867e207SSatish Balay     }
5529867e207SSatish Balay   }
5539867e207SSatish Balay 
5549867e207SSatish Balay   if (rr) {
5559867e207SSatish Balay     VecGetArray(rr,&r); VecGetSize(rr,&rn);
5569867e207SSatish Balay     if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length");
5579867e207SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
5589867e207SSatish Balay       M  = ai[i+1] - ai[i];
5599867e207SSatish Balay       v  = aa + bs*bs*ai[i];
5609867e207SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
5619867e207SSatish Balay         ri = r + bs*aj[ai[i]+j];
5629867e207SSatish Balay         for ( k=0; k<bs; k++ ) {
5639867e207SSatish Balay           x = ri[k];
5649867e207SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
5659867e207SSatish Balay         }
5669867e207SSatish Balay       }
5679867e207SSatish Balay     }
5689867e207SSatish Balay   }
5699867e207SSatish Balay   return 0;
5709867e207SSatish Balay }
5719867e207SSatish Balay 
5729867e207SSatish Balay 
573b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
574b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
575b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec);
5767fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
5777fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
5787fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
5797fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
5807fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
5817fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
5822593348eSBarry Smith 
583b6490206SBarry Smith static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
5842593348eSBarry Smith {
585b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5862593348eSBarry Smith   *m = a->m; *n = a->n;
5872593348eSBarry Smith   return 0;
5882593348eSBarry Smith }
5892593348eSBarry Smith 
590b6490206SBarry Smith static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
5912593348eSBarry Smith {
592b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
5932593348eSBarry Smith   *m = 0; *n = a->m;
5942593348eSBarry Smith   return 0;
5952593348eSBarry Smith }
596b6490206SBarry Smith 
597b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
5982593348eSBarry Smith {
599b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
6002593348eSBarry Smith   Scalar      *v = a->a;
6012593348eSBarry Smith   double      sum = 0.0;
602b6490206SBarry Smith   int         i;
6032593348eSBarry Smith 
6042593348eSBarry Smith   if (type == NORM_FROBENIUS) {
6052593348eSBarry Smith     for (i=0; i<a->nz; i++ ) {
6062593348eSBarry Smith #if defined(PETSC_COMPLEX)
6072593348eSBarry Smith       sum += real(conj(*v)*(*v)); v++;
6082593348eSBarry Smith #else
6092593348eSBarry Smith       sum += (*v)*(*v); v++;
6102593348eSBarry Smith #endif
6112593348eSBarry Smith     }
6122593348eSBarry Smith     *norm = sqrt(sum);
6132593348eSBarry Smith   }
6142593348eSBarry Smith   else {
615b6490206SBarry Smith     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
6162593348eSBarry Smith   }
6172593348eSBarry Smith   return 0;
6182593348eSBarry Smith }
6192593348eSBarry Smith 
6202593348eSBarry Smith /*
6212593348eSBarry Smith      note: This can only work for identity for row and col. It would
6222593348eSBarry Smith    be good to check this and otherwise generate an error.
6232593348eSBarry Smith */
624b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
6252593348eSBarry Smith {
626b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
6272593348eSBarry Smith   Mat         outA;
628de6a44a3SBarry Smith   int         ierr;
6292593348eSBarry Smith 
630b6490206SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
6312593348eSBarry Smith 
6322593348eSBarry Smith   outA          = inA;
6332593348eSBarry Smith   inA->factor   = FACTOR_LU;
6342593348eSBarry Smith   a->row        = row;
6352593348eSBarry Smith   a->col        = col;
6362593348eSBarry Smith 
637de6a44a3SBarry Smith   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
6382593348eSBarry Smith 
6392593348eSBarry Smith   if (!a->diag) {
640b6490206SBarry Smith     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
6412593348eSBarry Smith   }
6427fc0212eSBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
6432593348eSBarry Smith   return 0;
6442593348eSBarry Smith }
6452593348eSBarry Smith 
646b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
6472593348eSBarry Smith {
648b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
649b6490206SBarry Smith   int         one = 1, totalnz = a->bs*a->bs*a->nz;
650b6490206SBarry Smith   BLscal_( &totalnz, alpha, a->a, &one );
651b6490206SBarry Smith   PLogFlops(totalnz);
6522593348eSBarry Smith   return 0;
6532593348eSBarry Smith }
6542593348eSBarry Smith 
655b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A)
6562593348eSBarry Smith {
6572593348eSBarry Smith   static int called = 0;
6582593348eSBarry Smith 
6592593348eSBarry Smith   if (called) return 0; else called = 1;
6602593348eSBarry Smith   return 0;
6612593348eSBarry Smith }
6622593348eSBarry Smith /* -------------------------------------------------------------------*/
663b6490206SBarry Smith static struct _MatOps MatOps = {0,
6649867e207SSatish Balay        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
665b6490206SBarry Smith        MatMult_SeqBAIJ,0,
666b6490206SBarry Smith        0,0,
667de6a44a3SBarry Smith        MatSolve_SeqBAIJ,0,
668b6490206SBarry Smith        0,0,
669de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
670b6490206SBarry Smith        0,
671b6490206SBarry Smith        0,
67217e48fc4SSatish Balay        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
6739867e207SSatish Balay        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
674b6490206SBarry Smith        0,0,
675b6490206SBarry Smith        0,
676b6490206SBarry Smith        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
677b6490206SBarry Smith        MatGetReordering_SeqBAIJ,
6787fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
679b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
680de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
681b6490206SBarry Smith        0,0,/*  MatConvert_SeqBAIJ  */ 0,
682b6490206SBarry Smith        0,0,
683b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
684b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
685b6490206SBarry Smith        0,0,
686b6490206SBarry Smith        0,0,
687b6490206SBarry Smith        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
688b6490206SBarry Smith        0};
6892593348eSBarry Smith 
6902593348eSBarry Smith /*@C
691b6490206SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format
6922593348eSBarry Smith    (the default parallel PETSc format).  For good matrix assembly performance
6932593348eSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
6942593348eSBarry Smith    (or nzz).  By setting these parameters accurately, performance can be
6952593348eSBarry Smith    increased by more than a factor of 50.
6962593348eSBarry Smith 
6972593348eSBarry Smith    Input Parameters:
6982593348eSBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
699b6490206SBarry Smith .  bs - size of block
7002593348eSBarry Smith .  m - number of rows
7012593348eSBarry Smith .  n - number of columns
702b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
703b6490206SBarry Smith .  nzz - number of block nonzeros per block row or PETSC_NULL
704b6490206SBarry Smith          (possibly different for each row)
7052593348eSBarry Smith 
7062593348eSBarry Smith    Output Parameter:
7072593348eSBarry Smith .  A - the matrix
7082593348eSBarry Smith 
7092593348eSBarry Smith    Notes:
710b6490206SBarry Smith    The BAIJ format, is fully compatible with standard Fortran 77
7112593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
7122593348eSBarry Smith    either one (as in Fortran) or zero.  See the users manual for details.
7132593348eSBarry Smith 
7142593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
7152593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
7162593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
7172593348eSBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
7182593348eSBarry Smith 
7192593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
7202593348eSBarry Smith @*/
721b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
7222593348eSBarry Smith {
7232593348eSBarry Smith   Mat         B;
724b6490206SBarry Smith   Mat_SeqBAIJ *b;
725b6490206SBarry Smith   int         i,len,ierr, flg,mbs = m/bs;
726b6490206SBarry Smith 
727b6490206SBarry Smith   if (mbs*bs != m)
728b6490206SBarry Smith     SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize");
7292593348eSBarry Smith 
7302593348eSBarry Smith   *A      = 0;
731b6490206SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
7322593348eSBarry Smith   PLogObjectCreate(B);
733b6490206SBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
7342593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
7357fc0212eSBarry Smith   switch (bs) {
7367fc0212eSBarry Smith     case 1:
7377fc0212eSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
7387fc0212eSBarry Smith       break;
7394eeb42bcSBarry Smith     case 2:
7404eeb42bcSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
7416d84be18SBarry Smith       break;
742f327dd97SBarry Smith     case 3:
743f327dd97SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
7444eeb42bcSBarry Smith       break;
7456d84be18SBarry Smith     case 4:
7466d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
7476d84be18SBarry Smith       break;
7486d84be18SBarry Smith     case 5:
7496d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
7506d84be18SBarry Smith       break;
7517fc0212eSBarry Smith   }
752b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
753b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
7542593348eSBarry Smith   B->factor           = 0;
7552593348eSBarry Smith   B->lupivotthreshold = 1.0;
7562593348eSBarry Smith   b->row              = 0;
7572593348eSBarry Smith   b->col              = 0;
7582593348eSBarry Smith   b->reallocs         = 0;
7592593348eSBarry Smith 
7602593348eSBarry Smith   b->m       = m;
761b6490206SBarry Smith   b->mbs     = mbs;
7622593348eSBarry Smith   b->n       = n;
763b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
7642593348eSBarry Smith   if (nnz == PETSC_NULL) {
765de6a44a3SBarry Smith     if (nz == PETSC_DEFAULT) nz = 5;
7662593348eSBarry Smith     else if (nz <= 0)        nz = 1;
767b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
768b6490206SBarry Smith     nz = nz*mbs;
7692593348eSBarry Smith   }
7702593348eSBarry Smith   else {
7712593348eSBarry Smith     nz = 0;
772b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
7732593348eSBarry Smith   }
7742593348eSBarry Smith 
7752593348eSBarry Smith   /* allocate the matrix space */
776b6490206SBarry Smith   len   = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int);
7772593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
778b6490206SBarry Smith   PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar));
779b6490206SBarry Smith   b->j  = (int *) (b->a + nz*bs*bs);
7802593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
7812593348eSBarry Smith   b->i  = b->j + nz;
7822593348eSBarry Smith   b->singlemalloc = 1;
7832593348eSBarry Smith 
784de6a44a3SBarry Smith   b->i[0] = 0;
785b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
7862593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
7872593348eSBarry Smith   }
7882593348eSBarry Smith 
789b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
790b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
791b6490206SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
792b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
7932593348eSBarry Smith 
794b6490206SBarry Smith   b->bs               = bs;
795b6490206SBarry Smith   b->mbs              = mbs;
7962593348eSBarry Smith   b->nz               = 0;
7972593348eSBarry Smith   b->maxnz            = nz;
7982593348eSBarry Smith   b->sorted           = 0;
7992593348eSBarry Smith   b->roworiented      = 1;
8002593348eSBarry Smith   b->nonew            = 0;
8012593348eSBarry Smith   b->diag             = 0;
8022593348eSBarry Smith   b->solve_work       = 0;
803de6a44a3SBarry Smith   b->mult_work        = 0;
8042593348eSBarry Smith   b->spptr            = 0;
8052593348eSBarry Smith 
8062593348eSBarry Smith   *A = B;
8072593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
8082593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
8092593348eSBarry Smith   return 0;
8102593348eSBarry Smith }
8112593348eSBarry Smith 
812b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
8132593348eSBarry Smith {
8142593348eSBarry Smith   Mat         C;
815b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
816de6a44a3SBarry Smith   int         i,len, mbs = a->mbs, bs = a->bs,nz = a->nz;
817de6a44a3SBarry Smith 
818de6a44a3SBarry Smith   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
8192593348eSBarry Smith 
8202593348eSBarry Smith   *B = 0;
821b6490206SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
8222593348eSBarry Smith   PLogObjectCreate(C);
823b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
8242593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
825b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
826b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
8272593348eSBarry Smith   C->factor     = A->factor;
8282593348eSBarry Smith   c->row        = 0;
8292593348eSBarry Smith   c->col        = 0;
8302593348eSBarry Smith   C->assembled  = PETSC_TRUE;
8312593348eSBarry Smith 
8322593348eSBarry Smith   c->m          = a->m;
8332593348eSBarry Smith   c->n          = a->n;
834b6490206SBarry Smith   c->bs         = a->bs;
835b6490206SBarry Smith   c->mbs        = a->mbs;
836b6490206SBarry Smith   c->nbs        = a->nbs;
8372593348eSBarry Smith 
838b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
839b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
840b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
8412593348eSBarry Smith     c->imax[i] = a->imax[i];
8422593348eSBarry Smith     c->ilen[i] = a->ilen[i];
8432593348eSBarry Smith   }
8442593348eSBarry Smith 
8452593348eSBarry Smith   /* allocate the matrix space */
8462593348eSBarry Smith   c->singlemalloc = 1;
847de6a44a3SBarry Smith   len   = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int));
8482593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
849de6a44a3SBarry Smith   c->j  = (int *) (c->a + nz*bs*bs);
850de6a44a3SBarry Smith   c->i  = c->j + nz;
851b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
852b6490206SBarry Smith   if (mbs > 0) {
853de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
8542593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
855de6a44a3SBarry Smith       PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar));
8562593348eSBarry Smith     }
8572593348eSBarry Smith   }
8582593348eSBarry Smith 
859b6490206SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
8602593348eSBarry Smith   c->sorted      = a->sorted;
8612593348eSBarry Smith   c->roworiented = a->roworiented;
8622593348eSBarry Smith   c->nonew       = a->nonew;
8632593348eSBarry Smith 
8642593348eSBarry Smith   if (a->diag) {
865b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
866b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
867b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
8682593348eSBarry Smith       c->diag[i] = a->diag[i];
8692593348eSBarry Smith     }
8702593348eSBarry Smith   }
8712593348eSBarry Smith   else c->diag          = 0;
8722593348eSBarry Smith   c->nz                 = a->nz;
8732593348eSBarry Smith   c->maxnz              = a->maxnz;
8742593348eSBarry Smith   c->solve_work         = 0;
8752593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
8767fc0212eSBarry Smith   c->mult_work          = 0;
8772593348eSBarry Smith   *B = C;
8782593348eSBarry Smith   return 0;
8792593348eSBarry Smith }
8802593348eSBarry Smith 
88119bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
8822593348eSBarry Smith {
883b6490206SBarry Smith   Mat_SeqBAIJ  *a;
8842593348eSBarry Smith   Mat          B;
885de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
886b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
88735aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
888de6a44a3SBarry Smith   int          *masked, nmask,tmp,ishift,bs2;
889b6490206SBarry Smith   Scalar       *aa;
89019bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
8912593348eSBarry Smith 
89235aab85fSBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
893de6a44a3SBarry Smith   bs2  = bs*bs;
894b6490206SBarry Smith 
8952593348eSBarry Smith   MPI_Comm_size(comm,&size);
896b6490206SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
89790ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
89877c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
899de6a44a3SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
9002593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
9012593348eSBarry Smith 
902b6490206SBarry Smith   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
90335aab85fSBarry Smith 
90435aab85fSBarry Smith   /*
90535aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
90635aab85fSBarry Smith     divisible by the blocksize
90735aab85fSBarry Smith   */
908b6490206SBarry Smith   mbs        = M/bs;
90935aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
91035aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
91135aab85fSBarry Smith   else                  mbs++;
91235aab85fSBarry Smith   if (extra_rows) {
91335aab85fSBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize");
91435aab85fSBarry Smith   }
915b6490206SBarry Smith 
9162593348eSBarry Smith   /* read in row lengths */
91735aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
91877c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
91935aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
9202593348eSBarry Smith 
921b6490206SBarry Smith   /* read in column indices */
92235aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
92377c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
92435aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
925b6490206SBarry Smith 
926b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
927b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
928b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
92935aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
93035aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
93135aab85fSBarry Smith   masked = mask + mbs;
932b6490206SBarry Smith   rowcount = 0; nzcount = 0;
933b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
93435aab85fSBarry Smith     nmask = 0;
935b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
936b6490206SBarry Smith       kmax = rowlengths[rowcount];
937b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
93835aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
93935aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
940b6490206SBarry Smith       }
941b6490206SBarry Smith       rowcount++;
942b6490206SBarry Smith     }
94335aab85fSBarry Smith     browlengths[i] += nmask;
94435aab85fSBarry Smith     /* zero out the mask elements we set */
94535aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
946b6490206SBarry Smith   }
947b6490206SBarry Smith 
9482593348eSBarry Smith   /* create our matrix */
94935aab85fSBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
95035aab85fSBarry Smith          CHKERRQ(ierr);
9512593348eSBarry Smith   B = *A;
952b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
9532593348eSBarry Smith 
9542593348eSBarry Smith   /* set matrix "i" values */
955de6a44a3SBarry Smith   a->i[0] = 0;
956b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
957b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
958b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
9592593348eSBarry Smith   }
9609867e207SSatish Balay   a->indexshift = 0;
961b6490206SBarry Smith   a->nz         = 0;
962b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
9632593348eSBarry Smith 
964b6490206SBarry Smith   /* read in nonzero values */
96535aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
96677c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
96735aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
968b6490206SBarry Smith 
969b6490206SBarry Smith   /* set "a" and "j" values into matrix */
970b6490206SBarry Smith   nzcount = 0; jcount = 0;
971b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
972b6490206SBarry Smith     nzcountb = nzcount;
97335aab85fSBarry Smith     nmask    = 0;
974b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
975b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
976b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
97735aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
97835aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
979b6490206SBarry Smith       }
980b6490206SBarry Smith       rowcount++;
981b6490206SBarry Smith     }
982de6a44a3SBarry Smith     /* sort the masked values */
98377c4ece6SBarry Smith     PetscSortInt(nmask,masked);
984de6a44a3SBarry Smith 
985b6490206SBarry Smith     /* set "j" values into matrix */
986b6490206SBarry Smith     maskcount = 1;
98735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
98835aab85fSBarry Smith       a->j[jcount++]  = masked[j];
989de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
990b6490206SBarry Smith     }
991b6490206SBarry Smith     /* set "a" values into matrix */
992de6a44a3SBarry Smith     ishift = bs2*a->i[i];
993b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
994b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
995b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
996de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
997de6a44a3SBarry Smith         block  = mask[tmp] - 1;
998de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
999de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
1000b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
1001b6490206SBarry Smith       }
1002b6490206SBarry Smith     }
100335aab85fSBarry Smith     /* zero out the mask elements we set */
100435aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1005b6490206SBarry Smith   }
100635aab85fSBarry Smith   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
1007b6490206SBarry Smith 
1008b6490206SBarry Smith   PetscFree(rowlengths);
1009b6490206SBarry Smith   PetscFree(browlengths);
1010b6490206SBarry Smith   PetscFree(aa);
1011b6490206SBarry Smith   PetscFree(jj);
1012b6490206SBarry Smith   PetscFree(mask);
1013b6490206SBarry Smith 
1014b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1015de6a44a3SBarry Smith 
1016de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
1017de6a44a3SBarry Smith   if (flg) {
101819bcc07fSBarry Smith     Viewer tviewer;
101919bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
102090ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr);
102119bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
102219bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1023de6a44a3SBarry Smith   }
1024de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
1025de6a44a3SBarry Smith   if (flg) {
102619bcc07fSBarry Smith     Viewer tviewer;
102719bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
102890ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
102919bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
103019bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1031de6a44a3SBarry Smith   }
1032de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
1033de6a44a3SBarry Smith   if (flg) {
103419bcc07fSBarry Smith     Viewer tviewer;
103519bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
103619bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
103719bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1038de6a44a3SBarry Smith   }
1039de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
1040de6a44a3SBarry Smith   if (flg) {
104119bcc07fSBarry Smith     Viewer tviewer;
104219bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
104390ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr);
104419bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
104519bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1046de6a44a3SBarry Smith   }
1047de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
1048de6a44a3SBarry Smith   if (flg) {
104919bcc07fSBarry Smith     Viewer tviewer;
105019bcc07fSBarry Smith     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
105119bcc07fSBarry Smith     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
105219bcc07fSBarry Smith     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
105319bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1054de6a44a3SBarry Smith   }
10552593348eSBarry Smith   return 0;
10562593348eSBarry Smith }
10572593348eSBarry Smith 
10582593348eSBarry Smith 
10592593348eSBarry Smith 
1060