xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 0a8549d7fee0bb766707fd704f758681ba25da76)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: baij.c,v 1.30 1996/04/07 22:45:39 curfman Exp balay $";
4 #endif
5 
6 /*
7     Defines the basic matrix operations for the BAIJ (compressed row)
8   matrix storage format.
9 */
10 #include "baij.h"
11 #include "vec/vecimpl.h"
12 #include "inline/spops.h"
13 #include "petsc.h"
14 
15 extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
16 
17 static int MatGetReordering_SeqBAIJ(Mat A,MatOrdering type,IS *rperm,IS *cperm)
18 {
19   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
20   int         ierr, *ia, *ja,n = a->mbs,*idx,i,ishift,oshift;
21 
22   /*
23      this is tacky: In the future when we have written special factorization
24      and solve routines for the identity permutation we should use a
25      stride index set instead of the general one.
26   */
27   if (type  == ORDER_NATURAL) {
28     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
29     for ( i=0; i<n; i++ ) idx[i] = i;
30     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
31     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
32     PetscFree(idx);
33     ISSetPermutation(*rperm);
34     ISSetPermutation(*cperm);
35     ISSetIdentity(*rperm);
36     ISSetIdentity(*cperm);
37     return 0;
38   }
39 
40   MatReorderingRegisterAll();
41   ishift = a->indexshift;
42   oshift = -MatReorderingIndexShift[(int)type];
43   if (MatReorderingRequiresSymmetric[(int)type]) {
44     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja);
45     CHKERRQ(ierr);
46     ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
47     PetscFree(ia); PetscFree(ja);
48   } else {
49     if (ishift == oshift) {
50       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
51     }
52     else if (ishift == -1) {
53       /* temporarily subtract 1 from i and j indices */
54       int nz = a->i[a->n] - 1;
55       for ( i=0; i<nz; i++ ) a->j[i]--;
56       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
57       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
58       for ( i=0; i<nz; i++ ) a->j[i]++;
59       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
60     } else {
61       /* temporarily add 1 to i and j indices */
62       int nz = a->i[a->n] - 1;
63       for ( i=0; i<nz; i++ ) a->j[i]++;
64       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
65       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
66       for ( i=0; i<nz; i++ ) a->j[i]--;
67       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
68     }
69   }
70   return 0;
71 }
72 
73 /*
74      Adds diagonal pointers to sparse matrix structure.
75 */
76 
77 int MatMarkDiag_SeqBAIJ(Mat A)
78 {
79   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
80   int         i,j, *diag, m = a->mbs;
81 
82   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
83   PLogObjectMemory(A,(m+1)*sizeof(int));
84   for ( i=0; i<m; i++ ) {
85     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
86       if (a->j[j] == i) {
87         diag[i] = j;
88         break;
89       }
90     }
91   }
92   a->diag = diag;
93   return 0;
94 }
95 
96 #include "draw.h"
97 #include "pinclude/pviewer.h"
98 #include "sys.h"
99 
100 static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
101 {
102   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
103   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l,bs2=bs*bs;
104   Scalar      *aa;
105 
106   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
107   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
108   col_lens[0] = MAT_COOKIE;
109   col_lens[1] = a->m;
110   col_lens[2] = a->n;
111   col_lens[3] = a->nz*bs2;
112 
113   /* store lengths of each row and write (including header) to file */
114   count = 0;
115   for ( i=0; i<a->mbs; i++ ) {
116     for ( j=0; j<bs; j++ ) {
117       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
118     }
119   }
120   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
121   PetscFree(col_lens);
122 
123   /* store column indices (zero start index) */
124   jj = (int *) PetscMalloc( a->nz*bs2*sizeof(int) ); CHKPTRQ(jj);
125   count = 0;
126   for ( i=0; i<a->mbs; i++ ) {
127     for ( j=0; j<bs; j++ ) {
128       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
129         for ( l=0; l<bs; l++ ) {
130           jj[count++] = bs*a->j[k] + l;
131         }
132       }
133     }
134   }
135   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,BINARY_INT,0); CHKERRQ(ierr);
136   PetscFree(jj);
137 
138   /* store nonzero values */
139   aa = (Scalar *) PetscMalloc(a->nz*bs2*sizeof(Scalar)); CHKPTRQ(aa);
140   count = 0;
141   for ( i=0; i<a->mbs; i++ ) {
142     for ( j=0; j<bs; j++ ) {
143       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
144         for ( l=0; l<bs; l++ ) {
145           aa[count++] = a->a[bs2*k + l*bs + j];
146         }
147       }
148     }
149   }
150   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
151   PetscFree(aa);
152   return 0;
153 }
154 
155 static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
156 {
157   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
158   int         ierr, i,j,format,bs = a->bs,k,l,bs2=bs*bs;
159   FILE        *fd;
160   char        *outputname;
161 
162   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
163   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
164   ierr = ViewerGetFormat(viewer,&format);
165   if (format == ASCII_FORMAT_INFO) {
166     /* no need to print additional information */ ;
167   }
168   else if (format == ASCII_FORMAT_MATLAB) {
169     SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported");
170   }
171   else if (format == ASCII_FORMAT_COMMON) {
172     for ( i=0; i<a->mbs; i++ ) {
173       for ( j=0; j<bs; j++ ) {
174         fprintf(fd,"row %d:",i*bs+j);
175         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
176           for ( l=0; l<bs; l++ ) {
177 #if defined(PETSC_COMPLEX)
178           if (imag(a->a[bs2*k + l*bs + j]) != 0.0 && real(a->a[bs2*k + l*bs + j]) != 0.0)
179             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
180               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
181           else if (real(a->a[bs2*k + l*bs + j]) != 0.0)
182             fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
183 #else
184           if (a->a[bs2*k + l*bs + j] != 0.0)
185             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
186 #endif
187           }
188         }
189         fprintf(fd,"\n");
190       }
191     }
192   }
193   else {
194     for ( i=0; i<a->mbs; i++ ) {
195       for ( j=0; j<bs; j++ ) {
196         fprintf(fd,"row %d:",i*bs+j);
197         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
198           for ( l=0; l<bs; l++ ) {
199 #if defined(PETSC_COMPLEX)
200           if (imag(a->a[bs2*k + l*bs + j]) != 0.0) {
201             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
202               real(a->a[bs2*k + l*bs + j]),imag(a->a[bs2*k + l*bs + j]));
203           }
204           else {
205             fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs2*k + l*bs + j]));
206           }
207 #else
208             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);
209 #endif
210           }
211         }
212         fprintf(fd,"\n");
213       }
214     }
215   }
216   fflush(fd);
217   return 0;
218 }
219 
220 static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
221 {
222   Mat         A = (Mat) obj;
223   ViewerType  vtype;
224   int         ierr;
225 
226   if (!viewer) {
227     viewer = STDOUT_VIEWER_SELF;
228   }
229 
230   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
231   if (vtype == MATLAB_VIEWER) {
232     SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported");
233   }
234   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
235     return MatView_SeqBAIJ_ASCII(A,viewer);
236   }
237   else if (vtype == BINARY_FILE_VIEWER) {
238     return MatView_SeqBAIJ_Binary(A,viewer);
239   }
240   else if (vtype == DRAW_VIEWER) {
241     SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported");
242   }
243   return 0;
244 }
245 
246 #define CHUNKSIZE  10
247 
248 /* This version has row oriented v  */
249 static int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
250 {
251   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
252   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
253   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
254   int         *aj=a->j,nonew=a->nonew,shift=a->indexshift,bs=a->bs,brow,bcol;
255   int          ridx,cidx,bs2=bs*bs;
256   Scalar      *ap,value,*aa=a->a,*bap;
257 
258   for ( k=0; k<m; k++ ) { /* loop over added rows */
259     row  = im[k]; brow = row/bs;
260     if (row < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative row");
261     if (row >= a->m) SETERRQ(1,"MatSetValues_SeqBAIJ:Row too large");
262     rp   = aj + ai[brow] + shift; ap = aa + bs2*ai[brow] + shift;
263     rmax = imax[brow]; nrow = ailen[brow];
264     low = 0;
265     for ( l=0; l<n; l++ ) { /* loop over added columns */
266       if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqBAIJ:Negative column");
267       if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqBAIJ:Column too large");
268       col = in[l] - shift; bcol = col/bs;
269       ridx = row % bs; cidx = col % bs;
270       if (roworiented) {
271         value = *v++;
272       }
273       else {
274         value = v[k + l*m];
275       }
276       if (!sorted) low = 0; high = nrow;
277       while (high-low > 5) {
278         t = (low+high)/2;
279         if (rp[t] > bcol) high = t;
280         else             low  = t;
281       }
282       for ( i=low; i<high; i++ ) {
283         if (rp[i] > bcol) break;
284         if (rp[i] == bcol) {
285           bap  = ap +  bs2*i + bs*cidx + ridx;
286           if (is == ADD_VALUES) *bap += value;
287           else                  *bap  = value;
288           goto noinsert;
289         }
290       }
291       if (nonew) goto noinsert;
292       if (nrow >= rmax) {
293         /* there is no extra room in row, therefore enlarge */
294         int    new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
295         Scalar *new_a;
296 
297         /* malloc new storage space */
298         len     = new_nz*(sizeof(int)+bs2*sizeof(Scalar))+(a->mbs+1)*sizeof(int);
299         new_a   = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a);
300         new_j   = (int *) (new_a + bs2*new_nz);
301         new_i   = new_j + new_nz;
302 
303         /* copy over old data into new slots */
304         for ( ii=0; ii<brow+1; ii++ ) {new_i[ii] = ai[ii];}
305         for ( ii=brow+1; ii<a->mbs+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;}
306         PetscMemcpy(new_j,aj,(ai[brow]+nrow+shift)*sizeof(int));
307         len = (new_nz - CHUNKSIZE - ai[brow] - nrow - shift);
308         PetscMemcpy(new_j+ai[brow]+shift+nrow+CHUNKSIZE,aj+ai[brow]+shift+nrow,
309                                                            len*sizeof(int));
310         PetscMemcpy(new_a,aa,(ai[brow]+nrow+shift)*bs2*sizeof(Scalar));
311         PetscMemzero(new_a+bs2*(ai[brow]+shift+nrow),bs2*CHUNKSIZE*sizeof(Scalar));
312         PetscMemcpy(new_a+bs2*(ai[brow]+shift+nrow+CHUNKSIZE),
313                     aa+bs2*(ai[brow]+shift+nrow),bs2*len*sizeof(Scalar));
314         /* free up old matrix storage */
315         PetscFree(a->a);
316         if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);}
317         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
318         a->singlemalloc = 1;
319 
320         rp   = aj + ai[brow] + shift; ap = aa + bs2*ai[brow] + shift;
321         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
322         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(Scalar)));
323         a->maxnz += CHUNKSIZE;
324         a->reallocs++;
325         a->nz++;
326       }
327       N = nrow++ - 1;
328       /* shift up all the later entries in this row */
329       for ( ii=N; ii>=i; ii-- ) {
330         rp[ii+1] = rp[ii];
331          PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(Scalar));
332       }
333       if (N>=i) PetscMemzero(ap+bs2*i,bs2*sizeof(Scalar));
334       rp[i] = bcol;
335       ap[bs2*i + bs*cidx + ridx] = value;
336       noinsert:;
337       low = i;
338     }
339     ailen[brow] = nrow;
340   }
341   return 0;
342 }
343 
344 static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
345 {
346   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
347   *m = a->m; *n = a->n;
348   return 0;
349 }
350 
351 static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
352 {
353   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
354   *m = 0; *n = a->m;
355   return 0;
356 }
357 
358 int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
359 {
360   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
361   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2;
362   Scalar      *aa,*v_i,*aa_i;
363 
364   bs  = a->bs;
365   ai  = a->i;
366   aj  = a->j;
367   aa  = a->a;
368   bs2 = bs*bs;
369 
370   if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqBAIJ:Row out of range");
371 
372   bn  = row/bs;   /* Block number */
373   bp  = row % bs; /* Block Position */
374   M   = ai[bn+1] - ai[bn];
375   *nz = bs*M;
376 
377   if (v) {
378     *v = 0;
379     if (*nz) {
380       *v = (Scalar *) PetscMalloc( (*nz)*sizeof(Scalar) ); CHKPTRQ(*v);
381       for ( i=0; i<M; i++ ) { /* for each block in the block row */
382         v_i  = *v + i*bs;
383         aa_i = aa + bs2*(ai[bn] + i);
384         for ( j=bp,k=0; j<bs2; j+=bs,k++ ) {v_i[k] = aa_i[j];}
385       }
386     }
387   }
388 
389   if (idx) {
390     *idx = 0;
391     if (*nz) {
392       *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx);
393       for ( i=0; i<M; i++ ) { /* for each block in the block row */
394         idx_i = *idx + i*bs;
395         itmp  = bs*aj[ai[bn] + i];
396         for ( j=0; j<bs; j++ ) {idx_i[j] = itmp++;}
397       }
398     }
399   }
400   return 0;
401 }
402 
403 int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
404 {
405   if (idx) {if (*idx) PetscFree(*idx);}
406   if (v)   {if (*v)   PetscFree(*v);}
407   return 0;
408 }
409 
410 static int MatTranspose_SeqBAIJ(Mat A,Mat *B)
411 {
412   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
413   Mat         C;
414   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
415   int         shift=a->indexshift,*rows,*cols,bs2=bs*bs;
416   Scalar      *array=a->a;
417 
418   if (B==PETSC_NULL && mbs!=nbs)
419     SETERRQ(1,"MatTranspose_SeqBAIJ:Square matrix only for in-place");
420   col = (int *) PetscMalloc((1+nbs)*sizeof(int)); CHKPTRQ(col);
421   PetscMemzero(col,(1+nbs)*sizeof(int));
422   if (shift) {
423     for ( i=0; i<ai[mbs]-1; i++ ) aj[i] -= 1;
424   }
425   for ( i=0; i<ai[mbs]+shift; i++ ) col[aj[i]] += 1;
426   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C); CHKERRQ(ierr);
427   PetscFree(col);
428   rows = (int *) PetscMalloc(2*bs*sizeof(int)); CHKPTRQ(rows);
429   cols = rows + bs;
430   for ( i=0; i<mbs; i++ ) {
431     cols[0] = i*bs;
432     for (k=1; k<bs; k++ ) cols[k] = cols[k-1] + 1;
433     len = ai[i+1] - ai[i];
434     for ( j=0; j<len; j++ ) {
435       rows[0] = (*aj++)*bs;
436       for (k=1; k<bs; k++ ) rows[k] = rows[k-1] + 1;
437       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES); CHKERRQ(ierr);
438       array += bs2;
439     }
440   }
441   PetscFree(rows);
442   if (shift) {
443     for ( i=0; i<ai[mbs]-1; i++ ) aj[i] += 1;
444   }
445 
446   ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
447   ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr);
448 
449   if (B != PETSC_NULL) {
450     *B = C;
451   } else {
452     /* This isn't really an in-place transpose */
453     PetscFree(a->a);
454     if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);}
455     if (a->diag) PetscFree(a->diag);
456     if (a->ilen) PetscFree(a->ilen);
457     if (a->imax) PetscFree(a->imax);
458     if (a->solve_work) PetscFree(a->solve_work);
459     PetscFree(a);
460     PetscMemcpy(A,C,sizeof(struct _Mat));
461     PetscHeaderDestroy(C);
462   }
463   return 0;
464 }
465 
466 
467 static int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
468 {
469   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
470   int        fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax;
471   int        m = a->m,*ip, N, *ailen = a->ilen,shift = a->indexshift;
472   int        bs = a->bs,  mbs = a->mbs, bs2 = bs*bs;
473   Scalar     *aa = a->a, *ap;
474 
475   if (mode == FLUSH_ASSEMBLY) return 0;
476 
477   for ( i=1; i<mbs; i++ ) {
478     /* move each row back by the amount of empty slots (fshift) before it*/
479     fshift += imax[i-1] - ailen[i-1];
480     if (fshift) {
481       ip = aj + ai[i] + shift; ap = aa + bs2*ai[i] + shift;
482       N = ailen[i];
483       for ( j=0; j<N; j++ ) {
484         ip[j-fshift] = ip[j];
485         PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(Scalar));
486       }
487     }
488     ai[i] = ai[i-1] + ailen[i-1];
489   }
490   if (mbs) {
491     fshift += imax[mbs-1] - ailen[mbs-1];
492     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
493   }
494   /* reset ilen and imax for each row */
495   for ( i=0; i<mbs; i++ ) {
496     ailen[i] = imax[i] = ai[i+1] - ai[i];
497   }
498   a->nz = ai[mbs] + shift;
499 
500   /* diagonals may have moved, so kill the diagonal pointers */
501   if (fshift && a->diag) {
502     PetscFree(a->diag);
503     PLogObjectMemory(A,-(m+1)*sizeof(int));
504     a->diag = 0;
505   }
506   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Unneeded storage space %d used %d rows %d\n",
507            fshift,a->nz,m);
508   PLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues %d\n",
509            a->reallocs);
510   return 0;
511 }
512 
513 static int MatZeroEntries_SeqBAIJ(Mat A)
514 {
515   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
516   PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar));
517   return 0;
518 }
519 
520 int MatDestroy_SeqBAIJ(PetscObject obj)
521 {
522   Mat         A  = (Mat) obj;
523   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
524 
525 #if defined(PETSC_LOG)
526   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
527 #endif
528   PetscFree(a->a);
529   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
530   if (a->diag) PetscFree(a->diag);
531   if (a->ilen) PetscFree(a->ilen);
532   if (a->imax) PetscFree(a->imax);
533   if (a->solve_work) PetscFree(a->solve_work);
534   if (a->mult_work) PetscFree(a->mult_work);
535   PetscFree(a);
536   PLogObjectDestroy(A);
537   PetscHeaderDestroy(A);
538   return 0;
539 }
540 
541 static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
542 {
543   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
544   if      (op == ROW_ORIENTED)              a->roworiented = 1;
545   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
546   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
547   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
548   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
549   else if (op == ROWS_SORTED ||
550            op == SYMMETRIC_MATRIX ||
551            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
552            op == YES_NEW_DIAGONALS)
553     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
554   else if (op == NO_NEW_DIAGONALS)
555     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
556   else
557     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
558   return 0;
559 }
560 
561 
562 /* -------------------------------------------------------*/
563 /* Should check that shapes of vectors and matrices match */
564 /* -------------------------------------------------------*/
565 #include "pinclude/plapack.h"
566 
567 static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy)
568 {
569   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
570   Scalar          *xg,*yg;
571   register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5;
572   register Scalar x1,x2,x3,x4,x5;
573   int             mbs = a->mbs, m = a->m, i, *idx,*ii;
574   int             bs = a->bs,j,n,bs2 = bs*bs;
575 
576   VecGetArray(xx,&xg); x = xg;  VecGetArray(yy,&yg); y = yg;
577   PetscMemzero(y,m*sizeof(Scalar));
578   /*
579 x     = x; */
580   idx   = a->j;
581   v     = a->a;
582   ii    = a->i;
583 
584   switch (bs) {
585     case 1:
586      for ( i=0; i<mbs; i++ ) {
587        n    = ii[1] - ii[0]; ii++;
588        sum  = 0.0;
589        while (n--) sum += *v++ * x[*idx++];
590        y[i] = sum;
591       }
592       break;
593     case 2:
594       for ( i=0; i<mbs; i++ ) {
595         n  = ii[1] - ii[0]; ii++;
596         sum1 = 0.0; sum2 = 0.0;
597         for ( j=0; j<n; j++ ) {
598           xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
599           sum1 += v[0]*x1 + v[2]*x2;
600           sum2 += v[1]*x1 + v[3]*x2;
601           v += 4;
602         }
603         y[0] = sum1; y[1] = sum2;
604         y += 2;
605       }
606       break;
607     case 3:
608       for ( i=0; i<mbs; i++ ) {
609         n  = ii[1] - ii[0]; ii++;
610         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
611         for ( j=0; j<n; j++ ) {
612           xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
613           sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
614           sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
615           sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
616           v += 9;
617         }
618         y[0] = sum1; y[1] = sum2; y[2] = sum3;
619         y += 3;
620       }
621       break;
622     case 4:
623       for ( i=0; i<mbs; i++ ) {
624         n  = ii[1] - ii[0]; ii++;
625         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
626         for ( j=0; j<n; j++ ) {
627           xb = x + 4*(*idx++);
628           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
629           sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
630           sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
631           sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
632           sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
633           v += 16;
634         }
635         y[0] = sum1; y[1] = sum2; y[2] = sum3; y[3] = sum4;
636         y += 4;
637       }
638       break;
639     case 5:
640       for ( i=0; i<mbs; i++ ) {
641         n  = ii[1] - ii[0]; ii++;
642         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
643         for ( j=0; j<n; j++ ) {
644           xb = x + 5*(*idx++);
645           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
646           sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
647           sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
648           sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
649           sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
650           sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
651           v += 25;
652         }
653         y[0] = sum1; y[1] = sum2; y[2] = sum3; y[3] = sum4; y[4] = sum5;
654         y += 5;
655       }
656       break;
657       /* block sizes larger then 5 by 5 are handled by BLAS */
658     default: {
659       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
660       if (!a->mult_work) {
661         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
662         CHKPTRQ(a->mult_work);
663       }
664       work = a->mult_work;
665       for ( i=0; i<mbs; i++ ) {
666         n     = ii[1] - ii[0]; ii++;
667         ncols = n*bs;
668         workt = work;
669         for ( j=0; j<n; j++ ) {
670           xb = x + bs*(*idx++);
671           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
672           workt += bs;
673         }
674         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One);
675         v += n*bs2;
676         y += bs;
677       }
678     }
679   }
680   PLogFlops(2*bs2* (a->nz) - m);
681   return 0;
682 }
683 
684 static int MatMultAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
685 {
686   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
687   Scalar          *xg,*yg,*zg;
688   register Scalar *x,*y,*z,*v,sum,*xb,sum1,sum2,sum3,sum4,sum5;
689   register Scalar x1,x2,x3,x4,x5;
690   int             mbs=a->mbs,m=a->m,i,*idx,*ii;
691   int             bs=a->bs,j,n,bs2=bs*bs;
692 
693   VecGetArray(xx,&xg); x = xg;
694   VecGetArray(yy,&yg); y = yg;
695   VecGetArray(zz,&zg); z = zg;
696 
697   if (yy==PETSC_NULL) PetscMemzero(z,m*sizeof(Scalar));
698   else if (zz!=yy) PetscMemcpy(z,y,m*sizeof(Scalar));
699 
700   idx   = a->j;
701   v     = a->a;
702   ii    = a->i;
703 
704   switch (bs) {
705   case 1:
706     for ( i=0; i<mbs; i++ ) {
707       n    = ii[1] - ii[0]; ii++;
708       sum  = 0.0;
709       while (n--) sum += *v++ * x[*idx++];
710       z[i] += sum;
711     }
712     break;
713   case 2:
714     for ( i=0; i<mbs; i++ ) {
715       n  = ii[1] - ii[0]; ii++;
716       sum1 = 0.0; sum2 = 0.0;
717       for ( j=0; j<n; j++ ) {
718         xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
719         sum1 += v[0]*x1 + v[2]*x2;
720         sum2 += v[1]*x1 + v[3]*x2;
721         v += 4;
722       }
723       z[0] += sum1; z[1] += sum2;
724       z += 2;
725     }
726     break;
727   case 3:
728     for ( i=0; i<mbs; i++ ) {
729       n  = ii[1] - ii[0]; ii++;
730       sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
731       for ( j=0; j<n; j++ ) {
732         xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
733         sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
734         sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
735         sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
736         v += 9;
737       }
738       z[0] += sum1; z[1] += sum2; z[2] += sum3;
739       z += 3;
740     }
741     break;
742   case 4:
743     for ( i=0; i<mbs; i++ ) {
744       n  = ii[1] - ii[0]; ii++;
745       sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
746       for ( j=0; j<n; j++ ) {
747         xb = x + 4*(*idx++);
748         x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
749         sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
750         sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
751         sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
752         sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
753         v += 16;
754       }
755       z[0] += sum1; z[1] += sum2; z[2] += sum3; z[3] += sum4;
756       z += 4;
757     }
758     break;
759   case 5:
760     for ( i=0; i<mbs; i++ ) {
761       n  = ii[1] - ii[0]; ii++;
762       sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
763       for ( j=0; j<n; j++ ) {
764         xb = x + 5*(*idx++);
765         x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
766         sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
767         sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
768         sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
769         sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
770         sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
771         v += 25;
772       }
773       z[0] += sum1; z[1] += sum2; z[2] += sum3; z[3] += sum4; z[4] += sum5;
774       z += 5;
775     }
776     break;
777     /* block sizes larger then 5 by 5 are handled by BLAS */
778   default: {
779       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
780       if (!a->mult_work) {
781         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
782         CHKPTRQ(a->mult_work);
783       }
784       work = a->mult_work;
785       for ( i=0; i<mbs; i++ ) {
786         n     = ii[1] - ii[0]; ii++;
787         ncols = n*bs;
788         workt = work;
789         for ( j=0; j<n; j++ ) {
790           xb = x + bs*(*idx++);
791           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
792           workt += bs;
793         }
794         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One);
795         v += n*bs2;
796         z += bs;
797       }
798     }
799   }
800   PLogFlops(2*bs2*(a->nz));
801   return 0;
802 }
803 
804 static int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
805 {
806   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
807   Scalar          *xg,*yg,*zg,*zb;
808   register Scalar *x,*y,*z,*v,*xb,x1,x2,x3,x4,x5;
809   int             mbs=a->mbs,N=a->n,i,*idx,*ii,*ai=a->i,rval;
810   int             bs=a->bs,j,n,bs2=bs*bs,*ib;
811 
812   VecGetArray(xx,&xg); x = xg;
813   VecGetArray(yy,&yg); y = yg;
814   VecGetArray(zz,&zg); z = zg;
815 
816   if (yy==PETSC_NULL) PetscMemzero(z,N*sizeof(Scalar));
817   else if (zz!=yy) PetscMemcpy(z,y,N*sizeof(Scalar));
818 
819   idx   = a->j;
820   v     = a->a;
821   ii    = a->i;
822 
823   switch (bs) {
824   case 1:
825     for ( i=0; i<mbs; i++ ) {
826       n  = ii[1] - ii[0]; ii++;
827       xb = x + i; x1 = xb[0];
828       ib = idx + ai[i];
829       for ( j=0; j<n; j++ ) {
830         z[ib[j]] += *v++ * x1;
831       }
832     }
833     break;
834   case 2:
835     for ( i=0; i<mbs; i++ ) {
836       n  = ii[1] - ii[0]; ii++;
837       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
838       ib = idx + ai[i];
839       for ( j=0; j<n; j++ ) {
840         rval = ib[j]*2;
841         z[rval++] += v[0]*x1 + v[1]*x2;
842         z[rval++] += v[2]*x1 + v[3]*x2;
843         v += 4;
844       }
845     }
846     break;
847   case 3:
848     for ( i=0; i<mbs; i++ ) {
849       n  = ii[1] - ii[0]; ii++;
850       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
851       ib = idx + ai[i];
852       for ( j=0; j<n; j++ ) {
853         rval = ib[j]*3;
854         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
855         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
856         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
857         v += 9;
858       }
859     }
860     break;
861   case 4:
862     for ( i=0; i<mbs; i++ ) {
863       n  = ii[1] - ii[0]; ii++;
864       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
865       ib = idx + ai[i];
866       for ( j=0; j<n; j++ ) {
867         rval = ib[j]*4;
868         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
869         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
870         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
871         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
872         v += 16;
873       }
874     }
875     break;
876   case 5:
877     for ( i=0; i<mbs; i++ ) {
878       n  = ii[1] - ii[0]; ii++;
879       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
880       x4 = xb[3];   x5 = xb[4];
881       ib = idx + ai[i];
882       for ( j=0; j<n; j++ ) {
883         rval = ib[j]*5;
884         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
885         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
886         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
887         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
888         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
889         v += 25;
890       }
891     }
892     break;
893       /* block sizes larger then 5 by 5 are handled by BLAS */
894     default: {
895       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
896       if (!a->mult_work) {
897         /* must be max of m,n */
898         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
899         CHKPTRQ(a->mult_work);
900       }
901       work = a->mult_work;
902       for ( i=0; i<mbs; i++ ) {
903         n     = ii[1] - ii[0]; ii++;
904         ncols = n*bs;
905         PetscMemzero(work,ncols*sizeof(Scalar));
906         LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One);
907         v += n*bs2;
908         x += bs;
909         workt = work;
910         for ( j=0; j<n; j++ ) {
911           zb = z + bs*(*idx++);
912           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
913           workt += bs;
914         }
915       }
916     }
917   }
918   PLogFlops(2*bs2*(a->nz));
919   return 0;
920 }
921 
922 static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
923 {
924   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
925   if (nz)  *nz  = a->bs*a->bs*a->nz;
926   if (nza) *nza = a->maxnz;
927   if (mem) *mem = (int)A->mem;
928   return 0;
929 }
930 
931 static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
932 {
933   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
934 
935   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
936 
937   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
938   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
939       (a->nz != b->nz)||(a->indexshift != b->indexshift)) {
940     *flg = PETSC_FALSE; return 0;
941   }
942 
943   /* if the a->i are the same */
944   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
945     *flg = PETSC_FALSE; return 0;
946   }
947 
948   /* if a->j are the same */
949   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
950     *flg = PETSC_FALSE; return 0;
951   }
952 
953   /* if a->a are the same */
954   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
955     *flg = PETSC_FALSE; return 0;
956   }
957   *flg = PETSC_TRUE;
958   return 0;
959 
960 }
961 
962 static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
963 {
964   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
965   int         i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
966   Scalar      *x, zero = 0.0,*aa,*aa_j;
967 
968   bs   = a->bs;
969   aa   = a->a;
970   ai   = a->i;
971   aj   = a->j;
972   ambs = a->mbs;
973   bs2  = bs*bs;
974 
975   VecSet(&zero,v);
976   VecGetArray(v,&x); VecGetLocalSize(v,&n);
977   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqBAIJ:Nonconforming matrix and vector");
978   for ( i=0; i<ambs; i++ ) {
979     for ( j=ai[i]; j<ai[i+1]; j++ ) {
980       if (aj[j] == i) {
981         row  = i*bs;
982         aa_j = aa+j*bs2;
983         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
984         break;
985       }
986     }
987   }
988   return 0;
989 }
990 
991 static int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
992 {
993   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
994   Scalar      *l,*r,x,*v,*aa,*li,*ri;
995   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
996 
997   ai  = a->i;
998   aj  = a->j;
999   aa  = a->a;
1000   m   = a->m;
1001   n   = a->n;
1002   bs  = a->bs;
1003   mbs = a->mbs;
1004   bs2 = bs*bs;
1005   if (ll) {
1006     VecGetArray(ll,&l); VecGetSize(ll,&lm);
1007     if (lm != m) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Left scaling vector wrong length");
1008     for ( i=0; i<mbs; i++ ) { /* for each block row */
1009       M  = ai[i+1] - ai[i];
1010       li = l + i*bs;
1011       v  = aa + bs2*ai[i];
1012       for ( j=0; j<M; j++ ) { /* for each block */
1013         for ( k=0; k<bs2; k++ ) {
1014           (*v++) *= li[k%bs];
1015         }
1016       }
1017     }
1018   }
1019 
1020   if (rr) {
1021     VecGetArray(rr,&r); VecGetSize(rr,&rn);
1022     if (rn != n) SETERRQ(1,"MatDiagonalScale_SeqBAIJ:Right scaling vector wrong length");
1023     for ( i=0; i<mbs; i++ ) { /* for each block row */
1024       M  = ai[i+1] - ai[i];
1025       v  = aa + bs2*ai[i];
1026       for ( j=0; j<M; j++ ) { /* for each block */
1027         ri = r + bs*aj[ai[i]+j];
1028         for ( k=0; k<bs; k++ ) {
1029           x = ri[k];
1030           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
1031         }
1032       }
1033     }
1034   }
1035   return 0;
1036 }
1037 
1038 
1039 extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
1040 extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
1041 extern int MatSolve_SeqBAIJ(Mat,Vec,Vec);
1042 extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
1043 extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
1044 extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
1045 extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
1046 extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
1047 extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
1048 
1049 static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
1050 {
1051   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1052   Scalar      *v = a->a;
1053   double      sum = 0.0;
1054   int         i;
1055 
1056   if (type == NORM_FROBENIUS) {
1057     for (i=0; i<a->nz; i++ ) {
1058 #if defined(PETSC_COMPLEX)
1059       sum += real(conj(*v)*(*v)); v++;
1060 #else
1061       sum += (*v)*(*v); v++;
1062 #endif
1063     }
1064     *norm = sqrt(sum);
1065   }
1066   else {
1067     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
1068   }
1069   return 0;
1070 }
1071 
1072 /*
1073      note: This can only work for identity for row and col. It would
1074    be good to check this and otherwise generate an error.
1075 */
1076 static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
1077 {
1078   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
1079   Mat         outA;
1080   int         ierr;
1081 
1082   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
1083 
1084   outA          = inA;
1085   inA->factor   = FACTOR_LU;
1086   a->row        = row;
1087   a->col        = col;
1088 
1089   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
1090 
1091   if (!a->diag) {
1092     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
1093   }
1094   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
1095   return 0;
1096 }
1097 
1098 static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
1099 {
1100   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
1101   int         one = 1, totalnz = a->bs*a->bs*a->nz;
1102   BLscal_( &totalnz, alpha, a->a, &one );
1103   PLogFlops(totalnz);
1104   return 0;
1105 }
1106 
1107 static int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
1108 {
1109   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1110   int        *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j;
1111   int        *ai = a->i, *ailen = a->ilen, shift = a->indexshift;
1112   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=bs*bs;
1113   Scalar     *ap, *aa = a->a, zero = 0.0;
1114 
1115   for ( k=0; k<m; k++ ) { /* loop over rows */
1116     row  = im[k]; brow = row/bs;
1117     if (row < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative row");
1118     if (row >= a->m) SETERRQ(1,"MatGetValues_SeqBAIJ:Row too large");
1119     rp   = aj + ai[brow] + shift; ap = aa + bs2*ai[brow] + shift;
1120     nrow = ailen[brow];
1121     for ( l=0; l<n; l++ ) { /* loop over columns */
1122       if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqBAIJ:Negative column");
1123       if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqBAIJ:Column too large");
1124       col  = in[l] - shift;
1125       bcol = col/bs;
1126       cidx = col%bs;
1127       ridx = row%bs;
1128       high = nrow;
1129       low  = 0; /* assume unsorted */
1130       while (high-low > 5) {
1131         t = (low+high)/2;
1132         if (rp[t] > bcol) high = t;
1133         else             low  = t;
1134       }
1135       for ( i=low; i<high; i++ ) {
1136         if (rp[i] > bcol) break;
1137         if (rp[i] == bcol) {
1138           *v++ = ap[bs2*i+bs*cidx+ridx];
1139           goto finished;
1140         }
1141       }
1142       *v++ = zero;
1143       finished:;
1144     }
1145   }
1146   return 0;
1147 }
1148 
1149 int MatPrintHelp_SeqBAIJ(Mat A)
1150 {
1151   static int called = 0;
1152 
1153   if (called) return 0; else called = 1;
1154   return 0;
1155 }
1156 /* -------------------------------------------------------------------*/
1157 static struct _MatOps MatOps = {MatSetValues_SeqBAIJ,
1158        MatGetRow_SeqBAIJ,MatRestoreRow_SeqBAIJ,
1159        MatMult_SeqBAIJ,MatMultAdd_SeqBAIJ,
1160        0,MatMultTransAdd_SeqBAIJ,
1161        MatSolve_SeqBAIJ,0,
1162        0,0,
1163        MatLUFactor_SeqBAIJ,0,
1164        0,
1165        MatTranspose_SeqBAIJ,
1166        MatGetInfo_SeqBAIJ,MatEqual_SeqBAIJ,
1167        MatGetDiagonal_SeqBAIJ,MatDiagonalScale_SeqBAIJ,MatNorm_SeqBAIJ,
1168        0,MatAssemblyEnd_SeqBAIJ,
1169        0,
1170        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
1171        MatGetReordering_SeqBAIJ,
1172        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
1173        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
1174        MatILUFactorSymbolic_SeqBAIJ,0,
1175        0,0,/*  MatConvert_SeqBAIJ  */ 0,
1176        0,0,
1177        MatConvertSameType_SeqBAIJ,0,0,
1178        MatILUFactor_SeqBAIJ,0,0,
1179        0,0,
1180        MatGetValues_SeqBAIJ,0,
1181        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
1182        0};
1183 
1184 /*@C
1185    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1186    compressed row) format.  For good matrix assembly performance the
1187    user should preallocate the matrix storage by setting the parameter nz
1188    (or nzz).  By setting these parameters accurately, performance can be
1189    increased by more than a factor of 50.
1190 
1191    Input Parameters:
1192 .  comm - MPI communicator, set to MPI_COMM_SELF
1193 .  bs - size of block
1194 .  m - number of rows
1195 .  n - number of columns
1196 .  nz - number of block nonzeros per block row (same for all rows)
1197 .  nzz - number of block nonzeros per block row or PETSC_NULL
1198          (possibly different for each row)
1199 
1200    Output Parameter:
1201 .  A - the matrix
1202 
1203    Notes:
1204    The block AIJ format is fully compatible with standard Fortran 77
1205    storage.  That is, the stored row and column indices can begin at
1206    either one (as in Fortran) or zero.  See the users' manual for details.
1207 
1208    Specify the preallocated storage with either nz or nnz (not both).
1209    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1210    allocation.  For additional details, see the users manual chapter on
1211    matrices and the file $(PETSC_DIR)/Performance.
1212 
1213 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues()
1214 @*/
1215 int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
1216 {
1217   Mat         B;
1218   Mat_SeqBAIJ *b;
1219   int         i,len,ierr,flg,mbs=m/bs,nbs=n/bs,bs2=bs*bs;
1220 
1221   if (mbs*bs!=m || nbs*bs!=n)
1222     SETERRQ(1,"MatCreateSeqBAIJ:Number rows, cols must be divisible by blocksize");
1223 
1224   *A = 0;
1225   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
1226   PLogObjectCreate(B);
1227   B->data = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
1228   PetscMemzero(b,sizeof(Mat_SeqBAIJ));
1229   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
1230   switch (bs) {
1231     case 1:
1232       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1233       break;
1234     case 2:
1235       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1236       break;
1237     case 3:
1238       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1239       break;
1240     case 4:
1241       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1242       break;
1243     case 5:
1244       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1245       break;
1246   }
1247   B->destroy          = MatDestroy_SeqBAIJ;
1248   B->view             = MatView_SeqBAIJ;
1249   B->factor           = 0;
1250   B->lupivotthreshold = 1.0;
1251   b->row              = 0;
1252   b->col              = 0;
1253   b->reallocs         = 0;
1254 
1255   b->m       = m; B->m = m; B->M = m;
1256   b->n       = n; B->n = n; B->N = n;
1257   b->mbs     = mbs;
1258   b->nbs     = nbs;
1259   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
1260   if (nnz == PETSC_NULL) {
1261     if (nz == PETSC_DEFAULT) nz = 5;
1262     else if (nz <= 0)        nz = 1;
1263     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
1264     nz = nz*mbs;
1265   }
1266   else {
1267     nz = 0;
1268     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
1269   }
1270 
1271   /* allocate the matrix space */
1272   len   = nz*sizeof(int) + nz*bs2*sizeof(Scalar) + (b->m+1)*sizeof(int);
1273   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
1274   PetscMemzero(b->a,nz*bs2*sizeof(Scalar));
1275   b->j  = (int *) (b->a + nz*bs2);
1276   PetscMemzero(b->j,nz*sizeof(int));
1277   b->i  = b->j + nz;
1278   b->singlemalloc = 1;
1279 
1280   b->i[0] = 0;
1281   for (i=1; i<mbs+1; i++) {
1282     b->i[i] = b->i[i-1] + b->imax[i-1];
1283   }
1284 
1285   /* b->ilen will count nonzeros in each block row so far. */
1286   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
1287   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
1288   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
1289 
1290   b->bs               = bs;
1291   b->mbs              = mbs;
1292   b->nz               = 0;
1293   b->maxnz            = nz;
1294   b->sorted           = 0;
1295   b->roworiented      = 1;
1296   b->nonew            = 0;
1297   b->diag             = 0;
1298   b->solve_work       = 0;
1299   b->mult_work        = 0;
1300   b->spptr            = 0;
1301   b->indexshift       = 0;
1302   *A = B;
1303   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
1304   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
1305   return 0;
1306 }
1307 
1308 int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
1309 {
1310   Mat         C;
1311   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
1312   int         i,len, mbs = a->mbs, bs = a->bs,nz = a->nz,bs2 = bs*bs;
1313 
1314   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
1315 
1316   *B = 0;
1317   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
1318   PLogObjectCreate(C);
1319   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
1320   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
1321   C->destroy    = MatDestroy_SeqBAIJ;
1322   C->view       = MatView_SeqBAIJ;
1323   C->factor     = A->factor;
1324   c->row        = 0;
1325   c->col        = 0;
1326   C->assembled  = PETSC_TRUE;
1327 
1328   c->m = C->m   = a->m;
1329   c->n = C->n   = a->n;
1330   C->M          = a->m;
1331   C->N          = a->n;
1332 
1333   c->bs         = a->bs;
1334   c->mbs        = a->mbs;
1335   c->nbs        = a->nbs;
1336   c->indexshift = 0;
1337 
1338   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
1339   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
1340   for ( i=0; i<mbs; i++ ) {
1341     c->imax[i] = a->imax[i];
1342     c->ilen[i] = a->ilen[i];
1343   }
1344 
1345   /* allocate the matrix space */
1346   c->singlemalloc = 1;
1347   len   = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(Scalar) + sizeof(int));
1348   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
1349   c->j  = (int *) (c->a + nz*bs2);
1350   c->i  = c->j + nz;
1351   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
1352   if (mbs > 0) {
1353     PetscMemcpy(c->j,a->j,nz*sizeof(int));
1354     if (cpvalues == COPY_VALUES) {
1355       PetscMemcpy(c->a,a->a,bs2*nz*sizeof(Scalar));
1356     }
1357   }
1358 
1359   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
1360   c->sorted      = a->sorted;
1361   c->roworiented = a->roworiented;
1362   c->nonew       = a->nonew;
1363 
1364   if (a->diag) {
1365     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
1366     PLogObjectMemory(C,(mbs+1)*sizeof(int));
1367     for ( i=0; i<mbs; i++ ) {
1368       c->diag[i] = a->diag[i];
1369     }
1370   }
1371   else c->diag          = 0;
1372   c->nz                 = a->nz;
1373   c->maxnz              = a->maxnz;
1374   c->solve_work         = 0;
1375   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
1376   c->mult_work          = 0;
1377   *B = C;
1378   return 0;
1379 }
1380 
1381 int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
1382 {
1383   Mat_SeqBAIJ  *a;
1384   Mat          B;
1385   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
1386   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
1387   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1388   int          *masked, nmask,tmp,ishift,bs2;
1389   Scalar       *aa;
1390   MPI_Comm     comm = ((PetscObject) viewer)->comm;
1391 
1392   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
1393   bs2  = bs*bs;
1394 
1395   MPI_Comm_size(comm,&size);
1396   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
1397   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1398   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
1399   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
1400   M = header[1]; N = header[2]; nz = header[3];
1401 
1402   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
1403 
1404   /*
1405      This code adds extra rows to make sure the number of rows is
1406     divisible by the blocksize
1407   */
1408   mbs        = M/bs;
1409   extra_rows = bs - M + bs*(mbs);
1410   if (extra_rows == bs) extra_rows = 0;
1411   else                  mbs++;
1412   if (extra_rows) {
1413     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize");
1414   }
1415 
1416   /* read in row lengths */
1417   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
1418   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
1419   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
1420 
1421   /* read in column indices */
1422   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
1423   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
1424   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
1425 
1426   /* loop over row lengths determining block row lengths */
1427   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1428   PetscMemzero(browlengths,mbs*sizeof(int));
1429   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
1430   PetscMemzero(mask,mbs*sizeof(int));
1431   masked = mask + mbs;
1432   rowcount = 0; nzcount = 0;
1433   for ( i=0; i<mbs; i++ ) {
1434     nmask = 0;
1435     for ( j=0; j<bs; j++ ) {
1436       kmax = rowlengths[rowcount];
1437       for ( k=0; k<kmax; k++ ) {
1438         tmp = jj[nzcount++]/bs;
1439         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1440       }
1441       rowcount++;
1442     }
1443     browlengths[i] += nmask;
1444     /* zero out the mask elements we set */
1445     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1446   }
1447 
1448   /* create our matrix */
1449   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
1450          CHKERRQ(ierr);
1451   B = *A;
1452   a = (Mat_SeqBAIJ *) B->data;
1453 
1454   /* set matrix "i" values */
1455   a->i[0] = 0;
1456   for ( i=1; i<= mbs; i++ ) {
1457     a->i[i]      = a->i[i-1] + browlengths[i-1];
1458     a->ilen[i-1] = browlengths[i-1];
1459   }
1460   a->indexshift = 0;
1461   a->nz         = 0;
1462   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
1463 
1464   /* read in nonzero values */
1465   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
1466   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
1467   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
1468 
1469   /* set "a" and "j" values into matrix */
1470   nzcount = 0; jcount = 0;
1471   for ( i=0; i<mbs; i++ ) {
1472     nzcountb = nzcount;
1473     nmask    = 0;
1474     for ( j=0; j<bs; j++ ) {
1475       kmax = rowlengths[i*bs+j];
1476       for ( k=0; k<kmax; k++ ) {
1477         tmp = jj[nzcount++]/bs;
1478 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1479       }
1480       rowcount++;
1481     }
1482     /* sort the masked values */
1483     PetscSortInt(nmask,masked);
1484 
1485     /* set "j" values into matrix */
1486     maskcount = 1;
1487     for ( j=0; j<nmask; j++ ) {
1488       a->j[jcount++]  = masked[j];
1489       mask[masked[j]] = maskcount++;
1490     }
1491     /* set "a" values into matrix */
1492     ishift = bs2*a->i[i];
1493     for ( j=0; j<bs; j++ ) {
1494       kmax = rowlengths[i*bs+j];
1495       for ( k=0; k<kmax; k++ ) {
1496         tmp    = jj[nzcountb]/bs ;
1497         block  = mask[tmp] - 1;
1498         point  = jj[nzcountb] - bs*tmp;
1499         idx    = ishift + bs2*block + j + bs*point;
1500         a->a[idx] = aa[nzcountb++];
1501       }
1502     }
1503     /* zero out the mask elements we set */
1504     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
1505   }
1506   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
1507 
1508   PetscFree(rowlengths);
1509   PetscFree(browlengths);
1510   PetscFree(aa);
1511   PetscFree(jj);
1512   PetscFree(mask);
1513 
1514   B->assembled = PETSC_TRUE;
1515 
1516   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
1517   if (flg) {
1518     Viewer tviewer;
1519     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
1520     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr);
1521     ierr = MatView(B,tviewer); CHKERRQ(ierr);
1522     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1523   }
1524   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
1525   if (flg) {
1526     Viewer tviewer;
1527     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
1528     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
1529     ierr = MatView(B,tviewer); CHKERRQ(ierr);
1530     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1531   }
1532   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
1533   if (flg) {
1534     Viewer tviewer;
1535     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
1536     ierr = MatView(B,tviewer); CHKERRQ(ierr);
1537     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1538   }
1539   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
1540   if (flg) {
1541     Viewer tviewer;
1542     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
1543     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr);
1544     ierr = MatView(B,tviewer); CHKERRQ(ierr);
1545     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1546   }
1547   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
1548   if (flg) {
1549     Viewer tviewer;
1550     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
1551     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
1552     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
1553     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
1554   }
1555   return 0;
1556 }
1557 
1558 
1559 
1560