xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 6550863b4fa1924fb7aed57c971a3f7e3e1b3e91)
1bba1ac68SSatish Balay /*$Id: baij.c,v 1.245 2001/08/07 03:02:55 balay Exp bsmith $*/
22593348eSBarry Smith 
32593348eSBarry Smith /*
4b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
52593348eSBarry Smith   matrix storage format.
62593348eSBarry Smith */
770f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
81a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
91a6a6d98SBarry Smith #include "src/inline/spops.h"
10325e03aeSBarry Smith #include "petscsys.h"                     /*I "petscmat.h" I*/
113b547af2SSatish Balay 
12af674e45SBarry Smith /*
13af674e45SBarry Smith     Special version for Fun3d sequential benchmark
14af674e45SBarry Smith */
15af674e45SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
16af674e45SBarry Smith #define matsetvaluesblocked4_ MATSETVALUESBLOCKED4
17af674e45SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
18af674e45SBarry Smith #define matsetvaluesblocked4_ matsetvaluesblocked4
19af674e45SBarry Smith #endif
20af674e45SBarry Smith 
219c8c1041SBarry Smith EXTERN_C_BEGIN
22af674e45SBarry Smith #undef __FUNCT__
23af674e45SBarry Smith #define __FUNCT__ "matsetvaluesblocked4_"
244bb09213Spetsc void matsetvaluesblocked4_(Mat *AA,int *mm,int *im,int *nn,int *in,PetscScalar *v)
25af674e45SBarry Smith {
26af674e45SBarry Smith   Mat         A = *AA;
27af674e45SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
28a037b02bSBarry Smith   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,N,m = *mm,n = *nn;
2938fde467SHong Zhang   int         *ai=a->i,*ailen=a->ilen;
30a037b02bSBarry Smith   int         *aj=a->j,stepval;
314bb09213Spetsc   PetscScalar *value = v;
324bb09213Spetsc   MatScalar   *ap,*aa = a->a,*bap;
33af674e45SBarry Smith 
34af674e45SBarry Smith   PetscFunctionBegin;
35af674e45SBarry Smith   stepval = (n-1)*4;
36af674e45SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
37af674e45SBarry Smith     row  = im[k];
38af674e45SBarry Smith     rp   = aj + ai[row];
39af674e45SBarry Smith     ap   = aa + 16*ai[row];
40af674e45SBarry Smith     nrow = ailen[row];
41af674e45SBarry Smith     low  = 0;
42af674e45SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
43af674e45SBarry Smith       col = in[l];
44af674e45SBarry Smith       value = v + k*(stepval+4)*4 + l*4;
45af674e45SBarry Smith       low = 0; high = nrow;
46af674e45SBarry Smith       while (high-low > 7) {
47af674e45SBarry Smith         t = (low+high)/2;
48af674e45SBarry Smith         if (rp[t] > col) high = t;
49af674e45SBarry Smith         else             low  = t;
50af674e45SBarry Smith       }
51af674e45SBarry Smith       for (i=low; i<high; i++) {
52af674e45SBarry Smith         if (rp[i] > col) break;
53af674e45SBarry Smith         if (rp[i] == col) {
54af674e45SBarry Smith           bap  = ap +  16*i;
55af674e45SBarry Smith           for (ii=0; ii<4; ii++,value+=stepval) {
56af674e45SBarry Smith             for (jj=ii; jj<16; jj+=4) {
57af674e45SBarry Smith               bap[jj] += *value++;
58af674e45SBarry Smith             }
59af674e45SBarry Smith           }
60af674e45SBarry Smith           goto noinsert2;
61af674e45SBarry Smith         }
62af674e45SBarry Smith       }
63af674e45SBarry Smith       N = nrow++ - 1;
64af674e45SBarry Smith       /* shift up all the later entries in this row */
65af674e45SBarry Smith       for (ii=N; ii>=i; ii--) {
66af674e45SBarry Smith         rp[ii+1] = rp[ii];
67a037b02bSBarry Smith         PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar));
68af674e45SBarry Smith       }
69af674e45SBarry Smith       if (N >= i) {
70a037b02bSBarry Smith         PetscMemzero(ap+16*i,16*sizeof(MatScalar));
71af674e45SBarry Smith       }
72af674e45SBarry Smith       rp[i] = col;
73af674e45SBarry Smith       bap   = ap +  16*i;
74af674e45SBarry Smith       for (ii=0; ii<4; ii++,value+=stepval) {
75af674e45SBarry Smith         for (jj=ii; jj<16; jj+=4) {
76af674e45SBarry Smith           bap[jj] = *value++;
77af674e45SBarry Smith         }
78af674e45SBarry Smith       }
79af674e45SBarry Smith       noinsert2:;
80af674e45SBarry Smith       low = i;
81af674e45SBarry Smith     }
82af674e45SBarry Smith     ailen[row] = nrow;
83af674e45SBarry Smith   }
84af674e45SBarry Smith }
859c8c1041SBarry Smith EXTERN_C_END
86af674e45SBarry Smith 
87af674e45SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
88af674e45SBarry Smith #define matsetvalues4_ MATSETVALUES4
89af674e45SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
90af674e45SBarry Smith #define matsetvalues4_ matsetvalues4
91af674e45SBarry Smith #endif
92af674e45SBarry Smith 
939c8c1041SBarry Smith EXTERN_C_BEGIN
94af674e45SBarry Smith #undef __FUNCT__
95af674e45SBarry Smith #define __FUNCT__ "MatSetValues4_"
964bb09213Spetsc void matsetvalues4_(Mat *AA,int *mm,int *im,int *nn,int *in,PetscScalar *v)
97af674e45SBarry Smith {
98af674e45SBarry Smith   Mat         A = *AA;
99af674e45SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
100a037b02bSBarry Smith   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,N,n = *nn,m = *mm;
10138fde467SHong Zhang   int         *ai=a->i,*ailen=a->ilen;
102af674e45SBarry Smith   int         *aj=a->j,brow,bcol;
10338fde467SHong Zhang   int         ridx,cidx;
104af674e45SBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
105af674e45SBarry Smith 
106af674e45SBarry Smith   PetscFunctionBegin;
107af674e45SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
108af674e45SBarry Smith     row  = im[k]; brow = row/4;
109af674e45SBarry Smith     rp   = aj + ai[brow];
110af674e45SBarry Smith     ap   = aa + 16*ai[brow];
111af674e45SBarry Smith     nrow = ailen[brow];
112af674e45SBarry Smith     low  = 0;
113af674e45SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
114af674e45SBarry Smith       col = in[l]; bcol = col/4;
115af674e45SBarry Smith       ridx = row % 4; cidx = col % 4;
116af674e45SBarry Smith       value = v[l + k*n];
117af674e45SBarry Smith       low = 0; high = nrow;
118af674e45SBarry Smith       while (high-low > 7) {
119af674e45SBarry Smith         t = (low+high)/2;
120af674e45SBarry Smith         if (rp[t] > bcol) high = t;
121af674e45SBarry Smith         else              low  = t;
122af674e45SBarry Smith       }
123af674e45SBarry Smith       for (i=low; i<high; i++) {
124af674e45SBarry Smith         if (rp[i] > bcol) break;
125af674e45SBarry Smith         if (rp[i] == bcol) {
126af674e45SBarry Smith           bap  = ap +  16*i + 4*cidx + ridx;
127af674e45SBarry Smith           *bap += value;
128af674e45SBarry Smith           goto noinsert1;
129af674e45SBarry Smith         }
130af674e45SBarry Smith       }
131af674e45SBarry Smith       N = nrow++ - 1;
132af674e45SBarry Smith       /* shift up all the later entries in this row */
133af674e45SBarry Smith       for (ii=N; ii>=i; ii--) {
134af674e45SBarry Smith         rp[ii+1] = rp[ii];
135a037b02bSBarry Smith         PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar));
136af674e45SBarry Smith       }
137af674e45SBarry Smith       if (N>=i) {
138a037b02bSBarry Smith         PetscMemzero(ap+16*i,16*sizeof(MatScalar));
139af674e45SBarry Smith       }
140af674e45SBarry Smith       rp[i]                    = bcol;
141af674e45SBarry Smith       ap[16*i + 4*cidx + ridx] = value;
142af674e45SBarry Smith       noinsert1:;
143af674e45SBarry Smith       low = i;
144af674e45SBarry Smith     }
145af674e45SBarry Smith     ailen[brow] = nrow;
146af674e45SBarry Smith   }
147af674e45SBarry Smith }
1489c8c1041SBarry Smith EXTERN_C_END
149af674e45SBarry Smith 
15095d5f7c2SBarry Smith /*  UGLY, ugly, ugly
15187828ca2SBarry Smith    When MatScalar == PetscScalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does
1523477af42SKris Buschelman    not exist. Otherwise ..._MatScalar() takes matrix dlements in single precision and
15395d5f7c2SBarry Smith    inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ()
15495d5f7c2SBarry Smith    converts the entries into single precision and then calls ..._MatScalar() to put them
15595d5f7c2SBarry Smith    into the single precision data structures.
15695d5f7c2SBarry Smith */
15795d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
158ca44d042SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode);
15995d5f7c2SBarry Smith #else
16095d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ
16195d5f7c2SBarry Smith #endif
16204929863SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
16304929863SHong Zhang EXTERN int MatUseDSCPACK_MPIBAIJ(Mat);
16404929863SHong Zhang #endif
16595d5f7c2SBarry Smith 
1662d61bbb3SSatish Balay #define CHUNKSIZE  10
167de6a44a3SBarry Smith 
168be5855fcSBarry Smith /*
169be5855fcSBarry Smith      Checks for missing diagonals
170be5855fcSBarry Smith */
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ"
173c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A)
174be5855fcSBarry Smith {
175be5855fcSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
176883fce79SBarry Smith   int         *diag,*jj = a->j,i,ierr;
177be5855fcSBarry Smith 
178be5855fcSBarry Smith   PetscFunctionBegin;
179c4992f7dSBarry Smith   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
180883fce79SBarry Smith   diag = a->diag;
1810e8e8aceSBarry Smith   for (i=0; i<a->mbs; i++) {
182be5855fcSBarry Smith     if (jj[diag[i]] != i) {
18329bbc08cSBarry Smith       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
184be5855fcSBarry Smith     }
185be5855fcSBarry Smith   }
186be5855fcSBarry Smith   PetscFunctionReturn(0);
187be5855fcSBarry Smith }
188be5855fcSBarry Smith 
1894a2ae208SSatish Balay #undef __FUNCT__
1904a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ"
191c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A)
192de6a44a3SBarry Smith {
193de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
19482502324SSatish Balay   int         i,j,*diag,m = a->mbs,ierr;
195de6a44a3SBarry Smith 
1963a40ed3dSBarry Smith   PetscFunctionBegin;
197883fce79SBarry Smith   if (a->diag) PetscFunctionReturn(0);
198883fce79SBarry Smith 
199b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
200b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
2017fc0212eSBarry Smith   for (i=0; i<m; i++) {
202ecc615b2SBarry Smith     diag[i] = a->i[i+1];
203de6a44a3SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
204de6a44a3SBarry Smith       if (a->j[j] == i) {
205de6a44a3SBarry Smith         diag[i] = j;
206de6a44a3SBarry Smith         break;
207de6a44a3SBarry Smith       }
208de6a44a3SBarry Smith     }
209de6a44a3SBarry Smith   }
210de6a44a3SBarry Smith   a->diag = diag;
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
212de6a44a3SBarry Smith }
2132593348eSBarry Smith 
2142593348eSBarry Smith 
215ca44d042SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
2163b2fbd54SBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqBAIJ"
2190e6d2581SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
2203b2fbd54SBarry Smith {
2213b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
2223b2fbd54SBarry Smith   int         ierr,n = a->mbs,i;
2233b2fbd54SBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
2253b2fbd54SBarry Smith   *nn = n;
2263a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
2273b2fbd54SBarry Smith   if (symmetric) {
2283b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
2293b2fbd54SBarry Smith   } else if (oshift == 1) {
2303b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
231435da068SBarry Smith     int nz = a->i[n];
2323b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
2333b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]++;
2343b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
2353b2fbd54SBarry Smith   } else {
2363b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
2373b2fbd54SBarry Smith   }
2383b2fbd54SBarry Smith 
2393a40ed3dSBarry Smith   PetscFunctionReturn(0);
2403b2fbd54SBarry Smith }
2413b2fbd54SBarry Smith 
2424a2ae208SSatish Balay #undef __FUNCT__
2434a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ"
244435da068SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
2453b2fbd54SBarry Smith {
2463b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
247606d414cSSatish Balay   int         i,n = a->mbs,ierr;
2483b2fbd54SBarry Smith 
2493a40ed3dSBarry Smith   PetscFunctionBegin;
2503a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
2513b2fbd54SBarry Smith   if (symmetric) {
252606d414cSSatish Balay     ierr = PetscFree(*ia);CHKERRQ(ierr);
253606d414cSSatish Balay     ierr = PetscFree(*ja);CHKERRQ(ierr);
254af5da2bfSBarry Smith   } else if (oshift == 1) {
255435da068SBarry Smith     int nz = a->i[n]-1;
2563b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
2573b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]--;
2583b2fbd54SBarry Smith   }
2593a40ed3dSBarry Smith   PetscFunctionReturn(0);
2603b2fbd54SBarry Smith }
2613b2fbd54SBarry Smith 
2624a2ae208SSatish Balay #undef __FUNCT__
2634a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqBAIJ"
2642d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs)
2652d61bbb3SSatish Balay {
2662d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data;
2672d61bbb3SSatish Balay 
2682d61bbb3SSatish Balay   PetscFunctionBegin;
2692d61bbb3SSatish Balay   *bs = baij->bs;
2702d61bbb3SSatish Balay   PetscFunctionReturn(0);
2712d61bbb3SSatish Balay }
2722d61bbb3SSatish Balay 
2732d61bbb3SSatish Balay 
2744a2ae208SSatish Balay #undef __FUNCT__
2754a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqBAIJ"
276e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
2772d61bbb3SSatish Balay {
2782d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
279e51c0b9cSSatish Balay   int         ierr;
2802d61bbb3SSatish Balay 
281433994e6SBarry Smith   PetscFunctionBegin;
282aa482453SBarry Smith #if defined(PETSC_USE_LOG)
283b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
2842d61bbb3SSatish Balay #endif
285606d414cSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
286606d414cSSatish Balay   if (!a->singlemalloc) {
287606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
288606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
289606d414cSSatish Balay   }
290c38d4ed2SBarry Smith   if (a->row) {
291c38d4ed2SBarry Smith     ierr = ISDestroy(a->row);CHKERRQ(ierr);
292c38d4ed2SBarry Smith   }
293c38d4ed2SBarry Smith   if (a->col) {
294c38d4ed2SBarry Smith     ierr = ISDestroy(a->col);CHKERRQ(ierr);
295c38d4ed2SBarry Smith   }
296606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
297606d414cSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
298606d414cSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
299606d414cSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
300606d414cSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
301e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
302606d414cSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
303563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
304563b5814SBarry Smith   if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);}
305563b5814SBarry Smith #endif
306c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
307c4319e64SHong Zhang 
308606d414cSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
3092d61bbb3SSatish Balay   PetscFunctionReturn(0);
3102d61bbb3SSatish Balay }
3112d61bbb3SSatish Balay 
3124a2ae208SSatish Balay #undef __FUNCT__
3134a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqBAIJ"
3142d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
3152d61bbb3SSatish Balay {
3162d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
3172d61bbb3SSatish Balay 
3182d61bbb3SSatish Balay   PetscFunctionBegin;
319aa275fccSKris Buschelman   switch (op) {
320aa275fccSKris Buschelman   case MAT_ROW_ORIENTED:
321aa275fccSKris Buschelman     a->roworiented    = PETSC_TRUE;
322aa275fccSKris Buschelman     break;
323aa275fccSKris Buschelman   case MAT_COLUMN_ORIENTED:
324aa275fccSKris Buschelman     a->roworiented    = PETSC_FALSE;
325aa275fccSKris Buschelman     break;
326aa275fccSKris Buschelman   case MAT_COLUMNS_SORTED:
327aa275fccSKris Buschelman     a->sorted         = PETSC_TRUE;
328aa275fccSKris Buschelman     break;
329aa275fccSKris Buschelman   case MAT_COLUMNS_UNSORTED:
330aa275fccSKris Buschelman     a->sorted         = PETSC_FALSE;
331aa275fccSKris Buschelman     break;
332aa275fccSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
333aa275fccSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
334aa275fccSKris Buschelman     break;
335aa275fccSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
336aa275fccSKris Buschelman     a->nonew          = 1;
337aa275fccSKris Buschelman     break;
338aa275fccSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
339aa275fccSKris Buschelman     a->nonew          = -1;
340aa275fccSKris Buschelman     break;
341aa275fccSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
342aa275fccSKris Buschelman     a->nonew          = -2;
343aa275fccSKris Buschelman     break;
344aa275fccSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
345aa275fccSKris Buschelman     a->nonew          = 0;
346aa275fccSKris Buschelman     break;
347aa275fccSKris Buschelman   case MAT_ROWS_SORTED:
348aa275fccSKris Buschelman   case MAT_ROWS_UNSORTED:
349aa275fccSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
350aa275fccSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
351aa275fccSKris Buschelman   case MAT_USE_HASH_TABLE:
352b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
353aa275fccSKris Buschelman     break;
354aa275fccSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
35529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
356aa275fccSKris Buschelman   default:
35729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
3582d61bbb3SSatish Balay   }
3592d61bbb3SSatish Balay   PetscFunctionReturn(0);
3602d61bbb3SSatish Balay }
3612d61bbb3SSatish Balay 
3624a2ae208SSatish Balay #undef __FUNCT__
3634a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ"
36487828ca2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
3652d61bbb3SSatish Balay {
3662d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
36782502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr;
3683f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
36987828ca2SBarry Smith   PetscScalar  *v_i;
3702d61bbb3SSatish Balay 
3712d61bbb3SSatish Balay   PetscFunctionBegin;
3722d61bbb3SSatish Balay   bs  = a->bs;
3732d61bbb3SSatish Balay   ai  = a->i;
3742d61bbb3SSatish Balay   aj  = a->j;
3752d61bbb3SSatish Balay   aa  = a->a;
3762d61bbb3SSatish Balay   bs2 = a->bs2;
3772d61bbb3SSatish Balay 
378273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
3792d61bbb3SSatish Balay 
3802d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
3812d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
3822d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
3832d61bbb3SSatish Balay   *nz = bs*M;
3842d61bbb3SSatish Balay 
3852d61bbb3SSatish Balay   if (v) {
3862d61bbb3SSatish Balay     *v = 0;
3872d61bbb3SSatish Balay     if (*nz) {
38887828ca2SBarry Smith       ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr);
3892d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
3902d61bbb3SSatish Balay         v_i  = *v + i*bs;
3912d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
3922d61bbb3SSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
3932d61bbb3SSatish Balay       }
3942d61bbb3SSatish Balay     }
3952d61bbb3SSatish Balay   }
3962d61bbb3SSatish Balay 
3972d61bbb3SSatish Balay   if (idx) {
3982d61bbb3SSatish Balay     *idx = 0;
3992d61bbb3SSatish Balay     if (*nz) {
400b0a32e0cSBarry Smith       ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr);
4012d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
4022d61bbb3SSatish Balay         idx_i = *idx + i*bs;
4032d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
4042d61bbb3SSatish Balay         for (j=0; j<bs; j++) {idx_i[j] = itmp++;}
4052d61bbb3SSatish Balay       }
4062d61bbb3SSatish Balay     }
4072d61bbb3SSatish Balay   }
4082d61bbb3SSatish Balay   PetscFunctionReturn(0);
4092d61bbb3SSatish Balay }
4102d61bbb3SSatish Balay 
4114a2ae208SSatish Balay #undef __FUNCT__
4124a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ"
41387828ca2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
4142d61bbb3SSatish Balay {
415606d414cSSatish Balay   int ierr;
416606d414cSSatish Balay 
4172d61bbb3SSatish Balay   PetscFunctionBegin;
418606d414cSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
419606d414cSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
4202d61bbb3SSatish Balay   PetscFunctionReturn(0);
4212d61bbb3SSatish Balay }
4222d61bbb3SSatish Balay 
4234a2ae208SSatish Balay #undef __FUNCT__
4244a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ"
4252d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
4262d61bbb3SSatish Balay {
4272d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
4282d61bbb3SSatish Balay   Mat         C;
4292d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
430273d9f13SBarry Smith   int         *rows,*cols,bs2=a->bs2;
43187828ca2SBarry Smith   PetscScalar *array;
4322d61bbb3SSatish Balay 
4332d61bbb3SSatish Balay   PetscFunctionBegin;
43429bbc08cSBarry Smith   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place");
435b0a32e0cSBarry Smith   ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr);
436549d3d68SSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
4372d61bbb3SSatish Balay 
438375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
43987828ca2SBarry Smith   ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr);
44087828ca2SBarry Smith   for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i];
441375fe846SBarry Smith #else
4423eda8832SBarry Smith   array = a->a;
443375fe846SBarry Smith #endif
444375fe846SBarry Smith 
4452d61bbb3SSatish Balay   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
446273d9f13SBarry Smith   ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
447606d414cSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
448b0a32e0cSBarry Smith   ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr);
4492d61bbb3SSatish Balay   cols = rows + bs;
4502d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4512d61bbb3SSatish Balay     cols[0] = i*bs;
4522d61bbb3SSatish Balay     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
4532d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
4542d61bbb3SSatish Balay     for (j=0; j<len; j++) {
4552d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
4562d61bbb3SSatish Balay       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
4572d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
4582d61bbb3SSatish Balay       array += bs2;
4592d61bbb3SSatish Balay     }
4602d61bbb3SSatish Balay   }
461606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
462375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
463375fe846SBarry Smith   ierr = PetscFree(array);
464375fe846SBarry Smith #endif
4652d61bbb3SSatish Balay 
4662d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4672d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4682d61bbb3SSatish Balay 
469c4992f7dSBarry Smith   if (B) {
4702d61bbb3SSatish Balay     *B = C;
4712d61bbb3SSatish Balay   } else {
472273d9f13SBarry Smith     ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
4732d61bbb3SSatish Balay   }
4742d61bbb3SSatish Balay   PetscFunctionReturn(0);
4752d61bbb3SSatish Balay }
4762d61bbb3SSatish Balay 
4774a2ae208SSatish Balay #undef __FUNCT__
4784a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary"
479b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer)
4802593348eSBarry Smith {
481b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
4829df24120SSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
48387828ca2SBarry Smith   PetscScalar *aa;
484ce6f0cecSBarry Smith   FILE        *file;
4852593348eSBarry Smith 
4863a40ed3dSBarry Smith   PetscFunctionBegin;
487b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
488b0a32e0cSBarry Smith   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
489552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
4903b2fbd54SBarry Smith 
491273d9f13SBarry Smith   col_lens[1] = A->m;
492273d9f13SBarry Smith   col_lens[2] = A->n;
4937e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
4942593348eSBarry Smith 
4952593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
496b6490206SBarry Smith   count = 0;
497b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
498b6490206SBarry Smith     for (j=0; j<bs; j++) {
499b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
500b6490206SBarry Smith     }
5012593348eSBarry Smith   }
502273d9f13SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
503606d414cSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
5042593348eSBarry Smith 
5052593348eSBarry Smith   /* store column indices (zero start index) */
506b0a32e0cSBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
507b6490206SBarry Smith   count = 0;
508b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
509b6490206SBarry Smith     for (j=0; j<bs; j++) {
510b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
511b6490206SBarry Smith         for (l=0; l<bs; l++) {
512b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
5132593348eSBarry Smith         }
5142593348eSBarry Smith       }
515b6490206SBarry Smith     }
516b6490206SBarry Smith   }
5170752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
518606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
5192593348eSBarry Smith 
5202593348eSBarry Smith   /* store nonzero values */
52187828ca2SBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
522b6490206SBarry Smith   count = 0;
523b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
524b6490206SBarry Smith     for (j=0; j<bs; j++) {
525b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
526b6490206SBarry Smith         for (l=0; l<bs; l++) {
5277e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
528b6490206SBarry Smith         }
529b6490206SBarry Smith       }
530b6490206SBarry Smith     }
531b6490206SBarry Smith   }
5320752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
533606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
534ce6f0cecSBarry Smith 
535b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
536ce6f0cecSBarry Smith   if (file) {
537ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
538ce6f0cecSBarry Smith   }
5393a40ed3dSBarry Smith   PetscFunctionReturn(0);
5402593348eSBarry Smith }
5412593348eSBarry Smith 
54204929863SHong Zhang extern int MatMPIBAIJFactorInfo_DSCPACK(Mat,PetscViewer);
54304929863SHong Zhang 
5444a2ae208SSatish Balay #undef __FUNCT__
5454a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII"
546b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer)
5472593348eSBarry Smith {
548b6490206SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
549fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
550f3ef73ceSBarry Smith   PetscViewerFormat format;
5512593348eSBarry Smith 
5523a40ed3dSBarry Smith   PetscFunctionBegin;
553b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
554456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
555b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
556fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
557bcd9e38bSBarry Smith     Mat aij;
558bcd9e38bSBarry Smith     ierr = MatConvert(A,MATSEQAIJ,&aij);CHKERRQ(ierr);
559bcd9e38bSBarry Smith     ierr = MatView(aij,viewer);CHKERRQ(ierr);
560bcd9e38bSBarry Smith     ierr = MatDestroy(aij);CHKERRQ(ierr);
56104929863SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
56204929863SHong Zhang #if defined(PETSC_HAVE_DSCPACK) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
56304929863SHong Zhang      ierr = MatMPIBAIJFactorInfo_DSCPACK(A,viewer);CHKERRQ(ierr);
56404929863SHong Zhang #endif
56504929863SHong Zhang      PetscFunctionReturn(0);
566fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
567b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
56844cd7ae7SLois Curfman McInnes     for (i=0; i<a->mbs; i++) {
56944cd7ae7SLois Curfman McInnes       for (j=0; j<bs; j++) {
570b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
57144cd7ae7SLois Curfman McInnes         for (k=a->i[i]; k<a->i[i+1]; k++) {
57244cd7ae7SLois Curfman McInnes           for (l=0; l<bs; l++) {
573aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5740e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
57562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %gi) ",bs*a->j[k]+l,
5760e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
5770e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
57862b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %gi) ",bs*a->j[k]+l,
5790e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
5800e6d2581SBarry Smith             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
58162b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
5820ef38995SBarry Smith             }
58344cd7ae7SLois Curfman McInnes #else
5840ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
58562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
5860ef38995SBarry Smith             }
58744cd7ae7SLois Curfman McInnes #endif
58844cd7ae7SLois Curfman McInnes           }
58944cd7ae7SLois Curfman McInnes         }
590b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
59144cd7ae7SLois Curfman McInnes       }
59244cd7ae7SLois Curfman McInnes     }
593b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
5940ef38995SBarry Smith   } else {
595b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
596b6490206SBarry Smith     for (i=0; i<a->mbs; i++) {
597b6490206SBarry Smith       for (j=0; j<bs; j++) {
598b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
599b6490206SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
600b6490206SBarry Smith           for (l=0; l<bs; l++) {
601aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6020e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
60362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
6040e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
6050e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
60662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
6070e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
6080ef38995SBarry Smith             } else {
60962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
61088685aaeSLois Curfman McInnes             }
61188685aaeSLois Curfman McInnes #else
61262b941d6SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
61388685aaeSLois Curfman McInnes #endif
6142593348eSBarry Smith           }
6152593348eSBarry Smith         }
616b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6172593348eSBarry Smith       }
6182593348eSBarry Smith     }
619b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
620b6490206SBarry Smith   }
621b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6223a40ed3dSBarry Smith   PetscFunctionReturn(0);
6232593348eSBarry Smith }
6242593348eSBarry Smith 
6254a2ae208SSatish Balay #undef __FUNCT__
6264a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom"
627b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
6283270192aSSatish Balay {
62977ed5343SBarry Smith   Mat          A = (Mat) Aa;
6303270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
631b65db4caSBarry Smith   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2;
6320e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
6333f1db9ecSBarry Smith   MatScalar    *aa;
634b0a32e0cSBarry Smith   PetscViewer  viewer;
6353270192aSSatish Balay 
6363a40ed3dSBarry Smith   PetscFunctionBegin;
6373270192aSSatish Balay 
638b65db4caSBarry Smith   /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/
63977ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
64077ed5343SBarry Smith 
641b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
64277ed5343SBarry Smith 
6433270192aSSatish Balay   /* loop over matrix elements drawing boxes */
644b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
6453270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
6463270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
647273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
6483270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
6493270192aSSatish Balay       aa = a->a + j*bs2;
6503270192aSSatish Balay       for (k=0; k<bs; k++) {
6513270192aSSatish Balay         for (l=0; l<bs; l++) {
6520e6d2581SBarry Smith           if (PetscRealPart(*aa++) >=  0.) continue;
653b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
6543270192aSSatish Balay         }
6553270192aSSatish Balay       }
6563270192aSSatish Balay     }
6573270192aSSatish Balay   }
658b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
6593270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
6603270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
661273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
6623270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
6633270192aSSatish Balay       aa = a->a + j*bs2;
6643270192aSSatish Balay       for (k=0; k<bs; k++) {
6653270192aSSatish Balay         for (l=0; l<bs; l++) {
6660e6d2581SBarry Smith           if (PetscRealPart(*aa++) != 0.) continue;
667b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
6683270192aSSatish Balay         }
6693270192aSSatish Balay       }
6703270192aSSatish Balay     }
6713270192aSSatish Balay   }
6723270192aSSatish Balay 
673b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
6743270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
6753270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
676273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
6773270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
6783270192aSSatish Balay       aa = a->a + j*bs2;
6793270192aSSatish Balay       for (k=0; k<bs; k++) {
6803270192aSSatish Balay         for (l=0; l<bs; l++) {
6810e6d2581SBarry Smith           if (PetscRealPart(*aa++) <= 0.) continue;
682b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
6833270192aSSatish Balay         }
6843270192aSSatish Balay       }
6853270192aSSatish Balay     }
6863270192aSSatish Balay   }
68777ed5343SBarry Smith   PetscFunctionReturn(0);
68877ed5343SBarry Smith }
6893270192aSSatish Balay 
6904a2ae208SSatish Balay #undef __FUNCT__
6914a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw"
692b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer)
69377ed5343SBarry Smith {
6947dce120fSSatish Balay   int          ierr;
6950e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,w,h;
696b0a32e0cSBarry Smith   PetscDraw    draw;
69777ed5343SBarry Smith   PetscTruth   isnull;
6983270192aSSatish Balay 
69977ed5343SBarry Smith   PetscFunctionBegin;
70077ed5343SBarry Smith 
701b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
702b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
70377ed5343SBarry Smith 
70477ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
705273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
70677ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
707b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
708b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
70977ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
7103a40ed3dSBarry Smith   PetscFunctionReturn(0);
7113270192aSSatish Balay }
7123270192aSSatish Balay 
7134a2ae208SSatish Balay #undef __FUNCT__
7144a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ"
715b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer)
7162593348eSBarry Smith {
71719bcc07fSBarry Smith   int        ierr;
7186831982aSBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
7192593348eSBarry Smith 
7203a40ed3dSBarry Smith   PetscFunctionBegin;
721b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
722b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
723fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
724fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
7250f5bd95cSBarry Smith   if (issocket) {
72629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
7270f5bd95cSBarry Smith   } else if (isascii){
7283a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
7290f5bd95cSBarry Smith   } else if (isbinary) {
7303a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
7310f5bd95cSBarry Smith   } else if (isdraw) {
7323a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
7335cd90555SBarry Smith   } else {
73429bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
7352593348eSBarry Smith   }
7363a40ed3dSBarry Smith   PetscFunctionReturn(0);
7372593348eSBarry Smith }
738b6490206SBarry Smith 
739cd0e1443SSatish Balay 
7404a2ae208SSatish Balay #undef __FUNCT__
7414a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ"
74287828ca2SBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
743cd0e1443SSatish Balay {
744cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
7452d61bbb3SSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
7462d61bbb3SSatish Balay   int        *ai = a->i,*ailen = a->ilen;
7472d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
7483f1db9ecSBarry Smith   MatScalar  *ap,*aa = a->a,zero = 0.0;
749cd0e1443SSatish Balay 
7503a40ed3dSBarry Smith   PetscFunctionBegin;
7512d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
752cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
75329bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
754273d9f13SBarry Smith     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
7552d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
7562c3acbe9SBarry Smith     nrow = ailen[brow];
7572d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
75829bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
759273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
7602d61bbb3SSatish Balay       col  = in[l] ;
7612d61bbb3SSatish Balay       bcol = col/bs;
7622d61bbb3SSatish Balay       cidx = col%bs;
7632d61bbb3SSatish Balay       ridx = row%bs;
7642d61bbb3SSatish Balay       high = nrow;
7652d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
7662d61bbb3SSatish Balay       while (high-low > 5) {
767cd0e1443SSatish Balay         t = (low+high)/2;
768cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
769cd0e1443SSatish Balay         else             low  = t;
770cd0e1443SSatish Balay       }
771cd0e1443SSatish Balay       for (i=low; i<high; i++) {
772cd0e1443SSatish Balay         if (rp[i] > bcol) break;
773cd0e1443SSatish Balay         if (rp[i] == bcol) {
7742d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
7752d61bbb3SSatish Balay           goto finished;
776cd0e1443SSatish Balay         }
777cd0e1443SSatish Balay       }
7782d61bbb3SSatish Balay       *v++ = zero;
7792d61bbb3SSatish Balay       finished:;
780cd0e1443SSatish Balay     }
781cd0e1443SSatish Balay   }
7823a40ed3dSBarry Smith   PetscFunctionReturn(0);
783cd0e1443SSatish Balay }
784cd0e1443SSatish Balay 
78595d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
7864a2ae208SSatish Balay #undef __FUNCT__
7874a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
78887828ca2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,PetscScalar *v,InsertMode addv)
78995d5f7c2SBarry Smith {
790563b5814SBarry Smith   Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data;
791563b5814SBarry Smith   int         ierr,i,N = m*n*b->bs2;
792563b5814SBarry Smith   MatScalar   *vsingle;
79395d5f7c2SBarry Smith 
79495d5f7c2SBarry Smith   PetscFunctionBegin;
795563b5814SBarry Smith   if (N > b->setvalueslen) {
796563b5814SBarry Smith     if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);}
797b0a32e0cSBarry Smith     ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr);
798563b5814SBarry Smith     b->setvalueslen  = N;
799563b5814SBarry Smith   }
800563b5814SBarry Smith   vsingle = b->setvaluescopy;
80195d5f7c2SBarry Smith   for (i=0; i<N; i++) {
80295d5f7c2SBarry Smith     vsingle[i] = v[i];
80395d5f7c2SBarry Smith   }
80495d5f7c2SBarry Smith   ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr);
80595d5f7c2SBarry Smith   PetscFunctionReturn(0);
80695d5f7c2SBarry Smith }
80795d5f7c2SBarry Smith #endif
80895d5f7c2SBarry Smith 
8092d61bbb3SSatish Balay 
8104a2ae208SSatish Balay #undef __FUNCT__
8114a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
81295d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is)
81392c4ed94SBarry Smith {
81492c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
8158a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
816273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
817549d3d68SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
818273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
81995d5f7c2SBarry Smith   MatScalar   *value = v,*ap,*aa = a->a,*bap;
82092c4ed94SBarry Smith 
8213a40ed3dSBarry Smith   PetscFunctionBegin;
8220e324ae4SSatish Balay   if (roworiented) {
8230e324ae4SSatish Balay     stepval = (n-1)*bs;
8240e324ae4SSatish Balay   } else {
8250e324ae4SSatish Balay     stepval = (m-1)*bs;
8260e324ae4SSatish Balay   }
82792c4ed94SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
82892c4ed94SBarry Smith     row  = im[k];
8295ef9f2a5SBarry Smith     if (row < 0) continue;
830aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
83129bbc08cSBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
83292c4ed94SBarry Smith #endif
83392c4ed94SBarry Smith     rp   = aj + ai[row];
83492c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
83592c4ed94SBarry Smith     rmax = imax[row];
83692c4ed94SBarry Smith     nrow = ailen[row];
83792c4ed94SBarry Smith     low  = 0;
83892c4ed94SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
8395ef9f2a5SBarry Smith       if (in[l] < 0) continue;
840aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
84129bbc08cSBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
84292c4ed94SBarry Smith #endif
84392c4ed94SBarry Smith       col = in[l];
84492c4ed94SBarry Smith       if (roworiented) {
8450e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
8460e324ae4SSatish Balay       } else {
8470e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
84892c4ed94SBarry Smith       }
84992c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
85092c4ed94SBarry Smith       while (high-low > 7) {
85192c4ed94SBarry Smith         t = (low+high)/2;
85292c4ed94SBarry Smith         if (rp[t] > col) high = t;
85392c4ed94SBarry Smith         else             low  = t;
85492c4ed94SBarry Smith       }
85592c4ed94SBarry Smith       for (i=low; i<high; i++) {
85692c4ed94SBarry Smith         if (rp[i] > col) break;
85792c4ed94SBarry Smith         if (rp[i] == col) {
8588a84c255SSatish Balay           bap  = ap +  bs2*i;
8590e324ae4SSatish Balay           if (roworiented) {
8608a84c255SSatish Balay             if (is == ADD_VALUES) {
861dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
862dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
8638a84c255SSatish Balay                   bap[jj] += *value++;
864dd9472c6SBarry Smith                 }
865dd9472c6SBarry Smith               }
8660e324ae4SSatish Balay             } else {
867dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
868dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
8690e324ae4SSatish Balay                   bap[jj] = *value++;
8708a84c255SSatish Balay                 }
871dd9472c6SBarry Smith               }
872dd9472c6SBarry Smith             }
8730e324ae4SSatish Balay           } else {
8740e324ae4SSatish Balay             if (is == ADD_VALUES) {
875dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
876dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
8770e324ae4SSatish Balay                   *bap++ += *value++;
878dd9472c6SBarry Smith                 }
879dd9472c6SBarry Smith               }
8800e324ae4SSatish Balay             } else {
881dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
882dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
8830e324ae4SSatish Balay                   *bap++  = *value++;
8840e324ae4SSatish Balay                 }
8858a84c255SSatish Balay               }
886dd9472c6SBarry Smith             }
887dd9472c6SBarry Smith           }
888f1241b54SBarry Smith           goto noinsert2;
88992c4ed94SBarry Smith         }
89092c4ed94SBarry Smith       }
89189280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
89229bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
89392c4ed94SBarry Smith       if (nrow >= rmax) {
89492c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
89592c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
8963f1db9ecSBarry Smith         MatScalar *new_a;
89792c4ed94SBarry Smith 
89829bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
89989280ab3SLois Curfman McInnes 
90092c4ed94SBarry Smith         /* malloc new storage space */
9013f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
902b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
90392c4ed94SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
90492c4ed94SBarry Smith         new_i   = new_j + new_nz;
90592c4ed94SBarry Smith 
90692c4ed94SBarry Smith         /* copy over old data into new slots */
90792c4ed94SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
90892c4ed94SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
909549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
91092c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
911549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
912549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
913549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
914549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
91592c4ed94SBarry Smith         /* free up old matrix storage */
916606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
917606d414cSSatish Balay         if (!a->singlemalloc) {
918606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
919606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
920606d414cSSatish Balay         }
92192c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
922c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
92392c4ed94SBarry Smith 
92492c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
92592c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
926b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
92792c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
92892c4ed94SBarry Smith         a->reallocs++;
92992c4ed94SBarry Smith         a->nz++;
93092c4ed94SBarry Smith       }
93192c4ed94SBarry Smith       N = nrow++ - 1;
93292c4ed94SBarry Smith       /* shift up all the later entries in this row */
93392c4ed94SBarry Smith       for (ii=N; ii>=i; ii--) {
93492c4ed94SBarry Smith         rp[ii+1] = rp[ii];
935549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
93692c4ed94SBarry Smith       }
937549d3d68SSatish Balay       if (N >= i) {
938549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
939549d3d68SSatish Balay       }
94092c4ed94SBarry Smith       rp[i] = col;
9418a84c255SSatish Balay       bap   = ap +  bs2*i;
9420e324ae4SSatish Balay       if (roworiented) {
943dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
944dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
9450e324ae4SSatish Balay             bap[jj] = *value++;
946dd9472c6SBarry Smith           }
947dd9472c6SBarry Smith         }
9480e324ae4SSatish Balay       } else {
949dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
950dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++) {
9510e324ae4SSatish Balay             *bap++  = *value++;
9520e324ae4SSatish Balay           }
953dd9472c6SBarry Smith         }
954dd9472c6SBarry Smith       }
955f1241b54SBarry Smith       noinsert2:;
95692c4ed94SBarry Smith       low = i;
95792c4ed94SBarry Smith     }
95892c4ed94SBarry Smith     ailen[row] = nrow;
95992c4ed94SBarry Smith   }
9603a40ed3dSBarry Smith   PetscFunctionReturn(0);
96192c4ed94SBarry Smith }
96292c4ed94SBarry Smith 
9634a2ae208SSatish Balay #undef __FUNCT__
9644a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ"
9658f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
966584200bdSSatish Balay {
967584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
968584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
969273d9f13SBarry Smith   int        m = A->m,*ip,N,*ailen = a->ilen;
970549d3d68SSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
9713f1db9ecSBarry Smith   MatScalar  *aa = a->a,*ap;
972a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
973a56a16c8SHong Zhang   PetscTruth   flag;
974a56a16c8SHong Zhang #endif
975584200bdSSatish Balay 
9763a40ed3dSBarry Smith   PetscFunctionBegin;
9773a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
978584200bdSSatish Balay 
97943ee02c3SBarry Smith   if (m) rmax = ailen[0];
980584200bdSSatish Balay   for (i=1; i<mbs; i++) {
981584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
982584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
983d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
984584200bdSSatish Balay     if (fshift) {
985a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
986584200bdSSatish Balay       N = ailen[i];
987584200bdSSatish Balay       for (j=0; j<N; j++) {
988584200bdSSatish Balay         ip[j-fshift] = ip[j];
989549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
990584200bdSSatish Balay       }
991584200bdSSatish Balay     }
992584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
993584200bdSSatish Balay   }
994584200bdSSatish Balay   if (mbs) {
995584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
996584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
997584200bdSSatish Balay   }
998584200bdSSatish Balay   /* reset ilen and imax for each row */
999584200bdSSatish Balay   for (i=0; i<mbs; i++) {
1000584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
1001584200bdSSatish Balay   }
1002a7c10996SSatish Balay   a->nz = ai[mbs];
1003584200bdSSatish Balay 
1004584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
1005584200bdSSatish Balay   if (fshift && a->diag) {
1006606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
1007b424e231SHong Zhang     PetscLogObjectMemory(A,-(mbs+1)*sizeof(int));
1008584200bdSSatish Balay     a->diag = 0;
1009584200bdSSatish Balay   }
1010bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",m,A->n,a->bs,fshift*bs2,a->nz*bs2);
1011bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs);
1012b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
1013e2f3b5e9SSatish Balay   a->reallocs          = 0;
10140e6d2581SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
1015cf4441caSHong Zhang 
1016a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
10172c535e4dSHong Zhang   ierr = PetscOptionsHasName(A->prefix,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
1018a56a16c8SHong Zhang   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(A);CHKERRQ(ierr); }
1019a56a16c8SHong Zhang #endif
1020a56a16c8SHong Zhang 
10213a40ed3dSBarry Smith   PetscFunctionReturn(0);
1022584200bdSSatish Balay }
1023584200bdSSatish Balay 
10245a838052SSatish Balay 
1025bea157c4SSatish Balay 
1026bea157c4SSatish Balay /*
1027bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
1028bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
1029bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
1030bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
1031bea157c4SSatish Balay */
10324a2ae208SSatish Balay #undef __FUNCT__
10334a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks"
1034bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
1035d9b7c43dSSatish Balay {
1036bea157c4SSatish Balay   int        i,j,k,row;
1037bea157c4SSatish Balay   PetscTruth flg;
10383a40ed3dSBarry Smith 
1039433994e6SBarry Smith   PetscFunctionBegin;
1040bea157c4SSatish Balay   for (i=0,j=0; i<n; j++) {
1041bea157c4SSatish Balay     row = idx[i];
1042bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
1043bea157c4SSatish Balay       sizes[j] = 1;
1044bea157c4SSatish Balay       i++;
1045e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
1046bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
1047bea157c4SSatish Balay       i++;
1048bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
1049bea157c4SSatish Balay       flg = PETSC_TRUE;
1050bea157c4SSatish Balay       for (k=1; k<bs; k++) {
1051bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
1052bea157c4SSatish Balay           flg = PETSC_FALSE;
1053bea157c4SSatish Balay           break;
1054d9b7c43dSSatish Balay         }
1055bea157c4SSatish Balay       }
1056bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
1057bea157c4SSatish Balay         sizes[j] = bs;
1058bea157c4SSatish Balay         i+= bs;
1059bea157c4SSatish Balay       } else {
1060bea157c4SSatish Balay         sizes[j] = 1;
1061bea157c4SSatish Balay         i++;
1062bea157c4SSatish Balay       }
1063bea157c4SSatish Balay     }
1064bea157c4SSatish Balay   }
1065bea157c4SSatish Balay   *bs_max = j;
10663a40ed3dSBarry Smith   PetscFunctionReturn(0);
1067d9b7c43dSSatish Balay }
1068d9b7c43dSSatish Balay 
10694a2ae208SSatish Balay #undef __FUNCT__
10704a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ"
107187828ca2SBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,PetscScalar *diag)
1072d9b7c43dSSatish Balay {
1073d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
1074b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
1075bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
107687828ca2SBarry Smith   PetscScalar zero = 0.0;
10773f1db9ecSBarry Smith   MatScalar   *aa;
1078d9b7c43dSSatish Balay 
10793a40ed3dSBarry Smith   PetscFunctionBegin;
1080d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
1081b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
1082d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
1083d9b7c43dSSatish Balay 
1084bea157c4SSatish Balay   /* allocate memory for rows,sizes */
1085b0a32e0cSBarry Smith   ierr  = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
1086bea157c4SSatish Balay   sizes = rows + is_n;
1087bea157c4SSatish Balay 
1088563b5814SBarry Smith   /* copy IS values to rows, and sort them */
1089bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
1090bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
1091dffd3267SBarry Smith   if (baij->keepzeroedrows) {
1092dffd3267SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
1093dffd3267SBarry Smith     bs_max = is_n;
1094dffd3267SBarry Smith   } else {
1095bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
1096dffd3267SBarry Smith   }
1097b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
1098bea157c4SSatish Balay 
1099bea157c4SSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
1100bea157c4SSatish Balay     row   = rows[j];
1101273d9f13SBarry Smith     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
1102bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1103bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
1104dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
1105bea157c4SSatish Balay       if (diag) {
1106bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
1107bea157c4SSatish Balay           baij->ilen[row/bs]       = 1;
1108bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
1109bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
1110a07cd24cSSatish Balay         }
1111563b5814SBarry Smith         /* Now insert all the diagonal values for this bs */
1112bea157c4SSatish Balay         for (k=0; k<bs; k++) {
1113bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
1114bea157c4SSatish Balay         }
1115bea157c4SSatish Balay       } else { /* (!diag) */
1116bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
1117bea157c4SSatish Balay       } /* end (!diag) */
1118bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
1119aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
112029bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
1121bea157c4SSatish Balay #endif
1122bea157c4SSatish Balay       for (k=0; k<count; k++) {
1123d9b7c43dSSatish Balay         aa[0] =  zero;
1124d9b7c43dSSatish Balay         aa    += bs;
1125d9b7c43dSSatish Balay       }
1126d9b7c43dSSatish Balay       if (diag) {
1127f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1128d9b7c43dSSatish Balay       }
1129d9b7c43dSSatish Balay     }
1130bea157c4SSatish Balay   }
1131bea157c4SSatish Balay 
1132606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
11339a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1135d9b7c43dSSatish Balay }
11361c351548SSatish Balay 
11374a2ae208SSatish Balay #undef __FUNCT__
11384a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ"
113987828ca2SBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
11402d61bbb3SSatish Balay {
11412d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
11422d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1143273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
11442d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1145549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
1146273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
11473f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
11482d61bbb3SSatish Balay 
11492d61bbb3SSatish Balay   PetscFunctionBegin;
11502d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
11512d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
11525ef9f2a5SBarry Smith     if (row < 0) continue;
1153aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1154273d9f13SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
11552d61bbb3SSatish Balay #endif
11562d61bbb3SSatish Balay     rp   = aj + ai[brow];
11572d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
11582d61bbb3SSatish Balay     rmax = imax[brow];
11592d61bbb3SSatish Balay     nrow = ailen[brow];
11602d61bbb3SSatish Balay     low  = 0;
11612d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
11625ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1163aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1164273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
11652d61bbb3SSatish Balay #endif
11662d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
11672d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
11682d61bbb3SSatish Balay       if (roworiented) {
11695ef9f2a5SBarry Smith         value = v[l + k*n];
11702d61bbb3SSatish Balay       } else {
11712d61bbb3SSatish Balay         value = v[k + l*m];
11722d61bbb3SSatish Balay       }
11732d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
11742d61bbb3SSatish Balay       while (high-low > 7) {
11752d61bbb3SSatish Balay         t = (low+high)/2;
11762d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
11772d61bbb3SSatish Balay         else              low  = t;
11782d61bbb3SSatish Balay       }
11792d61bbb3SSatish Balay       for (i=low; i<high; i++) {
11802d61bbb3SSatish Balay         if (rp[i] > bcol) break;
11812d61bbb3SSatish Balay         if (rp[i] == bcol) {
11822d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
11832d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
11842d61bbb3SSatish Balay           else                  *bap  = value;
11852d61bbb3SSatish Balay           goto noinsert1;
11862d61bbb3SSatish Balay         }
11872d61bbb3SSatish Balay       }
11882d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
118929bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
11902d61bbb3SSatish Balay       if (nrow >= rmax) {
11912d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
11922d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
11933f1db9ecSBarry Smith         MatScalar *new_a;
11942d61bbb3SSatish Balay 
119529bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
11962d61bbb3SSatish Balay 
11972d61bbb3SSatish Balay         /* Malloc new storage space */
11983f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1199b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
12002d61bbb3SSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
12012d61bbb3SSatish Balay         new_i   = new_j + new_nz;
12022d61bbb3SSatish Balay 
12032d61bbb3SSatish Balay         /* copy over old data into new slots */
12042d61bbb3SSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
12052d61bbb3SSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1206549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
12072d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1208549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1209549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1210549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1211549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
12122d61bbb3SSatish Balay         /* free up old matrix storage */
1213606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1214606d414cSSatish Balay         if (!a->singlemalloc) {
1215606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1216606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1217606d414cSSatish Balay         }
12182d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1219c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
12202d61bbb3SSatish Balay 
12212d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
12222d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1223b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
12242d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
12252d61bbb3SSatish Balay         a->reallocs++;
12262d61bbb3SSatish Balay         a->nz++;
12272d61bbb3SSatish Balay       }
12282d61bbb3SSatish Balay       N = nrow++ - 1;
12292d61bbb3SSatish Balay       /* shift up all the later entries in this row */
12302d61bbb3SSatish Balay       for (ii=N; ii>=i; ii--) {
12312d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1232549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
12332d61bbb3SSatish Balay       }
1234549d3d68SSatish Balay       if (N>=i) {
1235549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1236549d3d68SSatish Balay       }
12372d61bbb3SSatish Balay       rp[i]                      = bcol;
12382d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
12392d61bbb3SSatish Balay       noinsert1:;
12402d61bbb3SSatish Balay       low = i;
12412d61bbb3SSatish Balay     }
12422d61bbb3SSatish Balay     ailen[brow] = nrow;
12432d61bbb3SSatish Balay   }
12442d61bbb3SSatish Balay   PetscFunctionReturn(0);
12452d61bbb3SSatish Balay }
12462d61bbb3SSatish Balay 
12472d61bbb3SSatish Balay 
12484a2ae208SSatish Balay #undef __FUNCT__
12494a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ"
1250b380c88cSHong Zhang int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
12512d61bbb3SSatish Balay {
12522d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
12532d61bbb3SSatish Balay   Mat         outA;
12542d61bbb3SSatish Balay   int         ierr;
1255667159a5SBarry Smith   PetscTruth  row_identity,col_identity;
12562d61bbb3SSatish Balay 
12572d61bbb3SSatish Balay   PetscFunctionBegin;
1258d3d32019SBarry Smith   if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
1259667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1260667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1261667159a5SBarry Smith   if (!row_identity || !col_identity) {
126229bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1263667159a5SBarry Smith   }
12642d61bbb3SSatish Balay 
12652d61bbb3SSatish Balay   outA          = inA;
12662d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
12672d61bbb3SSatish Balay 
12682d61bbb3SSatish Balay   if (!a->diag) {
1269c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
12702d61bbb3SSatish Balay   }
1271cf242676SKris Buschelman 
1272c38d4ed2SBarry Smith   a->row        = row;
1273c38d4ed2SBarry Smith   a->col        = col;
1274c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1275c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1276c38d4ed2SBarry Smith 
1277c38d4ed2SBarry Smith   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
12784c49b128SBarry Smith   ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1279b0a32e0cSBarry Smith   PetscLogObjectParent(inA,a->icol);
1280c38d4ed2SBarry Smith 
1281cf242676SKris Buschelman   /*
1282cf242676SKris Buschelman       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1283cf242676SKris Buschelman       for ILU(0) factorization with natural ordering
1284cf242676SKris Buschelman   */
1285cf242676SKris Buschelman   if (a->bs < 8) {
1286cf242676SKris Buschelman     ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr);
1287cf242676SKris Buschelman   } else {
1288c38d4ed2SBarry Smith     if (!a->solve_work) {
128987828ca2SBarry Smith       ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
129087828ca2SBarry Smith       PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar));
1291c38d4ed2SBarry Smith     }
12922d61bbb3SSatish Balay   }
12932d61bbb3SSatish Balay 
1294667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1295667159a5SBarry Smith 
12962d61bbb3SSatish Balay   PetscFunctionReturn(0);
12972d61bbb3SSatish Balay }
12984a2ae208SSatish Balay #undef __FUNCT__
12994a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ"
1300ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1301ba4ca20aSSatish Balay {
1302c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1303ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1304d132466eSBarry Smith   int               ierr;
1305ba4ca20aSSatish Balay 
13063a40ed3dSBarry Smith   PetscFunctionBegin;
1307c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1308d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1309d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
13103a40ed3dSBarry Smith   PetscFunctionReturn(0);
1311ba4ca20aSSatish Balay }
1312d9b7c43dSSatish Balay 
1313fb2e594dSBarry Smith EXTERN_C_BEGIN
13144a2ae208SSatish Balay #undef __FUNCT__
13154a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
131627a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
131727a8da17SBarry Smith {
131827a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
131914db4f74SSatish Balay   int         i,nz,nbs;
132027a8da17SBarry Smith 
132127a8da17SBarry Smith   PetscFunctionBegin;
132214db4f74SSatish Balay   nz  = baij->maxnz/baij->bs2;
132314db4f74SSatish Balay   nbs = baij->nbs;
132427a8da17SBarry Smith   for (i=0; i<nz; i++) {
132527a8da17SBarry Smith     baij->j[i] = indices[i];
132627a8da17SBarry Smith   }
132727a8da17SBarry Smith   baij->nz = nz;
132814db4f74SSatish Balay   for (i=0; i<nbs; i++) {
132927a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
133027a8da17SBarry Smith   }
133127a8da17SBarry Smith 
133227a8da17SBarry Smith   PetscFunctionReturn(0);
133327a8da17SBarry Smith }
1334fb2e594dSBarry Smith EXTERN_C_END
133527a8da17SBarry Smith 
13364a2ae208SSatish Balay #undef __FUNCT__
13374a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices"
133827a8da17SBarry Smith /*@
133927a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
134027a8da17SBarry Smith        in the matrix.
134127a8da17SBarry Smith 
134227a8da17SBarry Smith   Input Parameters:
134327a8da17SBarry Smith +  mat - the SeqBAIJ matrix
134427a8da17SBarry Smith -  indices - the column indices
134527a8da17SBarry Smith 
134615091d37SBarry Smith   Level: advanced
134715091d37SBarry Smith 
134827a8da17SBarry Smith   Notes:
134927a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
135027a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
135127a8da17SBarry Smith   of the MatSetValues() operation.
135227a8da17SBarry Smith 
135327a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
135427a8da17SBarry Smith   MatCreateSeqBAIJ().
135527a8da17SBarry Smith 
135627a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
135727a8da17SBarry Smith 
135827a8da17SBarry Smith @*/
135927a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
136027a8da17SBarry Smith {
136127a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
136227a8da17SBarry Smith 
136327a8da17SBarry Smith   PetscFunctionBegin;
136427a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1365c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
136627a8da17SBarry Smith   if (f) {
136727a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
136827a8da17SBarry Smith   } else {
136929bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
137027a8da17SBarry Smith   }
137127a8da17SBarry Smith   PetscFunctionReturn(0);
137227a8da17SBarry Smith }
137327a8da17SBarry Smith 
13744a2ae208SSatish Balay #undef __FUNCT__
13754a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ"
1376273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v)
1377273d9f13SBarry Smith {
1378273d9f13SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1379273d9f13SBarry Smith   int          ierr,i,j,n,row,bs,*ai,*aj,mbs;
1380273d9f13SBarry Smith   PetscReal    atmp;
138187828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
1382273d9f13SBarry Smith   MatScalar    *aa;
1383273d9f13SBarry Smith   int          ncols,brow,krow,kcol;
1384273d9f13SBarry Smith 
1385273d9f13SBarry Smith   PetscFunctionBegin;
1386273d9f13SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1387273d9f13SBarry Smith   bs   = a->bs;
1388273d9f13SBarry Smith   aa   = a->a;
1389273d9f13SBarry Smith   ai   = a->i;
1390273d9f13SBarry Smith   aj   = a->j;
1391273d9f13SBarry Smith   mbs = a->mbs;
1392273d9f13SBarry Smith 
1393273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1394273d9f13SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1395273d9f13SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1396273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1397273d9f13SBarry Smith   for (i=0; i<mbs; i++) {
1398273d9f13SBarry Smith     ncols = ai[1] - ai[0]; ai++;
1399273d9f13SBarry Smith     brow  = bs*i;
1400273d9f13SBarry Smith     for (j=0; j<ncols; j++){
1401273d9f13SBarry Smith       /* bcol = bs*(*aj); */
1402273d9f13SBarry Smith       for (kcol=0; kcol<bs; kcol++){
1403273d9f13SBarry Smith         for (krow=0; krow<bs; krow++){
1404273d9f13SBarry Smith           atmp = PetscAbsScalar(*aa); aa++;
1405273d9f13SBarry Smith           row = brow + krow;    /* row index */
1406273d9f13SBarry Smith           /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */
1407273d9f13SBarry Smith           if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp;
1408273d9f13SBarry Smith         }
1409273d9f13SBarry Smith       }
1410273d9f13SBarry Smith       aj++;
1411273d9f13SBarry Smith     }
1412273d9f13SBarry Smith   }
1413273d9f13SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1414273d9f13SBarry Smith   PetscFunctionReturn(0);
1415273d9f13SBarry Smith }
1416273d9f13SBarry Smith 
14174a2ae208SSatish Balay #undef __FUNCT__
14184a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ"
1419273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A)
1420273d9f13SBarry Smith {
1421273d9f13SBarry Smith   int        ierr;
1422273d9f13SBarry Smith 
1423273d9f13SBarry Smith   PetscFunctionBegin;
1424273d9f13SBarry Smith   ierr =  MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1425273d9f13SBarry Smith   PetscFunctionReturn(0);
1426273d9f13SBarry Smith }
1427273d9f13SBarry Smith 
14284a2ae208SSatish Balay #undef __FUNCT__
14294a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ"
143087828ca2SBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar **array)
1431f2a5309cSSatish Balay {
1432f2a5309cSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1433f2a5309cSSatish Balay   PetscFunctionBegin;
1434f2a5309cSSatish Balay   *array = a->a;
1435f2a5309cSSatish Balay   PetscFunctionReturn(0);
1436f2a5309cSSatish Balay }
1437f2a5309cSSatish Balay 
14384a2ae208SSatish Balay #undef __FUNCT__
14394a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ"
144087828ca2SBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar **array)
1441f2a5309cSSatish Balay {
1442f2a5309cSSatish Balay   PetscFunctionBegin;
1443f2a5309cSSatish Balay   PetscFunctionReturn(0);
1444f2a5309cSSatish Balay }
1445f2a5309cSSatish Balay 
144642ee4b1aSHong Zhang #include "petscblaslapack.h"
144742ee4b1aSHong Zhang #undef __FUNCT__
144842ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqBAIJ"
144942ee4b1aSHong Zhang int MatAXPY_SeqBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
145042ee4b1aSHong Zhang {
145142ee4b1aSHong Zhang   Mat_SeqBAIJ  *x  = (Mat_SeqBAIJ *)X->data,*y = (Mat_SeqBAIJ *)Y->data;
1452*6550863bSHong Zhang   int          ierr,one=1,i,bs=y->bs,j,bs2;
145342ee4b1aSHong Zhang 
145442ee4b1aSHong Zhang   PetscFunctionBegin;
145542ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
145642ee4b1aSHong Zhang     BLaxpy_(&x->nz,a,x->a,&one,y->a,&one);
1457c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1458c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1459c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1460c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1461c537a176SHong Zhang     }
1462c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1463c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1464c4319e64SHong Zhang       y->XtoY = X;
1465c537a176SHong Zhang     }
1466c4319e64SHong Zhang     bs2 = bs*bs;
1467c537a176SHong Zhang     for (i=0; i<x->nz; i++) {
1468c4319e64SHong Zhang       j = 0;
1469c4319e64SHong Zhang       while (j < bs2){
1470*6550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1471c4319e64SHong Zhang         j++;
1472c537a176SHong Zhang       }
1473c4319e64SHong Zhang     }
1474c4319e64SHong Zhang     PetscLogInfo(0,"MatAXPY_SeqBAIJ: ratio of nnz(X)/nnz(Y): %d/%d = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));
147542ee4b1aSHong Zhang   } else {
147642ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
147742ee4b1aSHong Zhang   }
147842ee4b1aSHong Zhang   PetscFunctionReturn(0);
147942ee4b1aSHong Zhang }
148042ee4b1aSHong Zhang 
14812593348eSBarry Smith /* -------------------------------------------------------------------*/
1482cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1483cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1484cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1485cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1486cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
14877c922b88SBarry Smith        MatMultTranspose_SeqBAIJ,
14887c922b88SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1489cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1490cc2dc46cSBarry Smith        0,
1491cc2dc46cSBarry Smith        0,
1492cc2dc46cSBarry Smith        0,
1493cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1494cc2dc46cSBarry Smith        0,
1495b6490206SBarry Smith        0,
1496f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1497cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1498cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1499cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1500cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1501cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1502b6490206SBarry Smith        0,
1503cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1504cc2dc46cSBarry Smith        0,
1505cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1506cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1507cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1508cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1509cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1510cc2dc46cSBarry Smith        0,
1511cc2dc46cSBarry Smith        0,
1512273d9f13SBarry Smith        MatSetUpPreallocation_SeqBAIJ,
1513cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1514cc2dc46cSBarry Smith        0,
1515f2a5309cSSatish Balay        MatGetArray_SeqBAIJ,
1516f2a5309cSSatish Balay        MatRestoreArray_SeqBAIJ,
15172e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1518cc2dc46cSBarry Smith        0,
1519cc2dc46cSBarry Smith        0,
1520cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1521cc2dc46cSBarry Smith        0,
152242ee4b1aSHong Zhang        MatAXPY_SeqBAIJ,
1523cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1524cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1525cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1526cc2dc46cSBarry Smith        0,
1527cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1528cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1529cc2dc46cSBarry Smith        0,
1530cc2dc46cSBarry Smith        0,
1531cc2dc46cSBarry Smith        0,
1532cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
15333b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
153492c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1535cc2dc46cSBarry Smith        0,
1536cc2dc46cSBarry Smith        0,
1537cc2dc46cSBarry Smith        0,
1538cc2dc46cSBarry Smith        0,
1539cc2dc46cSBarry Smith        0,
1540cc2dc46cSBarry Smith        0,
1541d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1542cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1543b9b97703SBarry Smith        MatDestroy_SeqBAIJ,
1544b9b97703SBarry Smith        MatView_SeqBAIJ,
15453a64cc32SBarry Smith        MatGetPetscMaps_Petsc,
1546273d9f13SBarry Smith        0,
1547273d9f13SBarry Smith        0,
1548273d9f13SBarry Smith        0,
1549273d9f13SBarry Smith        0,
1550273d9f13SBarry Smith        0,
1551273d9f13SBarry Smith        0,
1552273d9f13SBarry Smith        MatGetRowMax_SeqBAIJ,
1553273d9f13SBarry Smith        MatConvert_Basic};
15542593348eSBarry Smith 
15553e90b805SBarry Smith EXTERN_C_BEGIN
15564a2ae208SSatish Balay #undef __FUNCT__
15574a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ"
15583e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
15593e90b805SBarry Smith {
15603e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1561273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1562d132466eSBarry Smith   int         ierr;
15633e90b805SBarry Smith 
15643e90b805SBarry Smith   PetscFunctionBegin;
15653e90b805SBarry Smith   if (aij->nonew != 1) {
156629bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15673e90b805SBarry Smith   }
15683e90b805SBarry Smith 
15693e90b805SBarry Smith   /* allocate space for values if not already there */
15703e90b805SBarry Smith   if (!aij->saved_values) {
157187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
15723e90b805SBarry Smith   }
15733e90b805SBarry Smith 
15743e90b805SBarry Smith   /* copy values over */
157587828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
15763e90b805SBarry Smith   PetscFunctionReturn(0);
15773e90b805SBarry Smith }
15783e90b805SBarry Smith EXTERN_C_END
15793e90b805SBarry Smith 
15803e90b805SBarry Smith EXTERN_C_BEGIN
15814a2ae208SSatish Balay #undef __FUNCT__
15824a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ"
158332e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
15843e90b805SBarry Smith {
15853e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1586273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
15873e90b805SBarry Smith 
15883e90b805SBarry Smith   PetscFunctionBegin;
15893e90b805SBarry Smith   if (aij->nonew != 1) {
159029bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15913e90b805SBarry Smith   }
15923e90b805SBarry Smith   if (!aij->saved_values) {
159329bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
15943e90b805SBarry Smith   }
15953e90b805SBarry Smith 
15963e90b805SBarry Smith   /* copy values over */
159787828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
15983e90b805SBarry Smith   PetscFunctionReturn(0);
15993e90b805SBarry Smith }
16003e90b805SBarry Smith EXTERN_C_END
16013e90b805SBarry Smith 
1602273d9f13SBarry Smith EXTERN_C_BEGIN
1603273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*);
1604273d9f13SBarry Smith EXTERN_C_END
1605273d9f13SBarry Smith 
1606273d9f13SBarry Smith EXTERN_C_BEGIN
16074a2ae208SSatish Balay #undef __FUNCT__
16084a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ"
1609273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B)
16102593348eSBarry Smith {
1611273d9f13SBarry Smith   int         ierr,size;
1612b6490206SBarry Smith   Mat_SeqBAIJ *b;
16133b2fbd54SBarry Smith 
16143a40ed3dSBarry Smith   PetscFunctionBegin;
1615273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
161629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1617b6490206SBarry Smith 
1618273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1619273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1620b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr);
1621b0a32e0cSBarry Smith   B->data = (void*)b;
1622549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1623549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
16242593348eSBarry Smith   B->factor           = 0;
16252593348eSBarry Smith   B->lupivotthreshold = 1.0;
162690f02eecSBarry Smith   B->mapping          = 0;
16272593348eSBarry Smith   b->row              = 0;
16282593348eSBarry Smith   b->col              = 0;
1629e51c0b9cSSatish Balay   b->icol             = 0;
16302593348eSBarry Smith   b->reallocs         = 0;
16313e90b805SBarry Smith   b->saved_values     = 0;
1632cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
1633563b5814SBarry Smith   b->setvalueslen     = 0;
1634563b5814SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1635563b5814SBarry Smith #endif
16362593348eSBarry Smith 
16373a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
16383a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1639a5ae1ecdSBarry Smith 
1640c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1641c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
16422593348eSBarry Smith   b->nonew            = 0;
16432593348eSBarry Smith   b->diag             = 0;
16442593348eSBarry Smith   b->solve_work       = 0;
1645de6a44a3SBarry Smith   b->mult_work        = 0;
16462a1b7f2aSHong Zhang   B->spptr            = 0;
16470e6d2581SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1648883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
1649c4319e64SHong Zhang   b->xtoy              = 0;
1650c4319e64SHong Zhang   b->XtoY              = 0;
16514e220ebcSLois Curfman McInnes 
1652f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
16533e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1654bc4b532fSSatish Balay                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1655f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
16563e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1657bc4b532fSSatish Balay                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1658f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
165927a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1660bc4b532fSSatish Balay                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1661a6175056SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqaij_C",
1662273d9f13SBarry Smith                                      "MatConvert_SeqBAIJ_SeqAIJ",
1663273d9f13SBarry Smith                                       MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr);
16643a40ed3dSBarry Smith   PetscFunctionReturn(0);
16652593348eSBarry Smith }
1666273d9f13SBarry Smith EXTERN_C_END
16672593348eSBarry Smith 
16684a2ae208SSatish Balay #undef __FUNCT__
16694a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ"
16702e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
16712593348eSBarry Smith {
16722593348eSBarry Smith   Mat         C;
1673b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
167427a8da17SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1675de6a44a3SBarry Smith 
16763a40ed3dSBarry Smith   PetscFunctionBegin;
167729bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
16782593348eSBarry Smith 
16792593348eSBarry Smith   *B = 0;
1680273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1681273d9f13SBarry Smith   ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
1682273d9f13SBarry Smith   c    = (Mat_SeqBAIJ*)C->data;
168344cd7ae7SLois Curfman McInnes 
1684b6490206SBarry Smith   c->bs         = a->bs;
16859df24120SSatish Balay   c->bs2        = a->bs2;
1686b6490206SBarry Smith   c->mbs        = a->mbs;
1687b6490206SBarry Smith   c->nbs        = a->nbs;
168835d8aa7fSBarry Smith   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
16892593348eSBarry Smith 
1690b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1691b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
1692b6490206SBarry Smith   for (i=0; i<mbs; i++) {
16932593348eSBarry Smith     c->imax[i] = a->imax[i];
16942593348eSBarry Smith     c->ilen[i] = a->ilen[i];
16952593348eSBarry Smith   }
16962593348eSBarry Smith 
16972593348eSBarry Smith   /* allocate the matrix space */
1698c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
16993f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1700b0a32e0cSBarry Smith   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
17017e67e3f9SSatish Balay   c->j = (int*)(c->a + nz*bs2);
1702de6a44a3SBarry Smith   c->i = c->j + nz;
1703549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1704b6490206SBarry Smith   if (mbs > 0) {
1705549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
17062e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1707549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
17082e8a6d31SBarry Smith     } else {
1709549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
17102593348eSBarry Smith     }
17112593348eSBarry Smith   }
17122593348eSBarry Smith 
1713b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
17142593348eSBarry Smith   c->sorted      = a->sorted;
17152593348eSBarry Smith   c->roworiented = a->roworiented;
17162593348eSBarry Smith   c->nonew       = a->nonew;
17172593348eSBarry Smith 
17182593348eSBarry Smith   if (a->diag) {
1719b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1720b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1721b6490206SBarry Smith     for (i=0; i<mbs; i++) {
17222593348eSBarry Smith       c->diag[i] = a->diag[i];
17232593348eSBarry Smith     }
172498305bb5SBarry Smith   } else c->diag        = 0;
17252593348eSBarry Smith   c->nz                 = a->nz;
17262593348eSBarry Smith   c->maxnz              = a->maxnz;
17272593348eSBarry Smith   c->solve_work         = 0;
17282a1b7f2aSHong Zhang   C->spptr              = 0;     /* Dangerous -I'm throwing away a->spptr */
17297fc0212eSBarry Smith   c->mult_work          = 0;
1730273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
1731273d9f13SBarry Smith   C->assembled          = PETSC_TRUE;
17322593348eSBarry Smith   *B = C;
1733b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
17343a40ed3dSBarry Smith   PetscFunctionReturn(0);
17352593348eSBarry Smith }
17362593348eSBarry Smith 
1737273d9f13SBarry Smith EXTERN_C_BEGIN
17384a2ae208SSatish Balay #undef __FUNCT__
17394a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ"
1740b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A)
17412593348eSBarry Smith {
1742b6490206SBarry Smith   Mat_SeqBAIJ  *a;
17432593348eSBarry Smith   Mat          B;
1744f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1745b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
174635aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1747a7c10996SSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
174887828ca2SBarry Smith   PetscScalar  *aa;
174919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1750ecc4ab6dSBarry Smith #if defined(PETSC_HAVE_DSCPACK)
1751ecc4ab6dSBarry Smith   PetscTruth   flag;
1752ecc4ab6dSBarry Smith #endif
17532593348eSBarry Smith 
17543a40ed3dSBarry Smith   PetscFunctionBegin;
1755b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1756de6a44a3SBarry Smith   bs2  = bs*bs;
1757b6490206SBarry Smith 
1758d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
175929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1760b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
17610752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1762552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
17632593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
17642593348eSBarry Smith 
1765d64ed03dSBarry Smith   if (header[3] < 0) {
176629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ");
1767d64ed03dSBarry Smith   }
1768d64ed03dSBarry Smith 
176929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
177035aab85fSBarry Smith 
177135aab85fSBarry Smith   /*
177235aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
177335aab85fSBarry Smith     divisible by the blocksize
177435aab85fSBarry Smith   */
1775b6490206SBarry Smith   mbs        = M/bs;
177635aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
177735aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
177835aab85fSBarry Smith   else                  mbs++;
177935aab85fSBarry Smith   if (extra_rows) {
1780b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
178135aab85fSBarry Smith   }
1782b6490206SBarry Smith 
17832593348eSBarry Smith   /* read in row lengths */
1784b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
17850752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
178635aab85fSBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
17872593348eSBarry Smith 
1788b6490206SBarry Smith   /* read in column indices */
1789b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
17900752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
179135aab85fSBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1792b6490206SBarry Smith 
1793b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1794b0a32e0cSBarry Smith   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
1795549d3d68SSatish Balay   ierr     = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1796b0a32e0cSBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
1797549d3d68SSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
179835aab85fSBarry Smith   masked   = mask + mbs;
1799b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1800b6490206SBarry Smith   for (i=0; i<mbs; i++) {
180135aab85fSBarry Smith     nmask = 0;
1802b6490206SBarry Smith     for (j=0; j<bs; j++) {
1803b6490206SBarry Smith       kmax = rowlengths[rowcount];
1804b6490206SBarry Smith       for (k=0; k<kmax; k++) {
180535aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
180635aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1807b6490206SBarry Smith       }
1808b6490206SBarry Smith       rowcount++;
1809b6490206SBarry Smith     }
181035aab85fSBarry Smith     browlengths[i] += nmask;
181135aab85fSBarry Smith     /* zero out the mask elements we set */
181235aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1813b6490206SBarry Smith   }
1814b6490206SBarry Smith 
18152593348eSBarry Smith   /* create our matrix */
18163eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
18172593348eSBarry Smith   B = *A;
1818b6490206SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
18192593348eSBarry Smith 
18202593348eSBarry Smith   /* set matrix "i" values */
1821de6a44a3SBarry Smith   a->i[0] = 0;
1822b6490206SBarry Smith   for (i=1; i<= mbs; i++) {
1823b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1824b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
18252593348eSBarry Smith   }
1826b6490206SBarry Smith   a->nz         = 0;
1827b6490206SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
18282593348eSBarry Smith 
1829b6490206SBarry Smith   /* read in nonzero values */
183087828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
18310752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
183235aab85fSBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1833b6490206SBarry Smith 
1834b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1835b6490206SBarry Smith   nzcount = 0; jcount = 0;
1836b6490206SBarry Smith   for (i=0; i<mbs; i++) {
1837b6490206SBarry Smith     nzcountb = nzcount;
183835aab85fSBarry Smith     nmask    = 0;
1839b6490206SBarry Smith     for (j=0; j<bs; j++) {
1840b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1841b6490206SBarry Smith       for (k=0; k<kmax; k++) {
184235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
184335aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1844b6490206SBarry Smith       }
1845b6490206SBarry Smith     }
1846de6a44a3SBarry Smith     /* sort the masked values */
1847433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1848de6a44a3SBarry Smith 
1849b6490206SBarry Smith     /* set "j" values into matrix */
1850b6490206SBarry Smith     maskcount = 1;
185135aab85fSBarry Smith     for (j=0; j<nmask; j++) {
185235aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1853de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1854b6490206SBarry Smith     }
1855b6490206SBarry Smith     /* set "a" values into matrix */
1856de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1857b6490206SBarry Smith     for (j=0; j<bs; j++) {
1858b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1859b6490206SBarry Smith       for (k=0; k<kmax; k++) {
1860de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1861de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1862de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1863de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1864375fe846SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1865b6490206SBarry Smith       }
1866b6490206SBarry Smith     }
186735aab85fSBarry Smith     /* zero out the mask elements we set */
186835aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1869b6490206SBarry Smith   }
187029bbc08cSBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
1871b6490206SBarry Smith 
1872606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1873606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1874606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
1875606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
1876606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
1877b6490206SBarry Smith 
1878b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1879de6a44a3SBarry Smith 
1880ecc4ab6dSBarry Smith #if defined(PETSC_HAVE_DSCPACK)
1881ecc4ab6dSBarry Smith   ierr = PetscOptionsHasName(B->prefix,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
1882ecc4ab6dSBarry Smith   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(B);CHKERRQ(ierr); }
1883ecc4ab6dSBarry Smith #endif
1884ecc4ab6dSBarry Smith 
18859c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
18863a40ed3dSBarry Smith   PetscFunctionReturn(0);
18872593348eSBarry Smith }
1888273d9f13SBarry Smith EXTERN_C_END
18892593348eSBarry Smith 
18904a2ae208SSatish Balay #undef __FUNCT__
18914a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ"
1892273d9f13SBarry Smith /*@C
1893273d9f13SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1894273d9f13SBarry Smith    compressed row) format.  For good matrix assembly performance the
1895273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1896273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1897273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
18982593348eSBarry Smith 
1899273d9f13SBarry Smith    Collective on MPI_Comm
1900273d9f13SBarry Smith 
1901273d9f13SBarry Smith    Input Parameters:
1902273d9f13SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
1903273d9f13SBarry Smith .  bs - size of block
1904273d9f13SBarry Smith .  m - number of rows
1905273d9f13SBarry Smith .  n - number of columns
190635d8aa7fSBarry Smith .  nz - number of nonzero blocks  per block row (same for all rows)
190735d8aa7fSBarry Smith -  nnz - array containing the number of nonzero blocks in the various block rows
1908273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1909273d9f13SBarry Smith 
1910273d9f13SBarry Smith    Output Parameter:
1911273d9f13SBarry Smith .  A - the matrix
1912273d9f13SBarry Smith 
1913273d9f13SBarry Smith    Options Database Keys:
1914273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1915273d9f13SBarry Smith                      block calculations (much slower)
1916273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1917273d9f13SBarry Smith 
1918273d9f13SBarry Smith    Level: intermediate
1919273d9f13SBarry Smith 
1920273d9f13SBarry Smith    Notes:
192135d8aa7fSBarry Smith    A nonzero block is any block that as 1 or more nonzeros in it
192235d8aa7fSBarry Smith 
1923273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1924273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1925273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1926273d9f13SBarry Smith 
1927273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1928273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1929273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1930273d9f13SBarry Smith    matrices.
1931273d9f13SBarry Smith 
1932273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1933273d9f13SBarry Smith @*/
1934273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1935273d9f13SBarry Smith {
1936273d9f13SBarry Smith   int ierr;
1937273d9f13SBarry Smith 
1938273d9f13SBarry Smith   PetscFunctionBegin;
1939273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1940273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr);
1941273d9f13SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1942273d9f13SBarry Smith   PetscFunctionReturn(0);
1943273d9f13SBarry Smith }
1944273d9f13SBarry Smith 
19454a2ae208SSatish Balay #undef __FUNCT__
19464a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation"
1947273d9f13SBarry Smith /*@C
1948273d9f13SBarry Smith    MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros
1949273d9f13SBarry Smith    per row in the matrix. For good matrix assembly performance the
1950273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1951273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1952273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1953273d9f13SBarry Smith 
1954273d9f13SBarry Smith    Collective on MPI_Comm
1955273d9f13SBarry Smith 
1956273d9f13SBarry Smith    Input Parameters:
1957273d9f13SBarry Smith +  A - the matrix
1958273d9f13SBarry Smith .  bs - size of block
1959273d9f13SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1960273d9f13SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
1961273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1962273d9f13SBarry Smith 
1963273d9f13SBarry Smith    Options Database Keys:
1964273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1965273d9f13SBarry Smith                      block calculations (much slower)
1966273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1967273d9f13SBarry Smith 
1968273d9f13SBarry Smith    Level: intermediate
1969273d9f13SBarry Smith 
1970273d9f13SBarry Smith    Notes:
1971273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1972273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1973273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1974273d9f13SBarry Smith 
1975273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1976273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1977273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1978273d9f13SBarry Smith    matrices.
1979273d9f13SBarry Smith 
1980273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1981273d9f13SBarry Smith @*/
1982273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
1983273d9f13SBarry Smith {
1984273d9f13SBarry Smith   Mat_SeqBAIJ *b;
1985273d9f13SBarry Smith   int         i,len,ierr,mbs,nbs,bs2,newbs = bs;
1986273d9f13SBarry Smith   PetscTruth  flg;
1987273d9f13SBarry Smith 
1988273d9f13SBarry Smith   PetscFunctionBegin;
1989273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr);
1990273d9f13SBarry Smith   if (!flg) PetscFunctionReturn(0);
1991273d9f13SBarry Smith 
1992273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1993b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1994273d9f13SBarry Smith   if (nnz && newbs != bs) {
1995273d9f13SBarry Smith     SETERRQ(1,"Cannot change blocksize from command line if setting nnz");
1996273d9f13SBarry Smith   }
1997273d9f13SBarry Smith   bs   = newbs;
1998273d9f13SBarry Smith 
1999273d9f13SBarry Smith   mbs  = B->m/bs;
2000273d9f13SBarry Smith   nbs  = B->n/bs;
2001273d9f13SBarry Smith   bs2  = bs*bs;
2002273d9f13SBarry Smith 
2003273d9f13SBarry Smith   if (mbs*bs!=B->m || nbs*bs!=B->n) {
200435d8aa7fSBarry Smith     SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs);
2005273d9f13SBarry Smith   }
2006273d9f13SBarry Smith 
2007435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2008435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
2009273d9f13SBarry Smith   if (nnz) {
2010273d9f13SBarry Smith     for (i=0; i<mbs; i++) {
2011273d9f13SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
20123a7fca6bSBarry Smith       if (nnz[i] > nbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %d value %d rowlength %d",i,nnz[i],nbs);
2013273d9f13SBarry Smith     }
2014273d9f13SBarry Smith   }
2015273d9f13SBarry Smith 
2016273d9f13SBarry Smith   b       = (Mat_SeqBAIJ*)B->data;
2017b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
201854138f6bSKris Buschelman   B->ops->solve               = MatSolve_SeqBAIJ_Update;
201954138f6bSKris Buschelman   B->ops->solvetranspose      = MatSolveTranspose_SeqBAIJ_Update;
2020273d9f13SBarry Smith   if (!flg) {
2021273d9f13SBarry Smith     switch (bs) {
2022273d9f13SBarry Smith     case 1:
2023273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
2024273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
2025273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
2026273d9f13SBarry Smith       break;
2027273d9f13SBarry Smith     case 2:
2028273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
2029273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
2030273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
2031273d9f13SBarry Smith       break;
2032273d9f13SBarry Smith     case 3:
2033273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
2034273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
2035273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
2036273d9f13SBarry Smith       break;
2037273d9f13SBarry Smith     case 4:
2038273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
2039273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
2040273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
2041273d9f13SBarry Smith       break;
2042273d9f13SBarry Smith     case 5:
2043273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
2044273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
2045273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
2046273d9f13SBarry Smith       break;
2047273d9f13SBarry Smith     case 6:
2048273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
2049273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
2050273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
2051273d9f13SBarry Smith       break;
2052273d9f13SBarry Smith     case 7:
205354138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
2054273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
2055273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
2056273d9f13SBarry Smith       break;
2057273d9f13SBarry Smith     default:
205854138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
2059273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_N;
2060273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_N;
2061273d9f13SBarry Smith       break;
2062273d9f13SBarry Smith     }
2063273d9f13SBarry Smith   }
2064273d9f13SBarry Smith   b->bs      = bs;
2065273d9f13SBarry Smith   b->mbs     = mbs;
2066273d9f13SBarry Smith   b->nbs     = nbs;
2067b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
2068273d9f13SBarry Smith   if (!nnz) {
2069435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2070273d9f13SBarry Smith     else if (nz <= 0)        nz = 1;
2071273d9f13SBarry Smith     for (i=0; i<mbs; i++) b->imax[i] = nz;
2072273d9f13SBarry Smith     nz = nz*mbs;
2073273d9f13SBarry Smith   } else {
2074273d9f13SBarry Smith     nz = 0;
2075273d9f13SBarry Smith     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2076273d9f13SBarry Smith   }
2077273d9f13SBarry Smith 
2078273d9f13SBarry Smith   /* allocate the matrix space */
2079273d9f13SBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
2080b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&b->a);CHKERRQ(ierr);
2081273d9f13SBarry Smith   ierr  = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
2082273d9f13SBarry Smith   b->j  = (int*)(b->a + nz*bs2);
2083273d9f13SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
2084273d9f13SBarry Smith   b->i  = b->j + nz;
2085273d9f13SBarry Smith   b->singlemalloc = PETSC_TRUE;
2086273d9f13SBarry Smith 
2087273d9f13SBarry Smith   b->i[0] = 0;
2088273d9f13SBarry Smith   for (i=1; i<mbs+1; i++) {
2089273d9f13SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
2090273d9f13SBarry Smith   }
2091273d9f13SBarry Smith 
2092273d9f13SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
209316d6aa21SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
2094b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
2095273d9f13SBarry Smith   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
2096273d9f13SBarry Smith 
2097273d9f13SBarry Smith   b->bs               = bs;
2098273d9f13SBarry Smith   b->bs2              = bs2;
2099273d9f13SBarry Smith   b->mbs              = mbs;
2100273d9f13SBarry Smith   b->nz               = 0;
2101273d9f13SBarry Smith   b->maxnz            = nz*bs2;
2102273d9f13SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
2103273d9f13SBarry Smith   PetscFunctionReturn(0);
2104273d9f13SBarry Smith }
21052593348eSBarry Smith 
2106