xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 91d316f63ec8ac6450179e21ca8a756379bd7a2e)
1b6490206SBarry Smith 
22593348eSBarry Smith #ifndef lint
3*91d316f6SSatish Balay static char vcid[] = "$Id: baij.c,v 1.16 1996/03/23 20:43:02 bsmith Exp balay $";
42593348eSBarry Smith #endif
52593348eSBarry Smith 
62593348eSBarry Smith /*
7b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
82593348eSBarry Smith   matrix storage format.
92593348eSBarry Smith */
10b6490206SBarry Smith #include "baij.h"
112593348eSBarry Smith #include "vec/vecimpl.h"
122593348eSBarry Smith #include "inline/spops.h"
132593348eSBarry Smith #include "petsc.h"
142593348eSBarry Smith 
15bcd2baecSBarry Smith extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
162593348eSBarry Smith 
17b6490206SBarry Smith static int MatGetReordering_SeqBAIJ(Mat A,MatOrdering type,IS *rperm,IS *cperm)
182593348eSBarry Smith {
19b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
20bcd2baecSBarry Smith   int         ierr, *ia, *ja,n = a->mbs,*idx,i,ishift,oshift;
212593348eSBarry Smith 
222593348eSBarry Smith   /*
232593348eSBarry Smith      this is tacky: In the future when we have written special factorization
242593348eSBarry Smith      and solve routines for the identity permutation we should use a
252593348eSBarry Smith      stride index set instead of the general one.
262593348eSBarry Smith   */
272593348eSBarry Smith   if (type  == ORDER_NATURAL) {
282593348eSBarry Smith     idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx);
292593348eSBarry Smith     for ( i=0; i<n; i++ ) idx[i] = i;
302593348eSBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr);
312593348eSBarry Smith     ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr);
322593348eSBarry Smith     PetscFree(idx);
332593348eSBarry Smith     ISSetPermutation(*rperm);
342593348eSBarry Smith     ISSetPermutation(*cperm);
352593348eSBarry Smith     ISSetIdentity(*rperm);
362593348eSBarry Smith     ISSetIdentity(*cperm);
372593348eSBarry Smith     return 0;
382593348eSBarry Smith   }
392593348eSBarry Smith 
40bcd2baecSBarry Smith   MatReorderingRegisterAll();
41bcd2baecSBarry Smith   ishift = a->indexshift;
42bcd2baecSBarry Smith   oshift = -MatReorderingIndexShift[(int)type];
43bcd2baecSBarry Smith   if (MatReorderingRequiresSymmetric[(int)type]) {
44bcd2baecSBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,ishift,oshift,&ia,&ja);
45bcd2baecSBarry Smith     CHKERRQ(ierr);
46bcd2baecSBarry Smith     ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr);
472593348eSBarry Smith     PetscFree(ia); PetscFree(ja);
48bcd2baecSBarry Smith   } else {
49bcd2baecSBarry Smith     if (ishift == oshift) {
50bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
51bcd2baecSBarry Smith     }
52bcd2baecSBarry Smith     else if (ishift == -1) {
53bcd2baecSBarry Smith       /* temporarily subtract 1 from i and j indices */
54bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
55bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
56bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
57bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
58bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
59bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
60bcd2baecSBarry Smith     } else {
61bcd2baecSBarry Smith       /* temporarily add 1 to i and j indices */
62bcd2baecSBarry Smith       int nz = a->i[a->n] - 1;
63bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]++;
64bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]++;
65bcd2baecSBarry Smith       ierr = MatGetReordering_IJ(a->n,a->i,a->j,type,rperm,cperm);CHKERRQ(ierr);
66bcd2baecSBarry Smith       for ( i=0; i<nz; i++ ) a->j[i]--;
67bcd2baecSBarry Smith       for ( i=0; i<a->n+1; i++ ) a->i[i]--;
68bcd2baecSBarry Smith     }
69bcd2baecSBarry Smith   }
702593348eSBarry Smith   return 0;
712593348eSBarry Smith }
722593348eSBarry Smith 
73de6a44a3SBarry Smith /*
74de6a44a3SBarry Smith      Adds diagonal pointers to sparse matrix structure.
75de6a44a3SBarry Smith */
76de6a44a3SBarry Smith 
77de6a44a3SBarry Smith int MatMarkDiag_SeqBAIJ(Mat A)
78de6a44a3SBarry Smith {
79de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
807fc0212eSBarry Smith   int         i,j, *diag, m = a->mbs;
81de6a44a3SBarry Smith 
82de6a44a3SBarry Smith   diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag);
83de6a44a3SBarry Smith   PLogObjectMemory(A,(m+1)*sizeof(int));
847fc0212eSBarry Smith   for ( i=0; i<m; i++ ) {
85de6a44a3SBarry Smith     for ( j=a->i[i]; j<a->i[i+1]; j++ ) {
86de6a44a3SBarry Smith       if (a->j[j] == i) {
87de6a44a3SBarry Smith         diag[i] = j;
88de6a44a3SBarry Smith         break;
89de6a44a3SBarry Smith       }
90de6a44a3SBarry Smith     }
91de6a44a3SBarry Smith   }
92de6a44a3SBarry Smith   a->diag = diag;
93de6a44a3SBarry Smith   return 0;
94de6a44a3SBarry Smith }
952593348eSBarry Smith 
962593348eSBarry Smith #include "draw.h"
972593348eSBarry Smith #include "pinclude/pviewer.h"
9877c4ece6SBarry Smith #include "sys.h"
992593348eSBarry Smith 
100b6490206SBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,Viewer viewer)
1012593348eSBarry Smith {
102b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
103b6490206SBarry Smith   int         i, fd, *col_lens, ierr, bs = a->bs,count,*jj,j,k,l;
104b6490206SBarry Smith   Scalar      *aa;
1052593348eSBarry Smith 
10690ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
1072593348eSBarry Smith   col_lens = (int *) PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
1082593348eSBarry Smith   col_lens[0] = MAT_COOKIE;
1092593348eSBarry Smith   col_lens[1] = a->m;
1102593348eSBarry Smith   col_lens[2] = a->n;
111b6490206SBarry Smith   col_lens[3] = a->nz*bs*bs;
1122593348eSBarry Smith 
1132593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
114b6490206SBarry Smith   count = 0;
115b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
116b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
117b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
118b6490206SBarry Smith     }
1192593348eSBarry Smith   }
12077c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,BINARY_INT,1); CHKERRQ(ierr);
1212593348eSBarry Smith   PetscFree(col_lens);
1222593348eSBarry Smith 
1232593348eSBarry Smith   /* store column indices (zero start index) */
124b6490206SBarry Smith   jj = (int *) PetscMalloc( a->nz*bs*bs*sizeof(int) ); CHKPTRQ(jj);
125b6490206SBarry Smith   count = 0;
126b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
127b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
128b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
129b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
130b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
1312593348eSBarry Smith         }
1322593348eSBarry Smith       }
133b6490206SBarry Smith     }
134b6490206SBarry Smith   }
13577c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs*bs*a->nz,BINARY_INT,0); CHKERRQ(ierr);
136b6490206SBarry Smith   PetscFree(jj);
1372593348eSBarry Smith 
1382593348eSBarry Smith   /* store nonzero values */
139b6490206SBarry Smith   aa = (Scalar *) PetscMalloc(a->nz*bs*bs*sizeof(Scalar)); CHKPTRQ(aa);
140b6490206SBarry Smith   count = 0;
141b6490206SBarry Smith   for ( i=0; i<a->mbs; i++ ) {
142b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
143b6490206SBarry Smith       for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
144b6490206SBarry Smith         for ( l=0; l<bs; l++ ) {
145b6490206SBarry Smith           aa[count++] = a->a[bs*bs*k + l*bs + j];
146b6490206SBarry Smith         }
147b6490206SBarry Smith       }
148b6490206SBarry Smith     }
149b6490206SBarry Smith   }
15077c4ece6SBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs*bs*a->nz,BINARY_SCALAR,0); CHKERRQ(ierr);
151b6490206SBarry Smith   PetscFree(aa);
1522593348eSBarry Smith   return 0;
1532593348eSBarry Smith }
1542593348eSBarry Smith 
155b6490206SBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,Viewer viewer)
1562593348eSBarry Smith {
157b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
158de6a44a3SBarry Smith   int         ierr, i,j,format,bs = a->bs,k,l;
1592593348eSBarry Smith   FILE        *fd;
1602593348eSBarry Smith   char        *outputname;
1612593348eSBarry Smith 
16290ace30eSBarry Smith   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
1632593348eSBarry Smith   ierr = ViewerFileGetOutputname_Private(viewer,&outputname);CHKERRQ(ierr);
16490ace30eSBarry Smith   ierr = ViewerGetFormat(viewer,&format);
16590ace30eSBarry Smith   if (format == ASCII_FORMAT_INFO) {
1662593348eSBarry Smith     /* no need to print additional information */ ;
1672593348eSBarry Smith   }
16890ace30eSBarry Smith   else if (format == ASCII_FORMAT_MATLAB) {
169b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ_ASCII:Matlab format not supported");
1702593348eSBarry Smith   }
1712593348eSBarry Smith   else {
172b6490206SBarry Smith     for ( i=0; i<a->mbs; i++ ) {
173b6490206SBarry Smith       for ( j=0; j<bs; j++ ) {
174b6490206SBarry Smith         fprintf(fd,"row %d:",i*bs+j);
175b6490206SBarry Smith         for ( k=a->i[i]; k<a->i[i+1]; k++ ) {
176b6490206SBarry Smith           for ( l=0; l<bs; l++ ) {
17788685aaeSLois Curfman McInnes #if defined(PETSC_COMPLEX)
17888685aaeSLois Curfman McInnes           if (imag(a->a[bs*bs*k + l*bs + j]) != 0.0) {
17988685aaeSLois Curfman McInnes             fprintf(fd," %d %g + %g i",bs*a->j[k]+l,
18088685aaeSLois Curfman McInnes               real(a->a[bs*bs*k + l*bs + j]),imag(a->a[bs*bs*k + l*bs + j]));
18188685aaeSLois Curfman McInnes           }
18288685aaeSLois Curfman McInnes           else {
18388685aaeSLois Curfman McInnes             fprintf(fd," %d %g",bs*a->j[k]+l,real(a->a[bs*bs*k + l*bs + j]));
18488685aaeSLois Curfman McInnes           }
18588685aaeSLois Curfman McInnes #else
186de6a44a3SBarry Smith             fprintf(fd," %d %g",bs*a->j[k]+l,a->a[bs*bs*k + l*bs + j]);
18788685aaeSLois Curfman McInnes #endif
1882593348eSBarry Smith           }
1892593348eSBarry Smith         }
1902593348eSBarry Smith         fprintf(fd,"\n");
1912593348eSBarry Smith       }
1922593348eSBarry Smith     }
193b6490206SBarry Smith   }
1942593348eSBarry Smith   fflush(fd);
1952593348eSBarry Smith   return 0;
1962593348eSBarry Smith }
1972593348eSBarry Smith 
198b6490206SBarry Smith static int MatView_SeqBAIJ(PetscObject obj,Viewer viewer)
1992593348eSBarry Smith {
2002593348eSBarry Smith   Mat         A = (Mat) obj;
20119bcc07fSBarry Smith   ViewerType  vtype;
20219bcc07fSBarry Smith   int         ierr;
2032593348eSBarry Smith 
2042593348eSBarry Smith   if (!viewer) {
20519bcc07fSBarry Smith     viewer = STDOUT_VIEWER_SELF;
2062593348eSBarry Smith   }
20719bcc07fSBarry Smith 
20819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
20919bcc07fSBarry Smith   if (vtype == MATLAB_VIEWER) {
210b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Matlab viewer not supported");
2112593348eSBarry Smith   }
21219bcc07fSBarry Smith   else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER){
213b6490206SBarry Smith     return MatView_SeqBAIJ_ASCII(A,viewer);
2142593348eSBarry Smith   }
21519bcc07fSBarry Smith   else if (vtype == BINARY_FILE_VIEWER) {
216b6490206SBarry Smith     return MatView_SeqBAIJ_Binary(A,viewer);
2172593348eSBarry Smith   }
21819bcc07fSBarry Smith   else if (vtype == DRAW_VIEWER) {
219b6490206SBarry Smith     SETERRQ(1,"MatView_SeqBAIJ:Draw viewer not supported");
2202593348eSBarry Smith   }
2212593348eSBarry Smith   return 0;
2222593348eSBarry Smith }
223b6490206SBarry Smith 
224b6490206SBarry Smith 
225b6490206SBarry Smith static int MatZeroEntries_SeqBAIJ(Mat A)
2262593348eSBarry Smith {
227b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
228de6a44a3SBarry Smith   PetscMemzero(a->a,a->bs*a->bs*a->i[a->mbs]*sizeof(Scalar));
2292593348eSBarry Smith   return 0;
2302593348eSBarry Smith }
2312593348eSBarry Smith 
232b6490206SBarry Smith int MatDestroy_SeqBAIJ(PetscObject obj)
2332593348eSBarry Smith {
2342593348eSBarry Smith   Mat        A  = (Mat) obj;
235b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2362593348eSBarry Smith 
2372593348eSBarry Smith #if defined(PETSC_LOG)
2382593348eSBarry Smith   PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz);
2392593348eSBarry Smith #endif
2402593348eSBarry Smith   PetscFree(a->a);
2412593348eSBarry Smith   if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);}
2422593348eSBarry Smith   if (a->diag) PetscFree(a->diag);
2432593348eSBarry Smith   if (a->ilen) PetscFree(a->ilen);
2442593348eSBarry Smith   if (a->imax) PetscFree(a->imax);
2452593348eSBarry Smith   if (a->solve_work) PetscFree(a->solve_work);
246de6a44a3SBarry Smith   if (a->mult_work) PetscFree(a->mult_work);
2472593348eSBarry Smith   PetscFree(a);
2482593348eSBarry Smith   PLogObjectDestroy(A);
2492593348eSBarry Smith   PetscHeaderDestroy(A);
2502593348eSBarry Smith   return 0;
2512593348eSBarry Smith }
2522593348eSBarry Smith 
253b6490206SBarry Smith static int MatSetOption_SeqBAIJ(Mat A,MatOption op)
2542593348eSBarry Smith {
255b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
2562593348eSBarry Smith   if      (op == ROW_ORIENTED)              a->roworiented = 1;
2572593348eSBarry Smith   else if (op == COLUMN_ORIENTED)           a->roworiented = 0;
2582593348eSBarry Smith   else if (op == COLUMNS_SORTED)            a->sorted      = 1;
2592593348eSBarry Smith   else if (op == NO_NEW_NONZERO_LOCATIONS)  a->nonew       = 1;
2602593348eSBarry Smith   else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew       = 0;
2612593348eSBarry Smith   else if (op == ROWS_SORTED ||
2622593348eSBarry Smith            op == SYMMETRIC_MATRIX ||
2632593348eSBarry Smith            op == STRUCTURALLY_SYMMETRIC_MATRIX ||
2642593348eSBarry Smith            op == YES_NEW_DIAGONALS)
26594a424c1SBarry Smith     PLogInfo(A,"Info:MatSetOption_SeqBAIJ:Option ignored\n");
2662593348eSBarry Smith   else if (op == NO_NEW_DIAGONALS)
267b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:NO_NEW_DIAGONALS");}
2682593348eSBarry Smith   else
269b6490206SBarry Smith     {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqBAIJ:unknown option");}
2702593348eSBarry Smith   return 0;
2712593348eSBarry Smith }
2722593348eSBarry Smith 
2732593348eSBarry Smith 
2742593348eSBarry Smith /* -------------------------------------------------------*/
2752593348eSBarry Smith /* Should check that shapes of vectors and matrices match */
2762593348eSBarry Smith /* -------------------------------------------------------*/
277b6490206SBarry Smith #include "pinclude/plapack.h"
278b6490206SBarry Smith 
279b6490206SBarry Smith static int MatMult_SeqBAIJ(Mat A,Vec xx,Vec yy)
2802593348eSBarry Smith {
281b6490206SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
282b6490206SBarry Smith   Scalar          *xg,*yg;
283b6490206SBarry Smith   register Scalar *x, *y, *v, sum,*xb, sum1,sum2,sum3,sum4,sum5;
284b6490206SBarry Smith   register Scalar x1,x2,x3,x4,x5;
285b6490206SBarry Smith   int             mbs = a->mbs, m = a->m, i, *idx,*ii;
286b6490206SBarry Smith   int             bs = a->bs,j,n,bs2 = bs*bs;
2872593348eSBarry Smith 
288b6490206SBarry Smith   VecGetArray(xx,&xg); x = xg;  VecGetArray(yy,&yg); y = yg;
289b6490206SBarry Smith   PetscMemzero(y,m*sizeof(Scalar));
290b6490206SBarry Smith   x     = x;
2912593348eSBarry Smith   idx   = a->j;
2922593348eSBarry Smith   v     = a->a;
2932593348eSBarry Smith   ii    = a->i;
294b6490206SBarry Smith 
295b6490206SBarry Smith   switch (bs) {
296b6490206SBarry Smith     case 1:
2972593348eSBarry Smith      for ( i=0; i<m; i++ ) {
2982593348eSBarry Smith        n    = ii[1] - ii[0]; ii++;
2992593348eSBarry Smith        sum  = 0.0;
3002593348eSBarry Smith        while (n--) sum += *v++ * x[*idx++];
3012593348eSBarry Smith        y[i] = sum;
3022593348eSBarry Smith       }
3032593348eSBarry Smith       break;
304b6490206SBarry Smith     case 2:
305b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
306de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
307b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0;
308b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
309b6490206SBarry Smith           xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
310b6490206SBarry Smith           sum1 += v[0]*x1 + v[2]*x2;
311b6490206SBarry Smith           sum2 += v[1]*x1 + v[3]*x2;
312b6490206SBarry Smith           v += 4;
313b6490206SBarry Smith         }
314b6490206SBarry Smith         y[0] += sum1; y[1] += sum2;
315b6490206SBarry Smith         y += 2;
316b6490206SBarry Smith       }
317b6490206SBarry Smith       break;
318b6490206SBarry Smith     case 3:
319b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
320de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
321b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
322b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
323b6490206SBarry Smith           xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
324b6490206SBarry Smith           sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
325b6490206SBarry Smith           sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
326b6490206SBarry Smith           sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
327b6490206SBarry Smith           v += 9;
328b6490206SBarry Smith         }
329b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3;
330b6490206SBarry Smith         y += 3;
331b6490206SBarry Smith       }
332b6490206SBarry Smith       break;
333b6490206SBarry Smith     case 4:
334b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
335de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
336b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
337b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
338b6490206SBarry Smith           xb = x + 4*(*idx++);
339b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
340b6490206SBarry Smith           sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
341b6490206SBarry Smith           sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
342b6490206SBarry Smith           sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
343b6490206SBarry Smith           sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
344b6490206SBarry Smith           v += 16;
345b6490206SBarry Smith         }
346b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4;
347b6490206SBarry Smith         y += 4;
348b6490206SBarry Smith       }
349b6490206SBarry Smith       break;
350b6490206SBarry Smith     case 5:
351b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
352de6a44a3SBarry Smith         n  = ii[1] - ii[0]; ii++;
353b6490206SBarry Smith         sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
354b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
355b6490206SBarry Smith           xb = x + 5*(*idx++);
356b6490206SBarry Smith           x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
357b6490206SBarry Smith           sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
358b6490206SBarry Smith           sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
359b6490206SBarry Smith           sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
360b6490206SBarry Smith           sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
361b6490206SBarry Smith           sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
362b6490206SBarry Smith           v += 25;
363b6490206SBarry Smith         }
364b6490206SBarry Smith         y[0] += sum1; y[1] += sum2; y[2] += sum3; y[3] += sum4; y[4] += sum5;
365b6490206SBarry Smith         y += 5;
366b6490206SBarry Smith       }
367b6490206SBarry Smith       break;
368b6490206SBarry Smith       /* block sizes larger then 5 by 5 are handled by BLAS */
369b6490206SBarry Smith     default: {
370de6a44a3SBarry Smith       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
371de6a44a3SBarry Smith       if (!a->mult_work) {
372de6a44a3SBarry Smith         a->mult_work = (Scalar *) PetscMalloc(a->m*sizeof(Scalar));
373de6a44a3SBarry Smith         CHKPTRQ(a->mult_work);
374de6a44a3SBarry Smith       }
375de6a44a3SBarry Smith       work = a->mult_work;
376b6490206SBarry Smith       for ( i=0; i<mbs; i++ ) {
377de6a44a3SBarry Smith         n     = ii[1] - ii[0]; ii++;
378de6a44a3SBarry Smith         ncols = n*bs;
379de6a44a3SBarry Smith         workt = work;
380b6490206SBarry Smith         for ( j=0; j<n; j++ ) {
381b6490206SBarry Smith           xb = x + bs*(*idx++);
382de6a44a3SBarry Smith           for ( k=0; k<bs; k++ ) workt[k] = xb[k];
383de6a44a3SBarry Smith           workt += bs;
384b6490206SBarry Smith         }
385de6a44a3SBarry Smith         LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,y,&_One);
386de6a44a3SBarry Smith         v += n*bs2;
387b6490206SBarry Smith         y += bs;
3882593348eSBarry Smith       }
3892593348eSBarry Smith     }
3902593348eSBarry Smith   }
391b6490206SBarry Smith   PLogFlops(2*bs2*a->nz - m);
3922593348eSBarry Smith   return 0;
3932593348eSBarry Smith }
3942593348eSBarry Smith 
395de6a44a3SBarry Smith static int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,int *nz,int *nza,int *mem)
3962593348eSBarry Smith {
397b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
398bcd2baecSBarry Smith   if (nz)  *nz  = a->bs*a->bs*a->nz;
399bcd2baecSBarry Smith   if (nza) *nza = a->maxnz;
400bcd2baecSBarry Smith   if (mem) *mem = (int)A->mem;
4012593348eSBarry Smith   return 0;
4022593348eSBarry Smith }
4032593348eSBarry Smith 
404*91d316f6SSatish Balay static int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
405*91d316f6SSatish Balay {
406*91d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
407*91d316f6SSatish Balay 
408*91d316f6SSatish Balay   if (B->type !=MATSEQBAIJ)SETERRQ(1,"MatEqual_SeqBAIJ:Matrices must be same type");
409*91d316f6SSatish Balay 
410*91d316f6SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
411*91d316f6SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)||
412*91d316f6SSatish Balay       (a->nz != b->nz)||(a->indexshift != b->indexshift)) {
413*91d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
414*91d316f6SSatish Balay   }
415*91d316f6SSatish Balay 
416*91d316f6SSatish Balay   /* if the a->i are the same */
417*91d316f6SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
418*91d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
419*91d316f6SSatish Balay   }
420*91d316f6SSatish Balay 
421*91d316f6SSatish Balay   /* if a->j are the same */
422*91d316f6SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
423*91d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
424*91d316f6SSatish Balay   }
425*91d316f6SSatish Balay 
426*91d316f6SSatish Balay   /* if a->a are the same */
427*91d316f6SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
428*91d316f6SSatish Balay     *flg = PETSC_FALSE; return 0;
429*91d316f6SSatish Balay   }
430*91d316f6SSatish Balay   *flg = PETSC_TRUE;
431*91d316f6SSatish Balay   return 0;
432*91d316f6SSatish Balay 
433*91d316f6SSatish Balay }
434*91d316f6SSatish Balay 
435*91d316f6SSatish Balay static int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
436*91d316f6SSatish Balay {
437*91d316f6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
438*91d316f6SSatish Balay   int        i,j, n,shift = a->indexshift;
439*91d316f6SSatish Balay   Scalar     *x, zero = 0.0;
440*91d316f6SSatish Balay 
441*91d316f6SSatish Balay   VecSet(&zero,v);
442*91d316f6SSatish Balay   VecGetArray(v,&x); VecGetLocalSize(v,&n);
443*91d316f6SSatish Balay   if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector");
444*91d316f6SSatish Balay   for ( i=0; i<a->mbs; i++ ) {
445*91d316f6SSatish Balay     for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) {
446*91d316f6SSatish Balay       if (a->j[j]+shift == i) {
447*91d316f6SSatish Balay         x[i] = a->a[j];
448*91d316f6SSatish Balay         break;
449*91d316f6SSatish Balay       }
450*91d316f6SSatish Balay     }
451*91d316f6SSatish Balay   }
452*91d316f6SSatish Balay   return 0;
453*91d316f6SSatish Balay }
454*91d316f6SSatish Balay 
455b6490206SBarry Smith extern int MatLUFactorSymbolic_SeqBAIJ(Mat,IS,IS,double,Mat*);
456b6490206SBarry Smith extern int MatLUFactor_SeqBAIJ(Mat,IS,IS,double);
457b6490206SBarry Smith extern int MatSolve_SeqBAIJ(Mat,Vec,Vec);
4587fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_N(Mat,Mat*);
4597fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_1(Mat,Mat*);
4607fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_2(Mat,Mat*);
4617fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_3(Mat,Mat*);
4627fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_4(Mat,Mat*);
4637fc0212eSBarry Smith extern int MatLUFactorNumeric_SeqBAIJ_5(Mat,Mat*);
4642593348eSBarry Smith 
465b6490206SBarry Smith static int MatGetSize_SeqBAIJ(Mat A,int *m,int *n)
4662593348eSBarry Smith {
467b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4682593348eSBarry Smith   *m = a->m; *n = a->n;
4692593348eSBarry Smith   return 0;
4702593348eSBarry Smith }
4712593348eSBarry Smith 
472b6490206SBarry Smith static int MatGetOwnershipRange_SeqBAIJ(Mat A,int *m,int *n)
4732593348eSBarry Smith {
474b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4752593348eSBarry Smith   *m = 0; *n = a->m;
4762593348eSBarry Smith   return 0;
4772593348eSBarry Smith }
478b6490206SBarry Smith 
479b6490206SBarry Smith static int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
4802593348eSBarry Smith {
481b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
4822593348eSBarry Smith   Scalar      *v = a->a;
4832593348eSBarry Smith   double      sum = 0.0;
484b6490206SBarry Smith   int         i;
4852593348eSBarry Smith 
4862593348eSBarry Smith   if (type == NORM_FROBENIUS) {
4872593348eSBarry Smith     for (i=0; i<a->nz; i++ ) {
4882593348eSBarry Smith #if defined(PETSC_COMPLEX)
4892593348eSBarry Smith       sum += real(conj(*v)*(*v)); v++;
4902593348eSBarry Smith #else
4912593348eSBarry Smith       sum += (*v)*(*v); v++;
4922593348eSBarry Smith #endif
4932593348eSBarry Smith     }
4942593348eSBarry Smith     *norm = sqrt(sum);
4952593348eSBarry Smith   }
4962593348eSBarry Smith   else {
497b6490206SBarry Smith     SETERRQ(1,"MatNorm_SeqBAIJ:No support for this norm yet");
4982593348eSBarry Smith   }
4992593348eSBarry Smith   return 0;
5002593348eSBarry Smith }
5012593348eSBarry Smith 
5022593348eSBarry Smith /*
5032593348eSBarry Smith      note: This can only work for identity for row and col. It would
5042593348eSBarry Smith    be good to check this and otherwise generate an error.
5052593348eSBarry Smith */
506b6490206SBarry Smith static int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,double efill,int fill)
5072593348eSBarry Smith {
508b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
5092593348eSBarry Smith   Mat         outA;
510de6a44a3SBarry Smith   int         ierr;
5112593348eSBarry Smith 
512b6490206SBarry Smith   if (fill != 0) SETERRQ(1,"MatILUFactor_SeqBAIJ:Only fill=0 supported");
5132593348eSBarry Smith 
5142593348eSBarry Smith   outA          = inA;
5152593348eSBarry Smith   inA->factor   = FACTOR_LU;
5162593348eSBarry Smith   a->row        = row;
5172593348eSBarry Smith   a->col        = col;
5182593348eSBarry Smith 
519de6a44a3SBarry Smith   a->solve_work = (Scalar *) PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
5202593348eSBarry Smith 
5212593348eSBarry Smith   if (!a->diag) {
522b6490206SBarry Smith     ierr = MatMarkDiag_SeqBAIJ(inA); CHKERRQ(ierr);
5232593348eSBarry Smith   }
5247fc0212eSBarry Smith   ierr = MatLUFactorNumeric(inA,&outA); CHKERRQ(ierr);
5252593348eSBarry Smith   return 0;
5262593348eSBarry Smith }
5272593348eSBarry Smith 
528b6490206SBarry Smith static int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
5292593348eSBarry Smith {
530b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
531b6490206SBarry Smith   int         one = 1, totalnz = a->bs*a->bs*a->nz;
532b6490206SBarry Smith   BLscal_( &totalnz, alpha, a->a, &one );
533b6490206SBarry Smith   PLogFlops(totalnz);
5342593348eSBarry Smith   return 0;
5352593348eSBarry Smith }
5362593348eSBarry Smith 
537b6490206SBarry Smith int MatPrintHelp_SeqBAIJ(Mat A)
5382593348eSBarry Smith {
5392593348eSBarry Smith   static int called = 0;
5402593348eSBarry Smith 
5412593348eSBarry Smith   if (called) return 0; else called = 1;
5422593348eSBarry Smith   return 0;
5432593348eSBarry Smith }
5442593348eSBarry Smith /* -------------------------------------------------------------------*/
545b6490206SBarry Smith static struct _MatOps MatOps = {0,
546b6490206SBarry Smith        0,0,
547b6490206SBarry Smith        MatMult_SeqBAIJ,0,
548b6490206SBarry Smith        0,0,
549de6a44a3SBarry Smith        MatSolve_SeqBAIJ,0,
550b6490206SBarry Smith        0,0,
551de6a44a3SBarry Smith        MatLUFactor_SeqBAIJ,0,
552b6490206SBarry Smith        0,
553b6490206SBarry Smith        0,
554b6490206SBarry Smith        MatGetInfo_SeqBAIJ,0,
555b6490206SBarry Smith        0,0,MatNorm_SeqBAIJ,
556b6490206SBarry Smith        0,0,
557b6490206SBarry Smith        0,
558b6490206SBarry Smith        MatSetOption_SeqBAIJ,MatZeroEntries_SeqBAIJ,0,
559b6490206SBarry Smith        MatGetReordering_SeqBAIJ,
5607fc0212eSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,MatLUFactorNumeric_SeqBAIJ_N,0,0,
561b6490206SBarry Smith        MatGetSize_SeqBAIJ,MatGetSize_SeqBAIJ,MatGetOwnershipRange_SeqBAIJ,
562de6a44a3SBarry Smith        MatILUFactorSymbolic_SeqBAIJ,0,
563b6490206SBarry Smith        0,0,/*  MatConvert_SeqBAIJ  */ 0,
564b6490206SBarry Smith        0,0,
565b6490206SBarry Smith        MatConvertSameType_SeqBAIJ,0,0,
566b6490206SBarry Smith        MatILUFactor_SeqBAIJ,0,0,
567b6490206SBarry Smith        0,0,
568b6490206SBarry Smith        0,0,
569b6490206SBarry Smith        MatPrintHelp_SeqBAIJ,MatScale_SeqBAIJ,
570b6490206SBarry Smith        0};
5712593348eSBarry Smith 
5722593348eSBarry Smith /*@C
573b6490206SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in AIJ (compressed row) format
5742593348eSBarry Smith    (the default parallel PETSc format).  For good matrix assembly performance
5752593348eSBarry Smith    the user should preallocate the matrix storage by setting the parameter nz
5762593348eSBarry Smith    (or nzz).  By setting these parameters accurately, performance can be
5772593348eSBarry Smith    increased by more than a factor of 50.
5782593348eSBarry Smith 
5792593348eSBarry Smith    Input Parameters:
5802593348eSBarry Smith .  comm - MPI communicator, set to MPI_COMM_SELF
581b6490206SBarry Smith .  bs - size of block
5822593348eSBarry Smith .  m - number of rows
5832593348eSBarry Smith .  n - number of columns
584b6490206SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
585b6490206SBarry Smith .  nzz - number of block nonzeros per block row or PETSC_NULL
586b6490206SBarry Smith          (possibly different for each row)
5872593348eSBarry Smith 
5882593348eSBarry Smith    Output Parameter:
5892593348eSBarry Smith .  A - the matrix
5902593348eSBarry Smith 
5912593348eSBarry Smith    Notes:
592b6490206SBarry Smith    The BAIJ format, is fully compatible with standard Fortran 77
5932593348eSBarry Smith    storage.  That is, the stored row and column indices can begin at
5942593348eSBarry Smith    either one (as in Fortran) or zero.  See the users manual for details.
5952593348eSBarry Smith 
5962593348eSBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
5972593348eSBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
5982593348eSBarry Smith    allocation.  For additional details, see the users manual chapter on
5992593348eSBarry Smith    matrices and the file $(PETSC_DIR)/Performance.
6002593348eSBarry Smith 
6012593348eSBarry Smith .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues()
6022593348eSBarry Smith @*/
603b6490206SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz, Mat *A)
6042593348eSBarry Smith {
6052593348eSBarry Smith   Mat         B;
606b6490206SBarry Smith   Mat_SeqBAIJ *b;
607b6490206SBarry Smith   int         i,len,ierr, flg,mbs = m/bs;
608b6490206SBarry Smith 
609b6490206SBarry Smith   if (mbs*bs != m)
610b6490206SBarry Smith     SETERRQ(1,"MatCreateSeqBAIJ:Number rows must be divisible by blocksize");
6112593348eSBarry Smith 
6122593348eSBarry Smith   *A      = 0;
613b6490206SBarry Smith   PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQBAIJ,comm);
6142593348eSBarry Smith   PLogObjectCreate(B);
615b6490206SBarry Smith   B->data             = (void *) (b = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(b);
6162593348eSBarry Smith   PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps));
6177fc0212eSBarry Smith   switch (bs) {
6187fc0212eSBarry Smith     case 1:
6197fc0212eSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
6207fc0212eSBarry Smith       break;
6214eeb42bcSBarry Smith     case 2:
6224eeb42bcSBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
6236d84be18SBarry Smith       break;
624f327dd97SBarry Smith     case 3:
625f327dd97SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
6264eeb42bcSBarry Smith       break;
6276d84be18SBarry Smith     case 4:
6286d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
6296d84be18SBarry Smith       break;
6306d84be18SBarry Smith     case 5:
6316d84be18SBarry Smith       B->ops.lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
6326d84be18SBarry Smith       break;
6337fc0212eSBarry Smith   }
634b6490206SBarry Smith   B->destroy          = MatDestroy_SeqBAIJ;
635b6490206SBarry Smith   B->view             = MatView_SeqBAIJ;
6362593348eSBarry Smith   B->factor           = 0;
6372593348eSBarry Smith   B->lupivotthreshold = 1.0;
6382593348eSBarry Smith   b->row              = 0;
6392593348eSBarry Smith   b->col              = 0;
6402593348eSBarry Smith   b->reallocs         = 0;
6412593348eSBarry Smith 
6422593348eSBarry Smith   b->m       = m;
643b6490206SBarry Smith   b->mbs     = mbs;
6442593348eSBarry Smith   b->n       = n;
645b6490206SBarry Smith   b->imax    = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(b->imax);
6462593348eSBarry Smith   if (nnz == PETSC_NULL) {
647de6a44a3SBarry Smith     if (nz == PETSC_DEFAULT) nz = 5;
6482593348eSBarry Smith     else if (nz <= 0)        nz = 1;
649b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) b->imax[i] = nz;
650b6490206SBarry Smith     nz = nz*mbs;
6512593348eSBarry Smith   }
6522593348eSBarry Smith   else {
6532593348eSBarry Smith     nz = 0;
654b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];}
6552593348eSBarry Smith   }
6562593348eSBarry Smith 
6572593348eSBarry Smith   /* allocate the matrix space */
658b6490206SBarry Smith   len   = nz*sizeof(int) + nz*bs*bs*sizeof(Scalar) + (b->m+1)*sizeof(int);
6592593348eSBarry Smith   b->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a);
660b6490206SBarry Smith   PetscMemzero(b->a,nz*bs*bs*sizeof(Scalar));
661b6490206SBarry Smith   b->j  = (int *) (b->a + nz*bs*bs);
6622593348eSBarry Smith   PetscMemzero(b->j,nz*sizeof(int));
6632593348eSBarry Smith   b->i  = b->j + nz;
6642593348eSBarry Smith   b->singlemalloc = 1;
6652593348eSBarry Smith 
666de6a44a3SBarry Smith   b->i[0] = 0;
667b6490206SBarry Smith   for (i=1; i<mbs+1; i++) {
6682593348eSBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
6692593348eSBarry Smith   }
6702593348eSBarry Smith 
671b6490206SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
672b6490206SBarry Smith   b->ilen = (int *) PetscMalloc((mbs+1)*sizeof(int));
673b6490206SBarry Smith   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
674b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) { b->ilen[i] = 0;}
6752593348eSBarry Smith 
676b6490206SBarry Smith   b->bs               = bs;
677b6490206SBarry Smith   b->mbs              = mbs;
6782593348eSBarry Smith   b->nz               = 0;
6792593348eSBarry Smith   b->maxnz            = nz;
6802593348eSBarry Smith   b->sorted           = 0;
6812593348eSBarry Smith   b->roworiented      = 1;
6822593348eSBarry Smith   b->nonew            = 0;
6832593348eSBarry Smith   b->diag             = 0;
6842593348eSBarry Smith   b->solve_work       = 0;
685de6a44a3SBarry Smith   b->mult_work        = 0;
6862593348eSBarry Smith   b->spptr            = 0;
6872593348eSBarry Smith 
6882593348eSBarry Smith   *A = B;
6892593348eSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr);
6902593348eSBarry Smith   if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); }
6912593348eSBarry Smith   return 0;
6922593348eSBarry Smith }
6932593348eSBarry Smith 
694b6490206SBarry Smith int MatConvertSameType_SeqBAIJ(Mat A,Mat *B,int cpvalues)
6952593348eSBarry Smith {
6962593348eSBarry Smith   Mat         C;
697b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ *) A->data;
698de6a44a3SBarry Smith   int         i,len, mbs = a->mbs, bs = a->bs,nz = a->nz;
699de6a44a3SBarry Smith 
700de6a44a3SBarry Smith   if (a->i[mbs] != nz) SETERRQ(1,"MatConvertSameType_SeqBAIJ:Corrupt matrix");
7012593348eSBarry Smith 
7022593348eSBarry Smith   *B = 0;
703b6490206SBarry Smith   PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQBAIJ,A->comm);
7042593348eSBarry Smith   PLogObjectCreate(C);
705b6490206SBarry Smith   C->data       = (void *) (c = PetscNew(Mat_SeqBAIJ)); CHKPTRQ(c);
7062593348eSBarry Smith   PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps));
707b6490206SBarry Smith   C->destroy    = MatDestroy_SeqBAIJ;
708b6490206SBarry Smith   C->view       = MatView_SeqBAIJ;
7092593348eSBarry Smith   C->factor     = A->factor;
7102593348eSBarry Smith   c->row        = 0;
7112593348eSBarry Smith   c->col        = 0;
7122593348eSBarry Smith   C->assembled  = PETSC_TRUE;
7132593348eSBarry Smith 
7142593348eSBarry Smith   c->m          = a->m;
7152593348eSBarry Smith   c->n          = a->n;
716b6490206SBarry Smith   c->bs         = a->bs;
717b6490206SBarry Smith   c->mbs        = a->mbs;
718b6490206SBarry Smith   c->nbs        = a->nbs;
7192593348eSBarry Smith 
720b6490206SBarry Smith   c->imax       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->imax);
721b6490206SBarry Smith   c->ilen       = (int *) PetscMalloc((mbs+1)*sizeof(int)); CHKPTRQ(c->ilen);
722b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
7232593348eSBarry Smith     c->imax[i] = a->imax[i];
7242593348eSBarry Smith     c->ilen[i] = a->ilen[i];
7252593348eSBarry Smith   }
7262593348eSBarry Smith 
7272593348eSBarry Smith   /* allocate the matrix space */
7282593348eSBarry Smith   c->singlemalloc = 1;
729de6a44a3SBarry Smith   len   = (mbs+1)*sizeof(int) + nz*(bs*bs*sizeof(Scalar) + sizeof(int));
7302593348eSBarry Smith   c->a  = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a);
731de6a44a3SBarry Smith   c->j  = (int *) (c->a + nz*bs*bs);
732de6a44a3SBarry Smith   c->i  = c->j + nz;
733b6490206SBarry Smith   PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));
734b6490206SBarry Smith   if (mbs > 0) {
735de6a44a3SBarry Smith     PetscMemcpy(c->j,a->j,nz*sizeof(int));
7362593348eSBarry Smith     if (cpvalues == COPY_VALUES) {
737de6a44a3SBarry Smith       PetscMemcpy(c->a,a->a,bs*bs*nz*sizeof(Scalar));
7382593348eSBarry Smith     }
7392593348eSBarry Smith   }
7402593348eSBarry Smith 
741b6490206SBarry Smith   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqBAIJ));
7422593348eSBarry Smith   c->sorted      = a->sorted;
7432593348eSBarry Smith   c->roworiented = a->roworiented;
7442593348eSBarry Smith   c->nonew       = a->nonew;
7452593348eSBarry Smith 
7462593348eSBarry Smith   if (a->diag) {
747b6490206SBarry Smith     c->diag = (int *) PetscMalloc( (mbs+1)*sizeof(int) ); CHKPTRQ(c->diag);
748b6490206SBarry Smith     PLogObjectMemory(C,(mbs+1)*sizeof(int));
749b6490206SBarry Smith     for ( i=0; i<mbs; i++ ) {
7502593348eSBarry Smith       c->diag[i] = a->diag[i];
7512593348eSBarry Smith     }
7522593348eSBarry Smith   }
7532593348eSBarry Smith   else c->diag          = 0;
7542593348eSBarry Smith   c->nz                 = a->nz;
7552593348eSBarry Smith   c->maxnz              = a->maxnz;
7562593348eSBarry Smith   c->solve_work         = 0;
7572593348eSBarry Smith   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
7587fc0212eSBarry Smith   c->mult_work          = 0;
7592593348eSBarry Smith   *B = C;
7602593348eSBarry Smith   return 0;
7612593348eSBarry Smith }
7622593348eSBarry Smith 
76319bcc07fSBarry Smith int MatLoad_SeqBAIJ(Viewer viewer,MatType type,Mat *A)
7642593348eSBarry Smith {
765b6490206SBarry Smith   Mat_SeqBAIJ  *a;
7662593348eSBarry Smith   Mat          B;
767de6a44a3SBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1,flg;
768b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
76935aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
770de6a44a3SBarry Smith   int          *masked, nmask,tmp,ishift,bs2;
771b6490206SBarry Smith   Scalar       *aa;
77219bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject) viewer)->comm;
7732593348eSBarry Smith 
77435aab85fSBarry Smith   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,&flg);CHKERRQ(ierr);
775de6a44a3SBarry Smith   bs2  = bs*bs;
776b6490206SBarry Smith 
7772593348eSBarry Smith   MPI_Comm_size(comm,&size);
778b6490206SBarry Smith   if (size > 1) SETERRQ(1,"MatLoad_SeqBAIJ:view must have one processor");
77990ace30eSBarry Smith   ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr);
78077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr);
781de6a44a3SBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqBAIJ:not Mat object");
7822593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
7832593348eSBarry Smith 
784b6490206SBarry Smith   if (M != N) SETERRQ(1,"MatLoad_SeqBAIJ:Can only do square matrices");
78535aab85fSBarry Smith 
78635aab85fSBarry Smith   /*
78735aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
78835aab85fSBarry Smith     divisible by the blocksize
78935aab85fSBarry Smith   */
790b6490206SBarry Smith   mbs        = M/bs;
79135aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
79235aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
79335aab85fSBarry Smith   else                  mbs++;
79435aab85fSBarry Smith   if (extra_rows) {
79535aab85fSBarry Smith     PLogInfo(0,"MatLoad_SeqBAIJ:Padding loading matrix to match blocksize");
79635aab85fSBarry Smith   }
797b6490206SBarry Smith 
7982593348eSBarry Smith   /* read in row lengths */
79935aab85fSBarry Smith   rowlengths = (int*) PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
80077c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr);
80135aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) rowlengths[M+i] = 1;
8022593348eSBarry Smith 
803b6490206SBarry Smith   /* read in column indices */
80435aab85fSBarry Smith   jj = (int*) PetscMalloc( (nz+extra_rows)*sizeof(int) ); CHKPTRQ(jj);
80577c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,BINARY_INT); CHKERRQ(ierr);
80635aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) jj[nz+i] = M+i;
807b6490206SBarry Smith 
808b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
809b6490206SBarry Smith   browlengths = (int *) PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
810b6490206SBarry Smith   PetscMemzero(browlengths,mbs*sizeof(int));
81135aab85fSBarry Smith   mask   = (int *) PetscMalloc( 2*mbs*sizeof(int) ); CHKPTRQ(mask);
81235aab85fSBarry Smith   PetscMemzero(mask,mbs*sizeof(int));
81335aab85fSBarry Smith   masked = mask + mbs;
814b6490206SBarry Smith   rowcount = 0; nzcount = 0;
815b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
81635aab85fSBarry Smith     nmask = 0;
817b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
818b6490206SBarry Smith       kmax = rowlengths[rowcount];
819b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
82035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
82135aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
822b6490206SBarry Smith       }
823b6490206SBarry Smith       rowcount++;
824b6490206SBarry Smith     }
82535aab85fSBarry Smith     browlengths[i] += nmask;
82635aab85fSBarry Smith     /* zero out the mask elements we set */
82735aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
828b6490206SBarry Smith   }
829b6490206SBarry Smith 
8302593348eSBarry Smith   /* create our matrix */
83135aab85fSBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);
83235aab85fSBarry Smith          CHKERRQ(ierr);
8332593348eSBarry Smith   B = *A;
834b6490206SBarry Smith   a = (Mat_SeqBAIJ *) B->data;
8352593348eSBarry Smith 
8362593348eSBarry Smith   /* set matrix "i" values */
837de6a44a3SBarry Smith   a->i[0] = 0;
838b6490206SBarry Smith   for ( i=1; i<= mbs; i++ ) {
839b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
840b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
8412593348eSBarry Smith   }
842b6490206SBarry Smith   a->nz = 0;
843b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) a->nz += browlengths[i];
8442593348eSBarry Smith 
845b6490206SBarry Smith   /* read in nonzero values */
84635aab85fSBarry Smith   aa = (Scalar *) PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
84777c4ece6SBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,BINARY_SCALAR); CHKERRQ(ierr);
84835aab85fSBarry Smith   for ( i=0; i<extra_rows; i++ ) aa[nz+i] = 1.0;
849b6490206SBarry Smith 
850b6490206SBarry Smith   /* set "a" and "j" values into matrix */
851b6490206SBarry Smith   nzcount = 0; jcount = 0;
852b6490206SBarry Smith   for ( i=0; i<mbs; i++ ) {
853b6490206SBarry Smith     nzcountb = nzcount;
85435aab85fSBarry Smith     nmask    = 0;
855b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
856b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
857b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
85835aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
85935aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
860b6490206SBarry Smith       }
861b6490206SBarry Smith       rowcount++;
862b6490206SBarry Smith     }
863de6a44a3SBarry Smith     /* sort the masked values */
86477c4ece6SBarry Smith     PetscSortInt(nmask,masked);
865de6a44a3SBarry Smith 
866b6490206SBarry Smith     /* set "j" values into matrix */
867b6490206SBarry Smith     maskcount = 1;
86835aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) {
86935aab85fSBarry Smith       a->j[jcount++]  = masked[j];
870de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
871b6490206SBarry Smith     }
872b6490206SBarry Smith     /* set "a" values into matrix */
873de6a44a3SBarry Smith     ishift = bs2*a->i[i];
874b6490206SBarry Smith     for ( j=0; j<bs; j++ ) {
875b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
876b6490206SBarry Smith       for ( k=0; k<kmax; k++ ) {
877de6a44a3SBarry Smith         tmp    = jj[nzcountb]/bs ;
878de6a44a3SBarry Smith         block  = mask[tmp] - 1;
879de6a44a3SBarry Smith         point  = jj[nzcountb] - bs*tmp;
880de6a44a3SBarry Smith         idx    = ishift + bs2*block + j + bs*point;
881b6490206SBarry Smith         a->a[idx] = aa[nzcountb++];
882b6490206SBarry Smith       }
883b6490206SBarry Smith     }
88435aab85fSBarry Smith     /* zero out the mask elements we set */
88535aab85fSBarry Smith     for ( j=0; j<nmask; j++ ) mask[masked[j]] = 0;
886b6490206SBarry Smith   }
88735aab85fSBarry Smith   if (jcount != a->nz) SETERRQ(1,"MatLoad_SeqBAIJ:Error bad binary matrix");
888b6490206SBarry Smith 
889b6490206SBarry Smith   PetscFree(rowlengths);
890b6490206SBarry Smith   PetscFree(browlengths);
891b6490206SBarry Smith   PetscFree(aa);
892b6490206SBarry Smith   PetscFree(jj);
893b6490206SBarry Smith   PetscFree(mask);
894b6490206SBarry Smith 
895b6490206SBarry Smith   B->assembled = PETSC_TRUE;
896de6a44a3SBarry Smith 
897de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr);
898de6a44a3SBarry Smith   if (flg) {
89919bcc07fSBarry Smith     Viewer tviewer;
90019bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
90190ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO,0);CHKERRQ(ierr);
90219bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
90319bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
904de6a44a3SBarry Smith   }
905de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr);
906de6a44a3SBarry Smith   if (flg) {
90719bcc07fSBarry Smith     Viewer tviewer;
90819bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
90990ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr);
91019bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
91119bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
912de6a44a3SBarry Smith   }
913de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg); CHKERRQ(ierr);
914de6a44a3SBarry Smith   if (flg) {
91519bcc07fSBarry Smith     Viewer tviewer;
91619bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
91719bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
91819bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
919de6a44a3SBarry Smith   }
920de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg); CHKERRQ(ierr);
921de6a44a3SBarry Smith   if (flg) {
92219bcc07fSBarry Smith     Viewer tviewer;
92319bcc07fSBarry Smith     ierr = ViewerFileOpenASCII(B->comm,"stdout",&tviewer);CHKERRQ(ierr);
92490ace30eSBarry Smith     ierr = ViewerSetFormat(tviewer,ASCII_FORMAT_MATLAB,"M");CHKERRQ(ierr);
92519bcc07fSBarry Smith     ierr = MatView(B,tviewer); CHKERRQ(ierr);
92619bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
927de6a44a3SBarry Smith   }
928de6a44a3SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr);
929de6a44a3SBarry Smith   if (flg) {
93019bcc07fSBarry Smith     Viewer tviewer;
93119bcc07fSBarry Smith     ierr = ViewerDrawOpenX(B->comm,0,0,0,0,300,300,&tviewer); CHKERRQ(ierr);
93219bcc07fSBarry Smith     ierr = MatView(B,(Viewer)tviewer); CHKERRQ(ierr);
93319bcc07fSBarry Smith     ierr = ViewerFlush(tviewer); CHKERRQ(ierr);
93419bcc07fSBarry Smith     ierr = ViewerDestroy(tviewer); CHKERRQ(ierr);
935de6a44a3SBarry Smith   }
9362593348eSBarry Smith   return 0;
9372593348eSBarry Smith }
9382593348eSBarry Smith 
9392593348eSBarry Smith 
9402593348eSBarry Smith 
941