xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 4e7234bf172ce7019c4948834d73bd3002d026ba)
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_"
24f15d580aSBarry Smith void matsetvaluesblocked4_(Mat *AA,int *mm,const int im[],int *nn,const int in[],const 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;
31f15d580aSBarry Smith   const 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)
158f15d580aSBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,const int[],int,const int[],const 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"
219*4e7234bfSBarry 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"
244*4e7234bfSBarry 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"
742f15d580aSBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,const int im[],int n,const 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"
788f15d580aSBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,const int im[],int n,const int in[],const 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"
812f15d580aSBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,const int im[],int n,const int in[],const 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;
819f15d580aSBarry Smith   const MatScalar   *value = v;
820f15d580aSBarry Smith   MatScalar         *ap,*aa = a->a,*bap;
82192c4ed94SBarry Smith 
8223a40ed3dSBarry Smith   PetscFunctionBegin;
8230e324ae4SSatish Balay   if (roworiented) {
8240e324ae4SSatish Balay     stepval = (n-1)*bs;
8250e324ae4SSatish Balay   } else {
8260e324ae4SSatish Balay     stepval = (m-1)*bs;
8270e324ae4SSatish Balay   }
82892c4ed94SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
82992c4ed94SBarry Smith     row  = im[k];
8305ef9f2a5SBarry Smith     if (row < 0) continue;
831aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
832590ac198SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1);
83392c4ed94SBarry Smith #endif
83492c4ed94SBarry Smith     rp   = aj + ai[row];
83592c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
83692c4ed94SBarry Smith     rmax = imax[row];
83792c4ed94SBarry Smith     nrow = ailen[row];
83892c4ed94SBarry Smith     low  = 0;
83992c4ed94SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
8405ef9f2a5SBarry Smith       if (in[l] < 0) continue;
841aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
842590ac198SBarry Smith       if (in[l] >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],a->nbs-1);
84392c4ed94SBarry Smith #endif
84492c4ed94SBarry Smith       col = in[l];
84592c4ed94SBarry Smith       if (roworiented) {
8460e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
8470e324ae4SSatish Balay       } else {
8480e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
84992c4ed94SBarry Smith       }
85092c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
85192c4ed94SBarry Smith       while (high-low > 7) {
85292c4ed94SBarry Smith         t = (low+high)/2;
85392c4ed94SBarry Smith         if (rp[t] > col) high = t;
85492c4ed94SBarry Smith         else             low  = t;
85592c4ed94SBarry Smith       }
85692c4ed94SBarry Smith       for (i=low; i<high; i++) {
85792c4ed94SBarry Smith         if (rp[i] > col) break;
85892c4ed94SBarry Smith         if (rp[i] == col) {
8598a84c255SSatish Balay           bap  = ap +  bs2*i;
8600e324ae4SSatish Balay           if (roworiented) {
8618a84c255SSatish Balay             if (is == ADD_VALUES) {
862dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
863dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
8648a84c255SSatish Balay                   bap[jj] += *value++;
865dd9472c6SBarry Smith                 }
866dd9472c6SBarry Smith               }
8670e324ae4SSatish Balay             } else {
868dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
869dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
8700e324ae4SSatish Balay                   bap[jj] = *value++;
8718a84c255SSatish Balay                 }
872dd9472c6SBarry Smith               }
873dd9472c6SBarry Smith             }
8740e324ae4SSatish Balay           } else {
8750e324ae4SSatish Balay             if (is == ADD_VALUES) {
876dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
877dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
8780e324ae4SSatish Balay                   *bap++ += *value++;
879dd9472c6SBarry Smith                 }
880dd9472c6SBarry Smith               }
8810e324ae4SSatish Balay             } else {
882dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
883dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
8840e324ae4SSatish Balay                   *bap++  = *value++;
8850e324ae4SSatish Balay                 }
8868a84c255SSatish Balay               }
887dd9472c6SBarry Smith             }
888dd9472c6SBarry Smith           }
889f1241b54SBarry Smith           goto noinsert2;
89092c4ed94SBarry Smith         }
89192c4ed94SBarry Smith       }
89289280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
89329bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
89492c4ed94SBarry Smith       if (nrow >= rmax) {
89592c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
89692c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
8973f1db9ecSBarry Smith         MatScalar *new_a;
89892c4ed94SBarry Smith 
89929bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
90089280ab3SLois Curfman McInnes 
90192c4ed94SBarry Smith         /* malloc new storage space */
9023f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
903b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
90492c4ed94SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
90592c4ed94SBarry Smith         new_i   = new_j + new_nz;
90692c4ed94SBarry Smith 
90792c4ed94SBarry Smith         /* copy over old data into new slots */
90892c4ed94SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
90992c4ed94SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
910549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
91192c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
912549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
913549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
914549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
915549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
91692c4ed94SBarry Smith         /* free up old matrix storage */
917606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
918606d414cSSatish Balay         if (!a->singlemalloc) {
919606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
920606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
921606d414cSSatish Balay         }
92292c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
923c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
92492c4ed94SBarry Smith 
92592c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
92692c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
927b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
92892c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
92992c4ed94SBarry Smith         a->reallocs++;
93092c4ed94SBarry Smith         a->nz++;
93192c4ed94SBarry Smith       }
93292c4ed94SBarry Smith       N = nrow++ - 1;
93392c4ed94SBarry Smith       /* shift up all the later entries in this row */
93492c4ed94SBarry Smith       for (ii=N; ii>=i; ii--) {
93592c4ed94SBarry Smith         rp[ii+1] = rp[ii];
936549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
93792c4ed94SBarry Smith       }
938549d3d68SSatish Balay       if (N >= i) {
939549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
940549d3d68SSatish Balay       }
94192c4ed94SBarry Smith       rp[i] = col;
9428a84c255SSatish Balay       bap   = ap +  bs2*i;
9430e324ae4SSatish Balay       if (roworiented) {
944dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
945dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
9460e324ae4SSatish Balay             bap[jj] = *value++;
947dd9472c6SBarry Smith           }
948dd9472c6SBarry Smith         }
9490e324ae4SSatish Balay       } else {
950dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
951dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++) {
9520e324ae4SSatish Balay             *bap++  = *value++;
9530e324ae4SSatish Balay           }
954dd9472c6SBarry Smith         }
955dd9472c6SBarry Smith       }
956f1241b54SBarry Smith       noinsert2:;
95792c4ed94SBarry Smith       low = i;
95892c4ed94SBarry Smith     }
95992c4ed94SBarry Smith     ailen[row] = nrow;
96092c4ed94SBarry Smith   }
9613a40ed3dSBarry Smith   PetscFunctionReturn(0);
96292c4ed94SBarry Smith }
96392c4ed94SBarry Smith 
9644a2ae208SSatish Balay #undef __FUNCT__
9654a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ"
9668f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
967584200bdSSatish Balay {
968584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
969584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
970273d9f13SBarry Smith   int        m = A->m,*ip,N,*ailen = a->ilen;
971549d3d68SSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
9723f1db9ecSBarry Smith   MatScalar  *aa = a->a,*ap;
973a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
974a56a16c8SHong Zhang   PetscTruth   flag;
975a56a16c8SHong Zhang #endif
976584200bdSSatish Balay 
9773a40ed3dSBarry Smith   PetscFunctionBegin;
9783a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
979584200bdSSatish Balay 
98043ee02c3SBarry Smith   if (m) rmax = ailen[0];
981584200bdSSatish Balay   for (i=1; i<mbs; i++) {
982584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
983584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
984d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
985584200bdSSatish Balay     if (fshift) {
986a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
987584200bdSSatish Balay       N = ailen[i];
988584200bdSSatish Balay       for (j=0; j<N; j++) {
989584200bdSSatish Balay         ip[j-fshift] = ip[j];
990549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
991584200bdSSatish Balay       }
992584200bdSSatish Balay     }
993584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
994584200bdSSatish Balay   }
995584200bdSSatish Balay   if (mbs) {
996584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
997584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
998584200bdSSatish Balay   }
999584200bdSSatish Balay   /* reset ilen and imax for each row */
1000584200bdSSatish Balay   for (i=0; i<mbs; i++) {
1001584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
1002584200bdSSatish Balay   }
1003a7c10996SSatish Balay   a->nz = ai[mbs];
1004584200bdSSatish Balay 
1005584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
1006584200bdSSatish Balay   if (fshift && a->diag) {
1007606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
1008b424e231SHong Zhang     PetscLogObjectMemory(A,-(mbs+1)*sizeof(int));
1009584200bdSSatish Balay     a->diag = 0;
1010584200bdSSatish Balay   }
1011bba1ac68SSatish 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);
1012bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs);
1013b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
1014e2f3b5e9SSatish Balay   a->reallocs          = 0;
10150e6d2581SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
1016cf4441caSHong Zhang 
1017a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
10182c535e4dSHong Zhang   ierr = PetscOptionsHasName(A->prefix,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
1019a56a16c8SHong Zhang   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(A);CHKERRQ(ierr); }
1020a56a16c8SHong Zhang #endif
1021a56a16c8SHong Zhang 
10223a40ed3dSBarry Smith   PetscFunctionReturn(0);
1023584200bdSSatish Balay }
1024584200bdSSatish Balay 
10255a838052SSatish Balay 
1026bea157c4SSatish Balay 
1027bea157c4SSatish Balay /*
1028bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
1029bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
1030bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
1031bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
1032bea157c4SSatish Balay */
10334a2ae208SSatish Balay #undef __FUNCT__
10344a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks"
1035bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
1036d9b7c43dSSatish Balay {
1037bea157c4SSatish Balay   int        i,j,k,row;
1038bea157c4SSatish Balay   PetscTruth flg;
10393a40ed3dSBarry Smith 
1040433994e6SBarry Smith   PetscFunctionBegin;
1041bea157c4SSatish Balay   for (i=0,j=0; i<n; j++) {
1042bea157c4SSatish Balay     row = idx[i];
1043bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
1044bea157c4SSatish Balay       sizes[j] = 1;
1045bea157c4SSatish Balay       i++;
1046e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
1047bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
1048bea157c4SSatish Balay       i++;
1049bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
1050bea157c4SSatish Balay       flg = PETSC_TRUE;
1051bea157c4SSatish Balay       for (k=1; k<bs; k++) {
1052bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
1053bea157c4SSatish Balay           flg = PETSC_FALSE;
1054bea157c4SSatish Balay           break;
1055d9b7c43dSSatish Balay         }
1056bea157c4SSatish Balay       }
1057bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
1058bea157c4SSatish Balay         sizes[j] = bs;
1059bea157c4SSatish Balay         i+= bs;
1060bea157c4SSatish Balay       } else {
1061bea157c4SSatish Balay         sizes[j] = 1;
1062bea157c4SSatish Balay         i++;
1063bea157c4SSatish Balay       }
1064bea157c4SSatish Balay     }
1065bea157c4SSatish Balay   }
1066bea157c4SSatish Balay   *bs_max = j;
10673a40ed3dSBarry Smith   PetscFunctionReturn(0);
1068d9b7c43dSSatish Balay }
1069d9b7c43dSSatish Balay 
10704a2ae208SSatish Balay #undef __FUNCT__
10714a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ"
1072268466fbSBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,const PetscScalar *diag)
1073d9b7c43dSSatish Balay {
1074d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
1075b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
1076bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
107787828ca2SBarry Smith   PetscScalar zero = 0.0;
10783f1db9ecSBarry Smith   MatScalar   *aa;
1079d9b7c43dSSatish Balay 
10803a40ed3dSBarry Smith   PetscFunctionBegin;
1081d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
1082b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
1083d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
1084d9b7c43dSSatish Balay 
1085bea157c4SSatish Balay   /* allocate memory for rows,sizes */
1086b0a32e0cSBarry Smith   ierr  = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
1087bea157c4SSatish Balay   sizes = rows + is_n;
1088bea157c4SSatish Balay 
1089563b5814SBarry Smith   /* copy IS values to rows, and sort them */
1090bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
1091bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
1092dffd3267SBarry Smith   if (baij->keepzeroedrows) {
1093dffd3267SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
1094dffd3267SBarry Smith     bs_max = is_n;
1095dffd3267SBarry Smith   } else {
1096bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
1097dffd3267SBarry Smith   }
1098b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
1099bea157c4SSatish Balay 
1100bea157c4SSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
1101bea157c4SSatish Balay     row   = rows[j];
1102273d9f13SBarry Smith     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
1103bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1104bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
1105dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
1106bea157c4SSatish Balay       if (diag) {
1107bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
1108bea157c4SSatish Balay           baij->ilen[row/bs]       = 1;
1109bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
1110bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
1111a07cd24cSSatish Balay         }
1112563b5814SBarry Smith         /* Now insert all the diagonal values for this bs */
1113bea157c4SSatish Balay         for (k=0; k<bs; k++) {
1114bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
1115bea157c4SSatish Balay         }
1116bea157c4SSatish Balay       } else { /* (!diag) */
1117bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
1118bea157c4SSatish Balay       } /* end (!diag) */
1119bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
1120aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
112129bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
1122bea157c4SSatish Balay #endif
1123bea157c4SSatish Balay       for (k=0; k<count; k++) {
1124d9b7c43dSSatish Balay         aa[0] =  zero;
1125d9b7c43dSSatish Balay         aa    += bs;
1126d9b7c43dSSatish Balay       }
1127d9b7c43dSSatish Balay       if (diag) {
1128f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1129d9b7c43dSSatish Balay       }
1130d9b7c43dSSatish Balay     }
1131bea157c4SSatish Balay   }
1132bea157c4SSatish Balay 
1133606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
11349a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1136d9b7c43dSSatish Balay }
11371c351548SSatish Balay 
11384a2ae208SSatish Balay #undef __FUNCT__
11394a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ"
1140f15d580aSBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
11412d61bbb3SSatish Balay {
11422d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
11432d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1144273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
11452d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1146549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
1147273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
11483f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
11492d61bbb3SSatish Balay 
11502d61bbb3SSatish Balay   PetscFunctionBegin;
11512d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
11522d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
11535ef9f2a5SBarry Smith     if (row < 0) continue;
1154aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1155590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
11562d61bbb3SSatish Balay #endif
11572d61bbb3SSatish Balay     rp   = aj + ai[brow];
11582d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
11592d61bbb3SSatish Balay     rmax = imax[brow];
11602d61bbb3SSatish Balay     nrow = ailen[brow];
11612d61bbb3SSatish Balay     low  = 0;
11622d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
11635ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1164aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1165590ac198SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1);
11662d61bbb3SSatish Balay #endif
11672d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
11682d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
11692d61bbb3SSatish Balay       if (roworiented) {
11705ef9f2a5SBarry Smith         value = v[l + k*n];
11712d61bbb3SSatish Balay       } else {
11722d61bbb3SSatish Balay         value = v[k + l*m];
11732d61bbb3SSatish Balay       }
11742d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
11752d61bbb3SSatish Balay       while (high-low > 7) {
11762d61bbb3SSatish Balay         t = (low+high)/2;
11772d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
11782d61bbb3SSatish Balay         else              low  = t;
11792d61bbb3SSatish Balay       }
11802d61bbb3SSatish Balay       for (i=low; i<high; i++) {
11812d61bbb3SSatish Balay         if (rp[i] > bcol) break;
11822d61bbb3SSatish Balay         if (rp[i] == bcol) {
11832d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
11842d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
11852d61bbb3SSatish Balay           else                  *bap  = value;
11862d61bbb3SSatish Balay           goto noinsert1;
11872d61bbb3SSatish Balay         }
11882d61bbb3SSatish Balay       }
11892d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
119029bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
11912d61bbb3SSatish Balay       if (nrow >= rmax) {
11922d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
11932d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
11943f1db9ecSBarry Smith         MatScalar *new_a;
11952d61bbb3SSatish Balay 
119629bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
11972d61bbb3SSatish Balay 
11982d61bbb3SSatish Balay         /* Malloc new storage space */
11993f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1200b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
12012d61bbb3SSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
12022d61bbb3SSatish Balay         new_i   = new_j + new_nz;
12032d61bbb3SSatish Balay 
12042d61bbb3SSatish Balay         /* copy over old data into new slots */
12052d61bbb3SSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
12062d61bbb3SSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1207549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
12082d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1209549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1210549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1211549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1212549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
12132d61bbb3SSatish Balay         /* free up old matrix storage */
1214606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1215606d414cSSatish Balay         if (!a->singlemalloc) {
1216606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1217606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1218606d414cSSatish Balay         }
12192d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1220c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
12212d61bbb3SSatish Balay 
12222d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
12232d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1224b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
12252d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
12262d61bbb3SSatish Balay         a->reallocs++;
12272d61bbb3SSatish Balay         a->nz++;
12282d61bbb3SSatish Balay       }
12292d61bbb3SSatish Balay       N = nrow++ - 1;
12302d61bbb3SSatish Balay       /* shift up all the later entries in this row */
12312d61bbb3SSatish Balay       for (ii=N; ii>=i; ii--) {
12322d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1233549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
12342d61bbb3SSatish Balay       }
1235549d3d68SSatish Balay       if (N>=i) {
1236549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1237549d3d68SSatish Balay       }
12382d61bbb3SSatish Balay       rp[i]                      = bcol;
12392d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
12402d61bbb3SSatish Balay       noinsert1:;
12412d61bbb3SSatish Balay       low = i;
12422d61bbb3SSatish Balay     }
12432d61bbb3SSatish Balay     ailen[brow] = nrow;
12442d61bbb3SSatish Balay   }
12452d61bbb3SSatish Balay   PetscFunctionReturn(0);
12462d61bbb3SSatish Balay }
12472d61bbb3SSatish Balay 
12482d61bbb3SSatish Balay 
12494a2ae208SSatish Balay #undef __FUNCT__
12504a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ"
1251b380c88cSHong Zhang int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatFactorInfo *info)
12522d61bbb3SSatish Balay {
12532d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
12542d61bbb3SSatish Balay   Mat         outA;
12552d61bbb3SSatish Balay   int         ierr;
1256667159a5SBarry Smith   PetscTruth  row_identity,col_identity;
12572d61bbb3SSatish Balay 
12582d61bbb3SSatish Balay   PetscFunctionBegin;
1259d3d32019SBarry Smith   if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
1260667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1261667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1262667159a5SBarry Smith   if (!row_identity || !col_identity) {
126329bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1264667159a5SBarry Smith   }
12652d61bbb3SSatish Balay 
12662d61bbb3SSatish Balay   outA          = inA;
12672d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
12682d61bbb3SSatish Balay 
12692d61bbb3SSatish Balay   if (!a->diag) {
1270c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
12712d61bbb3SSatish Balay   }
1272cf242676SKris Buschelman 
1273c38d4ed2SBarry Smith   a->row        = row;
1274c38d4ed2SBarry Smith   a->col        = col;
1275c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1276c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1277c38d4ed2SBarry Smith 
1278c38d4ed2SBarry Smith   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
12794c49b128SBarry Smith   ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1280b0a32e0cSBarry Smith   PetscLogObjectParent(inA,a->icol);
1281c38d4ed2SBarry Smith 
1282cf242676SKris Buschelman   /*
1283cf242676SKris Buschelman       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1284cf242676SKris Buschelman       for ILU(0) factorization with natural ordering
1285cf242676SKris Buschelman   */
1286cf242676SKris Buschelman   if (a->bs < 8) {
1287cf242676SKris Buschelman     ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr);
1288cf242676SKris Buschelman   } else {
1289c38d4ed2SBarry Smith     if (!a->solve_work) {
129087828ca2SBarry Smith       ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
129187828ca2SBarry Smith       PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar));
1292c38d4ed2SBarry Smith     }
12932d61bbb3SSatish Balay   }
12942d61bbb3SSatish Balay 
1295667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1296667159a5SBarry Smith 
12972d61bbb3SSatish Balay   PetscFunctionReturn(0);
12982d61bbb3SSatish Balay }
12994a2ae208SSatish Balay #undef __FUNCT__
13004a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ"
1301ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1302ba4ca20aSSatish Balay {
1303c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1304ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1305d132466eSBarry Smith   int               ierr;
1306ba4ca20aSSatish Balay 
13073a40ed3dSBarry Smith   PetscFunctionBegin;
1308c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1309d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1310d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
13113a40ed3dSBarry Smith   PetscFunctionReturn(0);
1312ba4ca20aSSatish Balay }
1313d9b7c43dSSatish Balay 
1314fb2e594dSBarry Smith EXTERN_C_BEGIN
13154a2ae208SSatish Balay #undef __FUNCT__
13164a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
131727a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
131827a8da17SBarry Smith {
131927a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
132014db4f74SSatish Balay   int         i,nz,nbs;
132127a8da17SBarry Smith 
132227a8da17SBarry Smith   PetscFunctionBegin;
132314db4f74SSatish Balay   nz  = baij->maxnz/baij->bs2;
132414db4f74SSatish Balay   nbs = baij->nbs;
132527a8da17SBarry Smith   for (i=0; i<nz; i++) {
132627a8da17SBarry Smith     baij->j[i] = indices[i];
132727a8da17SBarry Smith   }
132827a8da17SBarry Smith   baij->nz = nz;
132914db4f74SSatish Balay   for (i=0; i<nbs; i++) {
133027a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
133127a8da17SBarry Smith   }
133227a8da17SBarry Smith 
133327a8da17SBarry Smith   PetscFunctionReturn(0);
133427a8da17SBarry Smith }
1335fb2e594dSBarry Smith EXTERN_C_END
133627a8da17SBarry Smith 
13374a2ae208SSatish Balay #undef __FUNCT__
13384a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices"
133927a8da17SBarry Smith /*@
134027a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
134127a8da17SBarry Smith        in the matrix.
134227a8da17SBarry Smith 
134327a8da17SBarry Smith   Input Parameters:
134427a8da17SBarry Smith +  mat - the SeqBAIJ matrix
134527a8da17SBarry Smith -  indices - the column indices
134627a8da17SBarry Smith 
134715091d37SBarry Smith   Level: advanced
134815091d37SBarry Smith 
134927a8da17SBarry Smith   Notes:
135027a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
135127a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
135227a8da17SBarry Smith   of the MatSetValues() operation.
135327a8da17SBarry Smith 
135427a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
135527a8da17SBarry Smith   MatCreateSeqBAIJ().
135627a8da17SBarry Smith 
135727a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
135827a8da17SBarry Smith 
135927a8da17SBarry Smith @*/
136027a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
136127a8da17SBarry Smith {
136227a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
136327a8da17SBarry Smith 
136427a8da17SBarry Smith   PetscFunctionBegin;
136527a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1366c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
136727a8da17SBarry Smith   if (f) {
136827a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
136927a8da17SBarry Smith   } else {
137029bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
137127a8da17SBarry Smith   }
137227a8da17SBarry Smith   PetscFunctionReturn(0);
137327a8da17SBarry Smith }
137427a8da17SBarry Smith 
13754a2ae208SSatish Balay #undef __FUNCT__
13764a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ"
1377273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v)
1378273d9f13SBarry Smith {
1379273d9f13SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1380273d9f13SBarry Smith   int          ierr,i,j,n,row,bs,*ai,*aj,mbs;
1381273d9f13SBarry Smith   PetscReal    atmp;
138287828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
1383273d9f13SBarry Smith   MatScalar    *aa;
1384273d9f13SBarry Smith   int          ncols,brow,krow,kcol;
1385273d9f13SBarry Smith 
1386273d9f13SBarry Smith   PetscFunctionBegin;
1387273d9f13SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1388273d9f13SBarry Smith   bs   = a->bs;
1389273d9f13SBarry Smith   aa   = a->a;
1390273d9f13SBarry Smith   ai   = a->i;
1391273d9f13SBarry Smith   aj   = a->j;
1392273d9f13SBarry Smith   mbs = a->mbs;
1393273d9f13SBarry Smith 
1394273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1395273d9f13SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1396273d9f13SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1397273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1398273d9f13SBarry Smith   for (i=0; i<mbs; i++) {
1399273d9f13SBarry Smith     ncols = ai[1] - ai[0]; ai++;
1400273d9f13SBarry Smith     brow  = bs*i;
1401273d9f13SBarry Smith     for (j=0; j<ncols; j++){
1402273d9f13SBarry Smith       /* bcol = bs*(*aj); */
1403273d9f13SBarry Smith       for (kcol=0; kcol<bs; kcol++){
1404273d9f13SBarry Smith         for (krow=0; krow<bs; krow++){
1405273d9f13SBarry Smith           atmp = PetscAbsScalar(*aa); aa++;
1406273d9f13SBarry Smith           row = brow + krow;    /* row index */
1407273d9f13SBarry Smith           /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */
1408273d9f13SBarry Smith           if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp;
1409273d9f13SBarry Smith         }
1410273d9f13SBarry Smith       }
1411273d9f13SBarry Smith       aj++;
1412273d9f13SBarry Smith     }
1413273d9f13SBarry Smith   }
1414273d9f13SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1415273d9f13SBarry Smith   PetscFunctionReturn(0);
1416273d9f13SBarry Smith }
1417273d9f13SBarry Smith 
14184a2ae208SSatish Balay #undef __FUNCT__
14194a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ"
1420273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A)
1421273d9f13SBarry Smith {
1422273d9f13SBarry Smith   int        ierr;
1423273d9f13SBarry Smith 
1424273d9f13SBarry Smith   PetscFunctionBegin;
1425273d9f13SBarry Smith   ierr =  MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1426273d9f13SBarry Smith   PetscFunctionReturn(0);
1427273d9f13SBarry Smith }
1428273d9f13SBarry Smith 
14294a2ae208SSatish Balay #undef __FUNCT__
14304a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ"
1431*4e7234bfSBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar *array[])
1432f2a5309cSSatish Balay {
1433f2a5309cSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1434f2a5309cSSatish Balay   PetscFunctionBegin;
1435f2a5309cSSatish Balay   *array = a->a;
1436f2a5309cSSatish Balay   PetscFunctionReturn(0);
1437f2a5309cSSatish Balay }
1438f2a5309cSSatish Balay 
14394a2ae208SSatish Balay #undef __FUNCT__
14404a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ"
1441*4e7234bfSBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar *array[])
1442f2a5309cSSatish Balay {
1443f2a5309cSSatish Balay   PetscFunctionBegin;
1444f2a5309cSSatish Balay   PetscFunctionReturn(0);
1445f2a5309cSSatish Balay }
1446f2a5309cSSatish Balay 
144742ee4b1aSHong Zhang #include "petscblaslapack.h"
144842ee4b1aSHong Zhang #undef __FUNCT__
144942ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqBAIJ"
1450268466fbSBarry Smith int MatAXPY_SeqBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
145142ee4b1aSHong Zhang {
145242ee4b1aSHong Zhang   Mat_SeqBAIJ  *x  = (Mat_SeqBAIJ *)X->data,*y = (Mat_SeqBAIJ *)Y->data;
14536550863bSHong Zhang   int          ierr,one=1,i,bs=y->bs,j,bs2;
145442ee4b1aSHong Zhang 
145542ee4b1aSHong Zhang   PetscFunctionBegin;
145642ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1457268466fbSBarry Smith     BLaxpy_(&x->nz,(PetscScalar*)a,x->a,&one,y->a,&one);
1458c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1459c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1460c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1461c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1462c537a176SHong Zhang     }
1463c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1464c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1465c4319e64SHong Zhang       y->XtoY = X;
1466c537a176SHong Zhang     }
1467c4319e64SHong Zhang     bs2 = bs*bs;
1468c537a176SHong Zhang     for (i=0; i<x->nz; i++) {
1469c4319e64SHong Zhang       j = 0;
1470c4319e64SHong Zhang       while (j < bs2){
14716550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1472c4319e64SHong Zhang         j++;
1473c537a176SHong Zhang       }
1474c4319e64SHong Zhang     }
1475c4319e64SHong 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));
147642ee4b1aSHong Zhang   } else {
147742ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
147842ee4b1aSHong Zhang   }
147942ee4b1aSHong Zhang   PetscFunctionReturn(0);
148042ee4b1aSHong Zhang }
148142ee4b1aSHong Zhang 
14822593348eSBarry Smith /* -------------------------------------------------------------------*/
1483cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1484cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1485cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1486cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1487cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
14887c922b88SBarry Smith        MatMultTranspose_SeqBAIJ,
14897c922b88SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1490cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1491cc2dc46cSBarry Smith        0,
1492cc2dc46cSBarry Smith        0,
1493cc2dc46cSBarry Smith        0,
1494cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1495cc2dc46cSBarry Smith        0,
1496b6490206SBarry Smith        0,
1497f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1498cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1499cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1500cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1501cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1502cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1503b6490206SBarry Smith        0,
1504cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1505cc2dc46cSBarry Smith        0,
1506cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1507cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1508cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1509cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1510cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1511cc2dc46cSBarry Smith        0,
1512cc2dc46cSBarry Smith        0,
1513273d9f13SBarry Smith        MatSetUpPreallocation_SeqBAIJ,
1514cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1515cc2dc46cSBarry Smith        0,
1516f2a5309cSSatish Balay        MatGetArray_SeqBAIJ,
1517f2a5309cSSatish Balay        MatRestoreArray_SeqBAIJ,
15182e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1519cc2dc46cSBarry Smith        0,
1520cc2dc46cSBarry Smith        0,
1521cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1522cc2dc46cSBarry Smith        0,
152342ee4b1aSHong Zhang        MatAXPY_SeqBAIJ,
1524cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1525cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1526cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1527cc2dc46cSBarry Smith        0,
1528cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1529cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1530cc2dc46cSBarry Smith        0,
1531cc2dc46cSBarry Smith        0,
1532cc2dc46cSBarry Smith        0,
1533cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
15343b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
153592c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1536cc2dc46cSBarry Smith        0,
1537cc2dc46cSBarry Smith        0,
1538cc2dc46cSBarry Smith        0,
1539cc2dc46cSBarry Smith        0,
1540cc2dc46cSBarry Smith        0,
1541cc2dc46cSBarry Smith        0,
1542d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1543cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1544b9b97703SBarry Smith        MatDestroy_SeqBAIJ,
1545b9b97703SBarry Smith        MatView_SeqBAIJ,
15463a64cc32SBarry Smith        MatGetPetscMaps_Petsc,
1547273d9f13SBarry Smith        0,
1548273d9f13SBarry Smith        0,
1549273d9f13SBarry Smith        0,
1550273d9f13SBarry Smith        0,
1551273d9f13SBarry Smith        0,
1552273d9f13SBarry Smith        0,
1553273d9f13SBarry Smith        MatGetRowMax_SeqBAIJ,
1554273d9f13SBarry Smith        MatConvert_Basic};
15552593348eSBarry Smith 
15563e90b805SBarry Smith EXTERN_C_BEGIN
15574a2ae208SSatish Balay #undef __FUNCT__
15584a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ"
15593e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
15603e90b805SBarry Smith {
15613e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1562273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1563d132466eSBarry Smith   int         ierr;
15643e90b805SBarry Smith 
15653e90b805SBarry Smith   PetscFunctionBegin;
15663e90b805SBarry Smith   if (aij->nonew != 1) {
156729bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15683e90b805SBarry Smith   }
15693e90b805SBarry Smith 
15703e90b805SBarry Smith   /* allocate space for values if not already there */
15713e90b805SBarry Smith   if (!aij->saved_values) {
157287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
15733e90b805SBarry Smith   }
15743e90b805SBarry Smith 
15753e90b805SBarry Smith   /* copy values over */
157687828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
15773e90b805SBarry Smith   PetscFunctionReturn(0);
15783e90b805SBarry Smith }
15793e90b805SBarry Smith EXTERN_C_END
15803e90b805SBarry Smith 
15813e90b805SBarry Smith EXTERN_C_BEGIN
15824a2ae208SSatish Balay #undef __FUNCT__
15834a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ"
158432e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
15853e90b805SBarry Smith {
15863e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1587273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
15883e90b805SBarry Smith 
15893e90b805SBarry Smith   PetscFunctionBegin;
15903e90b805SBarry Smith   if (aij->nonew != 1) {
159129bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15923e90b805SBarry Smith   }
15933e90b805SBarry Smith   if (!aij->saved_values) {
159429bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
15953e90b805SBarry Smith   }
15963e90b805SBarry Smith 
15973e90b805SBarry Smith   /* copy values over */
159887828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
15993e90b805SBarry Smith   PetscFunctionReturn(0);
16003e90b805SBarry Smith }
16013e90b805SBarry Smith EXTERN_C_END
16023e90b805SBarry Smith 
1603273d9f13SBarry Smith EXTERN_C_BEGIN
1604273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*);
1605273d9f13SBarry Smith EXTERN_C_END
1606273d9f13SBarry Smith 
1607273d9f13SBarry Smith EXTERN_C_BEGIN
16084a2ae208SSatish Balay #undef __FUNCT__
1609a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqBAIJSetPreallocation_SeqBAIJ"
1610a23d5eceSKris Buschelman int MatSeqBAIJSetPreallocation_SeqBAIJ(Mat B,int bs,int nz,int *nnz)
1611a23d5eceSKris Buschelman {
1612a23d5eceSKris Buschelman   Mat_SeqBAIJ *b;
1613a23d5eceSKris Buschelman   int         i,len,ierr,mbs,nbs,bs2,newbs = bs;
1614a23d5eceSKris Buschelman   PetscTruth  flg;
1615a23d5eceSKris Buschelman 
1616a23d5eceSKris Buschelman   PetscFunctionBegin;
1617a23d5eceSKris Buschelman 
1618a23d5eceSKris Buschelman   B->preallocated = PETSC_TRUE;
1619a23d5eceSKris Buschelman   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1620a23d5eceSKris Buschelman   if (nnz && newbs != bs) {
1621a23d5eceSKris Buschelman     SETERRQ(1,"Cannot change blocksize from command line if setting nnz");
1622a23d5eceSKris Buschelman   }
1623a23d5eceSKris Buschelman   bs   = newbs;
1624a23d5eceSKris Buschelman 
1625a23d5eceSKris Buschelman   mbs  = B->m/bs;
1626a23d5eceSKris Buschelman   nbs  = B->n/bs;
1627a23d5eceSKris Buschelman   bs2  = bs*bs;
1628a23d5eceSKris Buschelman 
1629a23d5eceSKris Buschelman   if (mbs*bs!=B->m || nbs*bs!=B->n) {
1630a23d5eceSKris Buschelman     SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs);
1631a23d5eceSKris Buschelman   }
1632a23d5eceSKris Buschelman 
1633a23d5eceSKris Buschelman   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1634a23d5eceSKris Buschelman   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
1635a23d5eceSKris Buschelman   if (nnz) {
1636a23d5eceSKris Buschelman     for (i=0; i<mbs; i++) {
1637a23d5eceSKris Buschelman       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
1638a23d5eceSKris Buschelman       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);
1639a23d5eceSKris Buschelman     }
1640a23d5eceSKris Buschelman   }
1641a23d5eceSKris Buschelman 
1642a23d5eceSKris Buschelman   b       = (Mat_SeqBAIJ*)B->data;
1643a23d5eceSKris Buschelman   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
1644a23d5eceSKris Buschelman   B->ops->solve               = MatSolve_SeqBAIJ_Update;
1645a23d5eceSKris Buschelman   B->ops->solvetranspose      = MatSolveTranspose_SeqBAIJ_Update;
1646a23d5eceSKris Buschelman   if (!flg) {
1647a23d5eceSKris Buschelman     switch (bs) {
1648a23d5eceSKris Buschelman     case 1:
1649a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1650a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_1;
1651a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
1652a23d5eceSKris Buschelman       break;
1653a23d5eceSKris Buschelman     case 2:
1654a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1655a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_2;
1656a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
1657a23d5eceSKris Buschelman       break;
1658a23d5eceSKris Buschelman     case 3:
1659a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1660a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_3;
1661a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
1662a23d5eceSKris Buschelman       break;
1663a23d5eceSKris Buschelman     case 4:
1664a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1665a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_4;
1666a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
1667a23d5eceSKris Buschelman       break;
1668a23d5eceSKris Buschelman     case 5:
1669a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1670a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_5;
1671a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
1672a23d5eceSKris Buschelman       break;
1673a23d5eceSKris Buschelman     case 6:
1674a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
1675a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_6;
1676a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
1677a23d5eceSKris Buschelman       break;
1678a23d5eceSKris Buschelman     case 7:
1679a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
1680a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_7;
1681a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
1682a23d5eceSKris Buschelman       break;
1683a23d5eceSKris Buschelman     default:
1684a23d5eceSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
1685a23d5eceSKris Buschelman       B->ops->mult            = MatMult_SeqBAIJ_N;
1686a23d5eceSKris Buschelman       B->ops->multadd         = MatMultAdd_SeqBAIJ_N;
1687a23d5eceSKris Buschelman       break;
1688a23d5eceSKris Buschelman     }
1689a23d5eceSKris Buschelman   }
1690a23d5eceSKris Buschelman   b->bs      = bs;
1691a23d5eceSKris Buschelman   b->mbs     = mbs;
1692a23d5eceSKris Buschelman   b->nbs     = nbs;
1693a23d5eceSKris Buschelman   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
1694a23d5eceSKris Buschelman   if (!nnz) {
1695a23d5eceSKris Buschelman     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1696a23d5eceSKris Buschelman     else if (nz <= 0)        nz = 1;
1697a23d5eceSKris Buschelman     for (i=0; i<mbs; i++) b->imax[i] = nz;
1698a23d5eceSKris Buschelman     nz = nz*mbs;
1699a23d5eceSKris Buschelman   } else {
1700a23d5eceSKris Buschelman     nz = 0;
1701a23d5eceSKris Buschelman     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1702a23d5eceSKris Buschelman   }
1703a23d5eceSKris Buschelman 
1704a23d5eceSKris Buschelman   /* allocate the matrix space */
1705a23d5eceSKris Buschelman   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
1706a23d5eceSKris Buschelman   ierr  = PetscMalloc(len,&b->a);CHKERRQ(ierr);
1707a23d5eceSKris Buschelman   ierr  = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1708a23d5eceSKris Buschelman   b->j  = (int*)(b->a + nz*bs2);
1709a23d5eceSKris Buschelman   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
1710a23d5eceSKris Buschelman   b->i  = b->j + nz;
1711a23d5eceSKris Buschelman   b->singlemalloc = PETSC_TRUE;
1712a23d5eceSKris Buschelman 
1713a23d5eceSKris Buschelman   b->i[0] = 0;
1714a23d5eceSKris Buschelman   for (i=1; i<mbs+1; i++) {
1715a23d5eceSKris Buschelman     b->i[i] = b->i[i-1] + b->imax[i-1];
1716a23d5eceSKris Buschelman   }
1717a23d5eceSKris Buschelman 
1718a23d5eceSKris Buschelman   /* b->ilen will count nonzeros in each block row so far. */
1719a23d5eceSKris Buschelman   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1720a23d5eceSKris Buschelman   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1721a23d5eceSKris Buschelman   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
1722a23d5eceSKris Buschelman 
1723a23d5eceSKris Buschelman   b->bs               = bs;
1724a23d5eceSKris Buschelman   b->bs2              = bs2;
1725a23d5eceSKris Buschelman   b->mbs              = mbs;
1726a23d5eceSKris Buschelman   b->nz               = 0;
1727a23d5eceSKris Buschelman   b->maxnz            = nz*bs2;
1728a23d5eceSKris Buschelman   B->info.nz_unneeded = (PetscReal)b->maxnz;
1729a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1730a23d5eceSKris Buschelman }
1731a23d5eceSKris Buschelman EXTERN_C_END
1732a23d5eceSKris Buschelman 
1733a23d5eceSKris Buschelman EXTERN_C_BEGIN
1734a23d5eceSKris Buschelman #undef __FUNCT__
17354a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ"
1736273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B)
17372593348eSBarry Smith {
1738273d9f13SBarry Smith   int         ierr,size;
1739b6490206SBarry Smith   Mat_SeqBAIJ *b;
17403b2fbd54SBarry Smith 
17413a40ed3dSBarry Smith   PetscFunctionBegin;
1742273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
174329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1744b6490206SBarry Smith 
1745273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1746273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1747b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr);
1748b0a32e0cSBarry Smith   B->data = (void*)b;
1749549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1750549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
17512593348eSBarry Smith   B->factor           = 0;
17522593348eSBarry Smith   B->lupivotthreshold = 1.0;
175390f02eecSBarry Smith   B->mapping          = 0;
17542593348eSBarry Smith   b->row              = 0;
17552593348eSBarry Smith   b->col              = 0;
1756e51c0b9cSSatish Balay   b->icol             = 0;
17572593348eSBarry Smith   b->reallocs         = 0;
17583e90b805SBarry Smith   b->saved_values     = 0;
1759cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
1760563b5814SBarry Smith   b->setvalueslen     = 0;
1761563b5814SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1762563b5814SBarry Smith #endif
17632593348eSBarry Smith 
17643a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
17653a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1766a5ae1ecdSBarry Smith 
1767c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1768c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
17692593348eSBarry Smith   b->nonew            = 0;
17702593348eSBarry Smith   b->diag             = 0;
17712593348eSBarry Smith   b->solve_work       = 0;
1772de6a44a3SBarry Smith   b->mult_work        = 0;
17732a1b7f2aSHong Zhang   B->spptr            = 0;
17740e6d2581SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1775883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
1776c4319e64SHong Zhang   b->xtoy              = 0;
1777c4319e64SHong Zhang   b->XtoY              = 0;
17784e220ebcSLois Curfman McInnes 
1779f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
17803e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1781bc4b532fSSatish Balay                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1782f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
17833e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1784bc4b532fSSatish Balay                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1785f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
178627a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1787bc4b532fSSatish Balay                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1788a6175056SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqaij_C",
1789273d9f13SBarry Smith                                      "MatConvert_SeqBAIJ_SeqAIJ",
1790273d9f13SBarry Smith                                       MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr);
1791a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetPreallocation_C",
1792a23d5eceSKris Buschelman                                      "MatSeqBAIJSetPreallocation_SeqBAIJ",
1793a23d5eceSKris Buschelman                                       MatSeqBAIJSetPreallocation_SeqBAIJ);CHKERRQ(ierr);
17943a40ed3dSBarry Smith   PetscFunctionReturn(0);
17952593348eSBarry Smith }
1796273d9f13SBarry Smith EXTERN_C_END
17972593348eSBarry Smith 
17984a2ae208SSatish Balay #undef __FUNCT__
17994a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ"
18002e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
18012593348eSBarry Smith {
18022593348eSBarry Smith   Mat         C;
1803b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
180427a8da17SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1805de6a44a3SBarry Smith 
18063a40ed3dSBarry Smith   PetscFunctionBegin;
180729bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
18082593348eSBarry Smith 
18092593348eSBarry Smith   *B = 0;
1810273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1811273d9f13SBarry Smith   ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
1812273d9f13SBarry Smith   c    = (Mat_SeqBAIJ*)C->data;
181344cd7ae7SLois Curfman McInnes 
1814b6490206SBarry Smith   c->bs         = a->bs;
18159df24120SSatish Balay   c->bs2        = a->bs2;
1816b6490206SBarry Smith   c->mbs        = a->mbs;
1817b6490206SBarry Smith   c->nbs        = a->nbs;
181835d8aa7fSBarry Smith   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
18192593348eSBarry Smith 
1820b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1821b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
1822b6490206SBarry Smith   for (i=0; i<mbs; i++) {
18232593348eSBarry Smith     c->imax[i] = a->imax[i];
18242593348eSBarry Smith     c->ilen[i] = a->ilen[i];
18252593348eSBarry Smith   }
18262593348eSBarry Smith 
18272593348eSBarry Smith   /* allocate the matrix space */
1828c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
18293f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1830b0a32e0cSBarry Smith   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
18317e67e3f9SSatish Balay   c->j = (int*)(c->a + nz*bs2);
1832de6a44a3SBarry Smith   c->i = c->j + nz;
1833549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1834b6490206SBarry Smith   if (mbs > 0) {
1835549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
18362e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1837549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
18382e8a6d31SBarry Smith     } else {
1839549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
18402593348eSBarry Smith     }
18412593348eSBarry Smith   }
18422593348eSBarry Smith 
1843b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
18442593348eSBarry Smith   c->sorted      = a->sorted;
18452593348eSBarry Smith   c->roworiented = a->roworiented;
18462593348eSBarry Smith   c->nonew       = a->nonew;
18472593348eSBarry Smith 
18482593348eSBarry Smith   if (a->diag) {
1849b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1850b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1851b6490206SBarry Smith     for (i=0; i<mbs; i++) {
18522593348eSBarry Smith       c->diag[i] = a->diag[i];
18532593348eSBarry Smith     }
185498305bb5SBarry Smith   } else c->diag        = 0;
18552593348eSBarry Smith   c->nz                 = a->nz;
18562593348eSBarry Smith   c->maxnz              = a->maxnz;
18572593348eSBarry Smith   c->solve_work         = 0;
18582a1b7f2aSHong Zhang   C->spptr              = 0;     /* Dangerous -I'm throwing away a->spptr */
18597fc0212eSBarry Smith   c->mult_work          = 0;
1860273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
1861273d9f13SBarry Smith   C->assembled          = PETSC_TRUE;
18622593348eSBarry Smith   *B = C;
1863b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
18643a40ed3dSBarry Smith   PetscFunctionReturn(0);
18652593348eSBarry Smith }
18662593348eSBarry Smith 
1867273d9f13SBarry Smith EXTERN_C_BEGIN
18684a2ae208SSatish Balay #undef __FUNCT__
18694a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ"
1870b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A)
18712593348eSBarry Smith {
1872b6490206SBarry Smith   Mat_SeqBAIJ  *a;
18732593348eSBarry Smith   Mat          B;
1874f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1875b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
187635aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1877a7c10996SSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
187887828ca2SBarry Smith   PetscScalar  *aa;
187919bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1880ecc4ab6dSBarry Smith #if defined(PETSC_HAVE_DSCPACK)
1881ecc4ab6dSBarry Smith   PetscTruth   flag;
1882ecc4ab6dSBarry Smith #endif
18832593348eSBarry Smith 
18843a40ed3dSBarry Smith   PetscFunctionBegin;
1885b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1886de6a44a3SBarry Smith   bs2  = bs*bs;
1887b6490206SBarry Smith 
1888d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
188929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1890b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
18910752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1892552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
18932593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
18942593348eSBarry Smith 
1895d64ed03dSBarry Smith   if (header[3] < 0) {
189629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ");
1897d64ed03dSBarry Smith   }
1898d64ed03dSBarry Smith 
189929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
190035aab85fSBarry Smith 
190135aab85fSBarry Smith   /*
190235aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
190335aab85fSBarry Smith     divisible by the blocksize
190435aab85fSBarry Smith   */
1905b6490206SBarry Smith   mbs        = M/bs;
190635aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
190735aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
190835aab85fSBarry Smith   else                  mbs++;
190935aab85fSBarry Smith   if (extra_rows) {
1910b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
191135aab85fSBarry Smith   }
1912b6490206SBarry Smith 
19132593348eSBarry Smith   /* read in row lengths */
1914b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
19150752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
191635aab85fSBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
19172593348eSBarry Smith 
1918b6490206SBarry Smith   /* read in column indices */
1919b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
19200752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
192135aab85fSBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1922b6490206SBarry Smith 
1923b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1924b0a32e0cSBarry Smith   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
1925549d3d68SSatish Balay   ierr     = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1926b0a32e0cSBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
1927549d3d68SSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
192835aab85fSBarry Smith   masked   = mask + mbs;
1929b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1930b6490206SBarry Smith   for (i=0; i<mbs; i++) {
193135aab85fSBarry Smith     nmask = 0;
1932b6490206SBarry Smith     for (j=0; j<bs; j++) {
1933b6490206SBarry Smith       kmax = rowlengths[rowcount];
1934b6490206SBarry Smith       for (k=0; k<kmax; k++) {
193535aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
193635aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1937b6490206SBarry Smith       }
1938b6490206SBarry Smith       rowcount++;
1939b6490206SBarry Smith     }
194035aab85fSBarry Smith     browlengths[i] += nmask;
194135aab85fSBarry Smith     /* zero out the mask elements we set */
194235aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1943b6490206SBarry Smith   }
1944b6490206SBarry Smith 
19452593348eSBarry Smith   /* create our matrix */
19463eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
19472593348eSBarry Smith   B = *A;
1948b6490206SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
19492593348eSBarry Smith 
19502593348eSBarry Smith   /* set matrix "i" values */
1951de6a44a3SBarry Smith   a->i[0] = 0;
1952b6490206SBarry Smith   for (i=1; i<= mbs; i++) {
1953b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1954b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
19552593348eSBarry Smith   }
1956b6490206SBarry Smith   a->nz         = 0;
1957b6490206SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
19582593348eSBarry Smith 
1959b6490206SBarry Smith   /* read in nonzero values */
196087828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
19610752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
196235aab85fSBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1963b6490206SBarry Smith 
1964b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1965b6490206SBarry Smith   nzcount = 0; jcount = 0;
1966b6490206SBarry Smith   for (i=0; i<mbs; i++) {
1967b6490206SBarry Smith     nzcountb = nzcount;
196835aab85fSBarry Smith     nmask    = 0;
1969b6490206SBarry Smith     for (j=0; j<bs; j++) {
1970b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1971b6490206SBarry Smith       for (k=0; k<kmax; k++) {
197235aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
197335aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1974b6490206SBarry Smith       }
1975b6490206SBarry Smith     }
1976de6a44a3SBarry Smith     /* sort the masked values */
1977433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1978de6a44a3SBarry Smith 
1979b6490206SBarry Smith     /* set "j" values into matrix */
1980b6490206SBarry Smith     maskcount = 1;
198135aab85fSBarry Smith     for (j=0; j<nmask; j++) {
198235aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1983de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1984b6490206SBarry Smith     }
1985b6490206SBarry Smith     /* set "a" values into matrix */
1986de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1987b6490206SBarry Smith     for (j=0; j<bs; j++) {
1988b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1989b6490206SBarry Smith       for (k=0; k<kmax; k++) {
1990de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1991de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1992de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1993de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1994375fe846SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1995b6490206SBarry Smith       }
1996b6490206SBarry Smith     }
199735aab85fSBarry Smith     /* zero out the mask elements we set */
199835aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1999b6490206SBarry Smith   }
200029bbc08cSBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
2001b6490206SBarry Smith 
2002606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
2003606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
2004606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
2005606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
2006606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
2007b6490206SBarry Smith 
2008b6490206SBarry Smith   B->assembled = PETSC_TRUE;
2009de6a44a3SBarry Smith 
2010ecc4ab6dSBarry Smith #if defined(PETSC_HAVE_DSCPACK)
2011ecc4ab6dSBarry Smith   ierr = PetscOptionsHasName(B->prefix,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
2012ecc4ab6dSBarry Smith   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(B);CHKERRQ(ierr); }
2013ecc4ab6dSBarry Smith #endif
2014ecc4ab6dSBarry Smith 
20159c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
20163a40ed3dSBarry Smith   PetscFunctionReturn(0);
20172593348eSBarry Smith }
2018273d9f13SBarry Smith EXTERN_C_END
20192593348eSBarry Smith 
20204a2ae208SSatish Balay #undef __FUNCT__
20214a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ"
2022273d9f13SBarry Smith /*@C
2023273d9f13SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
2024273d9f13SBarry Smith    compressed row) format.  For good matrix assembly performance the
2025273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
2026273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
2027273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
20282593348eSBarry Smith 
2029273d9f13SBarry Smith    Collective on MPI_Comm
2030273d9f13SBarry Smith 
2031273d9f13SBarry Smith    Input Parameters:
2032273d9f13SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
2033273d9f13SBarry Smith .  bs - size of block
2034273d9f13SBarry Smith .  m - number of rows
2035273d9f13SBarry Smith .  n - number of columns
203635d8aa7fSBarry Smith .  nz - number of nonzero blocks  per block row (same for all rows)
203735d8aa7fSBarry Smith -  nnz - array containing the number of nonzero blocks in the various block rows
2038273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
2039273d9f13SBarry Smith 
2040273d9f13SBarry Smith    Output Parameter:
2041273d9f13SBarry Smith .  A - the matrix
2042273d9f13SBarry Smith 
2043273d9f13SBarry Smith    Options Database Keys:
2044273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
2045273d9f13SBarry Smith                      block calculations (much slower)
2046273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
2047273d9f13SBarry Smith 
2048273d9f13SBarry Smith    Level: intermediate
2049273d9f13SBarry Smith 
2050273d9f13SBarry Smith    Notes:
205135d8aa7fSBarry Smith    A nonzero block is any block that as 1 or more nonzeros in it
205235d8aa7fSBarry Smith 
2053273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
2054273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
2055273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
2056273d9f13SBarry Smith 
2057273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
2058273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2059273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
2060273d9f13SBarry Smith    matrices.
2061273d9f13SBarry Smith 
2062273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
2063273d9f13SBarry Smith @*/
2064ca01db9bSBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A)
2065273d9f13SBarry Smith {
2066273d9f13SBarry Smith   int ierr;
2067273d9f13SBarry Smith 
2068273d9f13SBarry Smith   PetscFunctionBegin;
2069273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
2070273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr);
2071273d9f13SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
2072273d9f13SBarry Smith   PetscFunctionReturn(0);
2073273d9f13SBarry Smith }
2074273d9f13SBarry Smith 
20754a2ae208SSatish Balay #undef __FUNCT__
20764a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation"
2077273d9f13SBarry Smith /*@C
2078273d9f13SBarry Smith    MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros
2079273d9f13SBarry Smith    per row in the matrix. For good matrix assembly performance the
2080273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
2081273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
2082273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
2083273d9f13SBarry Smith 
2084273d9f13SBarry Smith    Collective on MPI_Comm
2085273d9f13SBarry Smith 
2086273d9f13SBarry Smith    Input Parameters:
2087273d9f13SBarry Smith +  A - the matrix
2088273d9f13SBarry Smith .  bs - size of block
2089273d9f13SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
2090273d9f13SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
2091273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
2092273d9f13SBarry Smith 
2093273d9f13SBarry Smith    Options Database Keys:
2094273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
2095273d9f13SBarry Smith                      block calculations (much slower)
2096273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
2097273d9f13SBarry Smith 
2098273d9f13SBarry Smith    Level: intermediate
2099273d9f13SBarry Smith 
2100273d9f13SBarry Smith    Notes:
2101273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
2102273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
2103273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
2104273d9f13SBarry Smith 
2105273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
2106273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
2107273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
2108273d9f13SBarry Smith    matrices.
2109273d9f13SBarry Smith 
2110273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
2111273d9f13SBarry Smith @*/
2112ca01db9bSBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[])
2113273d9f13SBarry Smith {
2114ca01db9bSBarry Smith   int ierr,(*f)(Mat,int,int,const int[]);
2115273d9f13SBarry Smith 
2116273d9f13SBarry Smith   PetscFunctionBegin;
2117a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
2118a23d5eceSKris Buschelman   if (f) {
2119a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
2120273d9f13SBarry Smith   }
2121273d9f13SBarry Smith   PetscFunctionReturn(0);
2122273d9f13SBarry Smith }
2123