xref: /petsc/src/mat/impls/baij/seq/baij.c (revision de6a44a3c49fb11edb0226bf70e60025dba29d78)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: baij.c,v 1.3 1996/02/13 23:29:47 bsmith Exp bsmith $";
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**);
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;
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   ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,&ia,&ja);CHKERRQ(ierr);
41   ierr = MatGetReordering_IJ(n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
42 
43   PetscFree(ia); PetscFree(ja);
44   return 0;
45 }
46 
47 /*
48      Adds diagonal pointers to sparse matrix structure.
49 */
50 
51 int MatMarkDiag_SeqBAIJ(Mat A)
52 {
53   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
54   int         i,j, *diag, m = a->m;
55 
56   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
57   PLogObjectMemory(A,(m+1)*sizeof(int));
58   for ( i=0; i<a->m; i++ ) {
59     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
60       if (a->j[j] == i) {
61         diag[i] = j;
62         break;
63       }
64     }
65   }
66   a->diag = diag;
67   return 0;
68 }
69 
70 #include "draw.h"
71 #include "pinclude/pviewer.h"
72 #include "sysio.h"
73 
74 static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
75 {
76   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
77   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l;
78   Scalar      *aa;
79 
80   ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr);
81   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
82   col_lens[0] = MAT_COOKIE;
83   col_lens[1] = a->m;
84   col_lens[2] = a->n;
85   col_lens[3] = a->nz*bs*bs;
86 
87   /* store lengths of each row and write (including header) to file */
88   count = 0;
89   for ( i=0; i<a->mbs; i++ ) {
90     for ( j=0; j<bs; j++ ) {
91       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
92     }
93   }
94   ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr);
95   PetscFree(col_lens);
96 
97   /* store column indices (zero start index) */
98   jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj);
99   count = 0;
100   for ( i=0; i<a->mbs; i++ ) {
101     for ( j=0; j<bs; j++ ) {
102       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
103         for ( l=0; l<bs; l++ ) {
104           jj[count++] = bs*a->j[k] + l;
105         }
106       }
107     }
108   }
109   ierr = SYWrite(fd,jj,bs*bs*a->nz,SYINT,0); CHKERRQ(ierr);
110   PetscFree(jj);
111 
112   /* store nonzero values */
113   aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa);
114   count = 0;
115   for ( i=0; i<a->mbs; i++ ) {
116     for ( j=0; j<bs; j++ ) {
117       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
118         for ( l=0; l<bs; l++ ) {
119           aa[count++] = a->a[bs*bs*k + l*bs + j];
120         }
121       }
122     }
123   }
124   ierr = SYWrite(fd,aa,bs*bs*a->nz,SYSCALAR,0); CHKERRQ(ierr);
125   PetscFree(aa);
126   return 0;
127 }
128 
129 static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
130 {
131   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
132   int         ierr, i,j,format,bs = a->bs,k,l;
133   FILE        *fd;
134   char        *outputname;
135 
136   ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr);
137   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
138   ierr = ViewerFileGetFormat_Private(viewer,&format);
139   if (format == FILE_FORMAT_INFO) {
140     /* no need to print additional information */ ;
141   }
142   else if (format == FILE_FORMAT_MATLAB) {
143     SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported");
144   }
145   else {
146     for ( i=0; i<a->mbs; i++ ) {
147       for ( j=0; j<bs; j++ ) {
148         fprintf(fd,"row %d:",i*bs+j);
149         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
150           for ( l=0; l<bs; l++ ) {
151             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs*bs*k + l*bs + j]);
152           }
153         }
154         fprintf(fd,"\n");
155       }
156     }
157   }
158   fflush(fd);
159   return 0;
160 }
161 
162 static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
163 {
164   Mat         A = (Mat) obj;
165   PetscObject vobj = (PetscObject) viewer;
166 
167   if (!viewer) {
168     viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer;
169   }
170   if (vobj->cookie == VIEWER_COOKIE) {
171     if (vobj->type == MATLAB_VIEWER) {
172       SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported");
173     }
174     else if (vobj->type == ASCII_FILE_VIEWER||vobj->type == ASCII_FILES_VIEWER){
175       return MatView_SeqBAIJ_ASCII(A,viewer);
176     }
177     else if (vobj->type == BINARY_FILE_VIEWER) {
178       return MatView_SeqBAIJ_Binary(A,viewer);
179     }
180   }
181   else if (vobj->cookie == DRAW_COOKIE) {
182     if (vobj->type == NULLWINDOW) return 0;
183     SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported");
184   }
185   return 0;
186 }
187 
188 
189 static int MatZeroEntries_SeqBAIJ(Mat A)
190 {
191   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
192   PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar));
193   return 0;
194 }
195 
196 int MatDestroy_SeqBAIJ(PetscObject obj)
197 {
198   Mat        A  = (Mat) obj;
199   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
200 
201 #if defined(PETSC_LOG)
202   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
203 #endif
204   PetscFree(a->a);
205   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
206   if (a->diag) PetscFree(a->diag);
207   if (a->ilen) PetscFree(a->ilen);
208   if (a->imax) PetscFree(a->imax);
209   if (a->solve_work) PetscFree(a->solve_work);
210   if (a->mult_work) PetscFree(a->mult_work);
211   PetscFree(a);
212   PLogObjectDestroy(A);
213   PetscHeaderDestroy(A);
214   return 0;
215 }
216 
217 static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
218 {
219   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
220   if      (op == ROW_ORIENTED)              a->roworiented = 1;
221   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
222   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
223   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
224   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
225   else if (op == ROWS_SORTED ||
226            op == SYMMETRIC_MATRIX ||
227            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
228            op == YES_NEW_DIAGONALS)
229     PLogInfo((PetscObject)A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
230   else if (op == NO_NEW_DIAGONALS)
231     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
232   else
233     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
234   return 0;
235 }
236 
237 
238 /* -------------------------------------------------------*/
239 /* Should check that shapes of vectors and matrices match */
240 /* -------------------------------------------------------*/
241 #include "pinclude/plapack.h"
242 
243 static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy)
244 {
245   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
246   Scalar          *xg,*yg;
247   register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5;
248   register Scalar x1,x2,x3,x4,x5;
249   int             mbs = a->mbs, m = a->m, i, *idx,*ii;
250   int             bs = a->bs,j,n,bs2 = bs*bs;
251 
252   VecGetArray(xx,&xg); x = xg;  VecGetArray(yy,&yg); y = yg;
253   PetscMemzero(y,m*sizeof(Scalar));
254   x     = x;
255   idx   = a->j;
256   v     = a->a;
257   ii    = a->i;
258 
259   switch (bs) {
260     case 1:
261      for ( i=0; i<m; i++ ) {
262        n    = ii[1] - ii[0]; ii++;
263        sum  = 0.0;
264        while (n--) sum += *v++ * x[*idx++];
265        y[i] = sum;
266       }
267       break;
268     case 2:
269       for ( i=0; i<mbs; i++ ) {
270         n  = ii[1] - ii[0]; ii++;
271         sum1 = 0.0; sum2 = 0.0;
272         for ( j=0; j<n; j++ ) {
273           xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
274           sum1 += v[0]*x1 + v[2]*x2;
275           sum2 += v[1]*x1 + v[3]*x2;
276           v += 4;
277         }
278         y[0] += sum1; y[1] += sum2;
279         y += 2;
280       }
281       break;
282     case 3:
283       for ( i=0; i<mbs; i++ ) {
284         n  = ii[1] - ii[0]; ii++;
285         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
286         for ( j=0; j<n; j++ ) {
287           xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
288           sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
289           sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
290           sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
291           v += 9;
292         }
293         y[0] += sum1; y[1] += sum2; y[2] += sum3;
294         y += 3;
295       }
296       break;
297     case 4:
298       for ( i=0; i<mbs; i++ ) {
299         n  = ii[1] - ii[0]; ii++;
300         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
301         for ( j=0; j<n; j++ ) {
302           xb = x + 4*(*idx++);
303           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
304           sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
305           sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
306           sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
307           sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
308           v += 16;
309         }
310         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4;
311         y += 4;
312       }
313       break;
314     case 5:
315       for ( i=0; i<mbs; i++ ) {
316         n  = ii[1] - ii[0]; ii++;
317         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
318         for ( j=0; j<n; j++ ) {
319           xb = x + 5*(*idx++);
320           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
321           sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
322           sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
323           sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
324           sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
325           sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
326           v += 25;
327         }
328         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5;
329         y += 5;
330       }
331       break;
332       /* block sizes larger then 5 by 5 are handled by BLAS */
333     default: {
334       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
335       if (!a->mult_work) {
336         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
337         CHKPTRQ(a->mult_work);
338       }
339       work = a->mult_work;
340       for ( i=0; i<mbs; i++ ) {
341         n     = ii[1] - ii[0]; ii++;
342         ncols = n*bs;
343         workt = work;
344         for ( j=0; j<n; j++ ) {
345           xb = x + bs*(*idx++);
346           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
347           workt += bs;
348         }
349         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One);
350         v += n*bs2;
351         y += bs;
352       }
353     }
354   }
355   PLogFlops(2*bs2*a->nz - m);
356   return 0;
357 }
358 
359 static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
360 {
361   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
362   *nz  = a->bs*a->bs*a->nz;
363   *nza = a->maxnz;
364   *mem = (int)A->mem;
365   return 0;
366 }
367 
368 extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
369 extern int MatLUFactorNumeric_SeqBAIJ(Mat,Mat*);
370 extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
371 extern int MatSolve_SeqBAIJ(Mat,Vec,Vec);
372 
373 static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
374 {
375   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
376   *m = a->m; *n = a->n;
377   return 0;
378 }
379 
380 static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
381 {
382   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
383   *m = 0; *n = a->m;
384   return 0;
385 }
386 
387 static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
388 {
389   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
390   Scalar      *v = a->a;
391   double      sum = 0.0;
392   int         i;
393 
394   if (type == NORM_FROBENIUS) {
395     for (i=0; i<a->nz; i++ ) {
396 #if defined(PETSC_COMPLEX)
397       sum += real(conj(*v)*(*v)); v++;
398 #else
399       sum += (*v)*(*v); v++;
400 #endif
401     }
402     *norm = sqrt(sum);
403   }
404   else {
405     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
406   }
407   return 0;
408 }
409 
410 /*
411      note: This can only work for identity for row and col. It would
412    be good to check this and otherwise generate an error.
413 */
414 static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
415 {
416   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
417   Mat         outA;
418   int         ierr;
419 
420   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
421 
422   outA          = inA;
423   inA->factor   = FACTOR_LU;
424   a->row        = row;
425   a->col        = col;
426 
427   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
428 
429   if (!a->diag) {
430     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
431   }
432   ierr = MatLUFactorNumeric_SeqBAIJ(inA,&outA); CHKERRQ(ierr);
433   return 0;
434 }
435 
436 static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
437 {
438   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
439   int         one = 1, totalnz = a->bs*a->bs*a->nz;
440   BLscal_( &totalnz, alpha, a->a, &one );
441   PLogFlops(totalnz);
442   return 0;
443 }
444 
445 int MatPrintHelp_SeqBAIJ(Mat A)
446 {
447   static int called = 0;
448   MPI_Comm   comm = A->comm;
449 
450   if (called) return 0; else called = 1;
451   return 0;
452 }
453 /* -------------------------------------------------------------------*/
454 static struct _MatOps MatOps = {0,
455        0,0,
456        MatMult_SeqBAIJ,0,
457        0,0,
458        MatSolve_SeqBAIJ,0,
459        0,0,
460        MatLUFactor_SeqBAIJ,0,
461        0,
462        0,
463        MatGetInfo_SeqBAIJ,0,
464        0,0,MatNorm_SeqBAIJ,
465        0,0,
466        0,
467        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
468        MatGetReordering_SeqBAIJ,
469        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ,0,0,
470        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
471        MatILUFactorSymbolic_SeqBAIJ,0,
472        0,0,/*MatConvert_SeqBAIJ*/ 0,
473        0,0,
474        MatConvertSameType_SeqBAIJ,0,0,
475        MatILUFactor_SeqBAIJ,0,0,
476        0,0,
477        0,0,
478        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
479        0};
480 
481 /*@C
482    MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format
483    (the default parallel PETSc format).  For good matrix assembly performance
484    the user should preallocate the matrix storage by setting the parameter nz
485    (or nzz).  By setting these parameters accurately, performance can be
486    increased by more than a factor of 50.
487 
488    Input Parameters:
489 .  comm - MPI communicator, set to MPI_COMM_SELF
490 .  bs - size of block
491 .  m - number of rows
492 .  n - number of columns
493 .  nz - number of block nonzeros per block row (same for all rows)
494 .  nzz - number of block nonzeros per block row or PETSC_NULL
495          (possibly different for each row)
496 
497    Output Parameter:
498 .  A - the matrix
499 
500    Notes:
501    The BAIJ format, is fully compatible with standard Fortran 77
502    storage.  That is, the stored row and column indices can begin at
503    either one (as in Fortran) or zero.  See the users manual for details.
504 
505    Specify the preallocated storage with either nz or nnz (not both).
506    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
507    allocation.  For additional details, see the users manual chapter on
508    matrices and the file $(PETSC_DIR)/Performance.
509 
510 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
511 @*/
512 int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
513 {
514   Mat         B;
515   Mat_SeqBAIJ *b;
516   int         i,len,ierr, flg,mbs = m/bs;
517 
518   if (mbs*bs != m)
519     SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize");
520 
521   *A      = 0;
522   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
523   PLogObjectCreate(B);
524   B->data             = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
525   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
526   B->destroy          = MatDestroy_SeqBAIJ;
527   B->view             = MatView_SeqBAIJ;
528   B->factor           = 0;
529   B->lupivotthreshold = 1.0;
530   b->row              = 0;
531   b->col              = 0;
532   b->reallocs         = 0;
533 
534   b->m       = m;
535   b->mbs     = mbs;
536   b->n       = n;
537   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
538   if (nnz == PETSC_NULL) {
539     if (nz == PETSC_DEFAULT) nz = 5;
540     else if (nz <= 0)        nz = 1;
541     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
542     nz = nz*mbs;
543   }
544   else {
545     nz = 0;
546     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
547   }
548 
549   /* allocate the matrix space */
550   len   = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int);
551   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
552   PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar));
553   b->j  = (int *) (b->a + nz*bs*bs);
554   PetscMemzero(b->j,nz*sizeof(int));
555   b->i  = b->j + nz;
556   b->singlemalloc = 1;
557 
558   b->i[0] = 0;
559   for (i=1; i<mbs+1; i++) {
560     b->i[i] = b->i[i-1] + b->imax[i-1];
561   }
562 
563   /* b->ilen will count nonzeros in each block row so far. */
564   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
565   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
566   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
567 
568   b->bs               = bs;
569   b->mbs              = mbs;
570   b->nz               = 0;
571   b->maxnz            = nz;
572   b->sorted           = 0;
573   b->roworiented      = 1;
574   b->nonew            = 0;
575   b->diag             = 0;
576   b->solve_work       = 0;
577   b->mult_work        = 0;
578   b->spptr            = 0;
579 
580   *A = B;
581   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
582   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
583   return 0;
584 }
585 
586 int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
587 {
588   Mat         C;
589   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
590   int         i,len, mbs = a->mbs, bs = a->bs,nz = a->nz;
591 
592   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
593 
594   *B = 0;
595   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
596   PLogObjectCreate(C);
597   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
598   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
599   C->destroy    = MatDestroy_SeqBAIJ;
600   C->view       = MatView_SeqBAIJ;
601   C->factor     = A->factor;
602   c->row        = 0;
603   c->col        = 0;
604   C->assembled  = PETSC_TRUE;
605 
606   c->m          = a->m;
607   c->n          = a->n;
608   c->bs         = a->bs;
609   c->mbs        = a->mbs;
610   c->nbs        = a->nbs;
611 
612   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
613   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
614   for ( i=0; i<mbs; i++ ) {
615     c->imax[i] = a->imax[i];
616     c->ilen[i] = a->ilen[i];
617   }
618 
619   /* allocate the matrix space */
620   c->singlemalloc = 1;
621   len   = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int));
622   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
623   c->j  = (int *) (c->a + nz*bs*bs);
624   c->i  = c->j + nz;
625   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
626   if (mbs > 0) {
627     PetscMemcpy(c->j,a->j,nz*sizeof(int));
628     if (cpvalues == COPY_VALUES) {
629       PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar));
630     }
631   }
632 
633   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
634   c->sorted      = a->sorted;
635   c->roworiented = a->roworiented;
636   c->nonew       = a->nonew;
637 
638   if (a->diag) {
639     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
640     PLogObjectMemory(C,(mbs+1)*sizeof(int));
641     for ( i=0; i<mbs; i++ ) {
642       c->diag[i] = a->diag[i];
643     }
644   }
645   else c->diag          = 0;
646   c->nz                 = a->nz;
647   c->maxnz              = a->maxnz;
648   c->solve_work         = 0;
649   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
650 
651   *B = C;
652   return 0;
653 }
654 
655 int MatLoad_SeqBAIJ(Viewer bview,MatType type,Mat *A)
656 {
657   Mat_SeqBAIJ  *a;
658   Mat          B;
659   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
660   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
661   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
662   int          *masked, nmask,tmp,ishift,bs2;
663   Scalar       *aa;
664   PetscObject  vobj = (PetscObject) bview;
665   MPI_Comm     comm = vobj->comm;
666 
667   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
668   bs2  = bs*bs;
669 
670   MPI_Comm_size(comm,&size);
671   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
672   ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr);
673   ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr);
674   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
675   M = header[1]; N = header[2]; nz = header[3];
676 
677   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
678 
679   /*
680      This code adds extra rows to make sure the number of rows is
681     divisible by the blocksize
682   */
683   mbs        = M/bs;
684   extra_rows = bs - M + bs*(mbs);
685   if (extra_rows == bs) extra_rows = 0;
686   else                  mbs++;
687   if (extra_rows) {
688     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize");
689   }
690 
691   /* read in row lengths */
692   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
693   ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr);
694   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
695 
696   /* read in column indices */
697   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
698   ierr = SYRead(fd,jj,nz,SYINT); CHKERRQ(ierr);
699   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
700 
701   /* loop over row lengths determining block row lengths */
702   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
703   PetscMemzero(browlengths,mbs*sizeof(int));
704   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
705   PetscMemzero(mask,mbs*sizeof(int));
706   masked = mask + mbs;
707   rowcount = 0; nzcount = 0;
708   for ( i=0; i<mbs; i++ ) {
709     nmask = 0;
710     for ( j=0; j<bs; j++ ) {
711       kmax = rowlengths[rowcount];
712       for ( k=0; k<kmax; k++ ) {
713         tmp = jj[nzcount++]/bs;
714         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
715       }
716       rowcount++;
717     }
718     browlengths[i] += nmask;
719     /* zero out the mask elements we set */
720     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
721   }
722 
723   /* create our matrix */
724   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
725          CHKERRQ(ierr);
726   B = *A;
727   a = (Mat_SeqBAIJ *) B->data;
728 
729   /* set matrix "i" values */
730   a->i[0] = 0;
731   for ( i=1; i<= mbs; i++ ) {
732     a->i[i]      = a->i[i-1] + browlengths[i-1];
733     a->ilen[i-1] = browlengths[i-1];
734   }
735   a->nz = 0;
736   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
737 
738   /* read in nonzero values */
739   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
740   ierr = SYRead(fd,aa,nz,SYSCALAR); CHKERRQ(ierr);
741   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
742 
743   /* set "a" and "j" values into matrix */
744   nzcount = 0; jcount = 0;
745   for ( i=0; i<mbs; i++ ) {
746     nzcountb = nzcount;
747     nmask    = 0;
748     for ( j=0; j<bs; j++ ) {
749       kmax = rowlengths[i*bs+j];
750       for ( k=0; k<kmax; k++ ) {
751         tmp = jj[nzcount++]/bs;
752 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
753       }
754       rowcount++;
755     }
756     /* sort the masked values */
757     SYIsort(nmask,masked);
758 
759     /* set "j" values into matrix */
760     maskcount = 1;
761     for ( j=0; j<nmask; j++ ) {
762       a->j[jcount++]  = masked[j];
763       mask[masked[j]] = maskcount++;
764     }
765     /* set "a" values into matrix */
766     ishift = bs2*a->i[i];
767     for ( j=0; j<bs; j++ ) {
768       kmax = rowlengths[i*bs+j];
769       for ( k=0; k<kmax; k++ ) {
770         tmp    = jj[nzcountb]/bs ;
771         block  = mask[tmp] - 1;
772         point  = jj[nzcountb] - bs*tmp;
773         idx    = ishift + bs2*block + j + bs*point;
774         a->a[idx] = aa[nzcountb++];
775       }
776     }
777     /* zero out the mask elements we set */
778     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
779   }
780   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
781 
782   PetscFree(rowlengths);
783   PetscFree(browlengths);
784   PetscFree(aa);
785   PetscFree(jj);
786   PetscFree(mask);
787 
788   B->assembled = PETSC_TRUE;
789 
790   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
791   if (flg) {
792     Viewer viewer;
793     ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr);
794     ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO,0);CHKERRQ(ierr);
795     ierr = MatView(B,viewer); CHKERRQ(ierr);
796     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
797   }
798   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
799   if (flg) {
800     Viewer viewer;
801     ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr);
802     ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
803     ierr = MatView(B,viewer); CHKERRQ(ierr);
804     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
805   }
806   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
807   if (flg) {
808     Viewer viewer;
809     ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr);
810     ierr = MatView(B,viewer); CHKERRQ(ierr);
811     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
812   }
813   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
814   if (flg) {
815     Viewer viewer;
816     ierr = ViewerFileOpenASCII(B->comm,"stdout",&viewer);CHKERRQ(ierr);
817     ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_MATLAB,"M");CHKERRQ(ierr);
818     ierr = MatView(B,viewer); CHKERRQ(ierr);
819     ierr = ViewerDestroy(viewer); CHKERRQ(ierr);
820   }
821   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
822   if (flg) {
823     Draw    win;
824     ierr = DrawOpenX(B->comm,0,0,0,0,300,300,&win); CHKERRQ(ierr);
825     ierr = MatView(B,(Viewer)win); CHKERRQ(ierr);
826     ierr = DrawSyncFlush(win); CHKERRQ(ierr);
827     ierr = DrawDestroy(win); CHKERRQ(ierr);
828   }
829   return 0;
830 }
831 
832 
833 
834