xref: /petsc/src/mat/impls/baij/seq/baij.c (revision ecc4ab6dbdf05484c6bb3f15e29a97309b3af462)
1bba1ac68SSatish Balay /*$Id: baij.c,v 1.245 2001/08/07 03:02:55 balay Exp bsmith $*/
22593348eSBarry Smith 
32593348eSBarry Smith /*
4b6490206SBarry Smith     Defines the basic matrix operations for the BAIJ (compressed row)
52593348eSBarry Smith   matrix storage format.
62593348eSBarry Smith */
770f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
81a6a6d98SBarry Smith #include "src/vec/vecimpl.h"
91a6a6d98SBarry Smith #include "src/inline/spops.h"
10325e03aeSBarry Smith #include "petscsys.h"                     /*I "petscmat.h" I*/
113b547af2SSatish Balay 
12af674e45SBarry Smith /*
13af674e45SBarry Smith     Special version for Fun3d sequential benchmark
14af674e45SBarry Smith */
15af674e45SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
16af674e45SBarry Smith #define matsetvaluesblocked4_ MATSETVALUESBLOCKED4
17af674e45SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
18af674e45SBarry Smith #define matsetvaluesblocked4_ matsetvaluesblocked4
19af674e45SBarry Smith #endif
20af674e45SBarry Smith 
219c8c1041SBarry Smith EXTERN_C_BEGIN
22af674e45SBarry Smith #undef __FUNCT__
23af674e45SBarry Smith #define __FUNCT__ "matsetvaluesblocked4_"
244bb09213Spetsc void matsetvaluesblocked4_(Mat *AA,int *mm,int *im,int *nn,int *in,PetscScalar *v)
25af674e45SBarry Smith {
26af674e45SBarry Smith   Mat         A = *AA;
27af674e45SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
28a037b02bSBarry Smith   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,N,m = *mm,n = *nn;
2938fde467SHong Zhang   int         *ai=a->i,*ailen=a->ilen;
30a037b02bSBarry Smith   int         *aj=a->j,stepval;
314bb09213Spetsc   PetscScalar *value = v;
324bb09213Spetsc   MatScalar   *ap,*aa = a->a,*bap;
33af674e45SBarry Smith 
34af674e45SBarry Smith   PetscFunctionBegin;
35af674e45SBarry Smith   stepval = (n-1)*4;
36af674e45SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
37af674e45SBarry Smith     row  = im[k];
38af674e45SBarry Smith     rp   = aj + ai[row];
39af674e45SBarry Smith     ap   = aa + 16*ai[row];
40af674e45SBarry Smith     nrow = ailen[row];
41af674e45SBarry Smith     low  = 0;
42af674e45SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
43af674e45SBarry Smith       col = in[l];
44af674e45SBarry Smith       value = v + k*(stepval+4)*4 + l*4;
45af674e45SBarry Smith       low = 0; high = nrow;
46af674e45SBarry Smith       while (high-low > 7) {
47af674e45SBarry Smith         t = (low+high)/2;
48af674e45SBarry Smith         if (rp[t] > col) high = t;
49af674e45SBarry Smith         else             low  = t;
50af674e45SBarry Smith       }
51af674e45SBarry Smith       for (i=low; i<high; i++) {
52af674e45SBarry Smith         if (rp[i] > col) break;
53af674e45SBarry Smith         if (rp[i] == col) {
54af674e45SBarry Smith           bap  = ap +  16*i;
55af674e45SBarry Smith           for (ii=0; ii<4; ii++,value+=stepval) {
56af674e45SBarry Smith             for (jj=ii; jj<16; jj+=4) {
57af674e45SBarry Smith               bap[jj] += *value++;
58af674e45SBarry Smith             }
59af674e45SBarry Smith           }
60af674e45SBarry Smith           goto noinsert2;
61af674e45SBarry Smith         }
62af674e45SBarry Smith       }
63af674e45SBarry Smith       N = nrow++ - 1;
64af674e45SBarry Smith       /* shift up all the later entries in this row */
65af674e45SBarry Smith       for (ii=N; ii>=i; ii--) {
66af674e45SBarry Smith         rp[ii+1] = rp[ii];
67a037b02bSBarry Smith         PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar));
68af674e45SBarry Smith       }
69af674e45SBarry Smith       if (N >= i) {
70a037b02bSBarry Smith         PetscMemzero(ap+16*i,16*sizeof(MatScalar));
71af674e45SBarry Smith       }
72af674e45SBarry Smith       rp[i] = col;
73af674e45SBarry Smith       bap   = ap +  16*i;
74af674e45SBarry Smith       for (ii=0; ii<4; ii++,value+=stepval) {
75af674e45SBarry Smith         for (jj=ii; jj<16; jj+=4) {
76af674e45SBarry Smith           bap[jj] = *value++;
77af674e45SBarry Smith         }
78af674e45SBarry Smith       }
79af674e45SBarry Smith       noinsert2:;
80af674e45SBarry Smith       low = i;
81af674e45SBarry Smith     }
82af674e45SBarry Smith     ailen[row] = nrow;
83af674e45SBarry Smith   }
84af674e45SBarry Smith }
859c8c1041SBarry Smith EXTERN_C_END
86af674e45SBarry Smith 
87af674e45SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
88af674e45SBarry Smith #define matsetvalues4_ MATSETVALUES4
89af674e45SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
90af674e45SBarry Smith #define matsetvalues4_ matsetvalues4
91af674e45SBarry Smith #endif
92af674e45SBarry Smith 
939c8c1041SBarry Smith EXTERN_C_BEGIN
94af674e45SBarry Smith #undef __FUNCT__
95af674e45SBarry Smith #define __FUNCT__ "MatSetValues4_"
964bb09213Spetsc void matsetvalues4_(Mat *AA,int *mm,int *im,int *nn,int *in,PetscScalar *v)
97af674e45SBarry Smith {
98af674e45SBarry Smith   Mat         A = *AA;
99af674e45SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
100a037b02bSBarry Smith   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,N,n = *nn,m = *mm;
10138fde467SHong Zhang   int         *ai=a->i,*ailen=a->ilen;
102af674e45SBarry Smith   int         *aj=a->j,brow,bcol;
10338fde467SHong Zhang   int         ridx,cidx;
104af674e45SBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
105af674e45SBarry Smith 
106af674e45SBarry Smith   PetscFunctionBegin;
107af674e45SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
108af674e45SBarry Smith     row  = im[k]; brow = row/4;
109af674e45SBarry Smith     rp   = aj + ai[brow];
110af674e45SBarry Smith     ap   = aa + 16*ai[brow];
111af674e45SBarry Smith     nrow = ailen[brow];
112af674e45SBarry Smith     low  = 0;
113af674e45SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
114af674e45SBarry Smith       col = in[l]; bcol = col/4;
115af674e45SBarry Smith       ridx = row % 4; cidx = col % 4;
116af674e45SBarry Smith       value = v[l + k*n];
117af674e45SBarry Smith       low = 0; high = nrow;
118af674e45SBarry Smith       while (high-low > 7) {
119af674e45SBarry Smith         t = (low+high)/2;
120af674e45SBarry Smith         if (rp[t] > bcol) high = t;
121af674e45SBarry Smith         else              low  = t;
122af674e45SBarry Smith       }
123af674e45SBarry Smith       for (i=low; i<high; i++) {
124af674e45SBarry Smith         if (rp[i] > bcol) break;
125af674e45SBarry Smith         if (rp[i] == bcol) {
126af674e45SBarry Smith           bap  = ap +  16*i + 4*cidx + ridx;
127af674e45SBarry Smith           *bap += value;
128af674e45SBarry Smith           goto noinsert1;
129af674e45SBarry Smith         }
130af674e45SBarry Smith       }
131af674e45SBarry Smith       N = nrow++ - 1;
132af674e45SBarry Smith       /* shift up all the later entries in this row */
133af674e45SBarry Smith       for (ii=N; ii>=i; ii--) {
134af674e45SBarry Smith         rp[ii+1] = rp[ii];
135a037b02bSBarry Smith         PetscMemcpy(ap+16*(ii+1),ap+16*(ii),16*sizeof(MatScalar));
136af674e45SBarry Smith       }
137af674e45SBarry Smith       if (N>=i) {
138a037b02bSBarry Smith         PetscMemzero(ap+16*i,16*sizeof(MatScalar));
139af674e45SBarry Smith       }
140af674e45SBarry Smith       rp[i]                    = bcol;
141af674e45SBarry Smith       ap[16*i + 4*cidx + ridx] = value;
142af674e45SBarry Smith       noinsert1:;
143af674e45SBarry Smith       low = i;
144af674e45SBarry Smith     }
145af674e45SBarry Smith     ailen[brow] = nrow;
146af674e45SBarry Smith   }
147af674e45SBarry Smith }
1489c8c1041SBarry Smith EXTERN_C_END
149af674e45SBarry Smith 
15095d5f7c2SBarry Smith /*  UGLY, ugly, ugly
15187828ca2SBarry Smith    When MatScalar == PetscScalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does
1523477af42SKris Buschelman    not exist. Otherwise ..._MatScalar() takes matrix dlements in single precision and
15395d5f7c2SBarry Smith    inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ()
15495d5f7c2SBarry Smith    converts the entries into single precision and then calls ..._MatScalar() to put them
15595d5f7c2SBarry Smith    into the single precision data structures.
15695d5f7c2SBarry Smith */
15795d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
158ca44d042SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode);
15995d5f7c2SBarry Smith #else
16095d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ
16195d5f7c2SBarry Smith #endif
16204929863SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
16304929863SHong Zhang EXTERN int MatUseDSCPACK_MPIBAIJ(Mat);
16404929863SHong Zhang #endif
16595d5f7c2SBarry Smith 
1662d61bbb3SSatish Balay #define CHUNKSIZE  10
167de6a44a3SBarry Smith 
168be5855fcSBarry Smith /*
169be5855fcSBarry Smith      Checks for missing diagonals
170be5855fcSBarry Smith */
1714a2ae208SSatish Balay #undef __FUNCT__
1724a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ"
173c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A)
174be5855fcSBarry Smith {
175be5855fcSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
176883fce79SBarry Smith   int         *diag,*jj = a->j,i,ierr;
177be5855fcSBarry Smith 
178be5855fcSBarry Smith   PetscFunctionBegin;
179c4992f7dSBarry Smith   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
180883fce79SBarry Smith   diag = a->diag;
1810e8e8aceSBarry Smith   for (i=0; i<a->mbs; i++) {
182be5855fcSBarry Smith     if (jj[diag[i]] != i) {
18329bbc08cSBarry Smith       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
184be5855fcSBarry Smith     }
185be5855fcSBarry Smith   }
186be5855fcSBarry Smith   PetscFunctionReturn(0);
187be5855fcSBarry Smith }
188be5855fcSBarry Smith 
1894a2ae208SSatish Balay #undef __FUNCT__
1904a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ"
191c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A)
192de6a44a3SBarry Smith {
193de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
19482502324SSatish Balay   int         i,j,*diag,m = a->mbs,ierr;
195de6a44a3SBarry Smith 
1963a40ed3dSBarry Smith   PetscFunctionBegin;
197883fce79SBarry Smith   if (a->diag) PetscFunctionReturn(0);
198883fce79SBarry Smith 
199b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
200b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
2017fc0212eSBarry Smith   for (i=0; i<m; i++) {
202ecc615b2SBarry Smith     diag[i] = a->i[i+1];
203de6a44a3SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
204de6a44a3SBarry Smith       if (a->j[j] == i) {
205de6a44a3SBarry Smith         diag[i] = j;
206de6a44a3SBarry Smith         break;
207de6a44a3SBarry Smith       }
208de6a44a3SBarry Smith     }
209de6a44a3SBarry Smith   }
210de6a44a3SBarry Smith   a->diag = diag;
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
212de6a44a3SBarry Smith }
2132593348eSBarry Smith 
2142593348eSBarry Smith 
215ca44d042SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
2163b2fbd54SBarry Smith 
2174a2ae208SSatish Balay #undef __FUNCT__
2184a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqBAIJ"
2190e6d2581SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
2203b2fbd54SBarry Smith {
2213b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
2223b2fbd54SBarry Smith   int         ierr,n = a->mbs,i;
2233b2fbd54SBarry Smith 
2243a40ed3dSBarry Smith   PetscFunctionBegin;
2253b2fbd54SBarry Smith   *nn = n;
2263a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
2273b2fbd54SBarry Smith   if (symmetric) {
2283b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
2293b2fbd54SBarry Smith   } else if (oshift == 1) {
2303b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
231435da068SBarry Smith     int nz = a->i[n];
2323b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
2333b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]++;
2343b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
2353b2fbd54SBarry Smith   } else {
2363b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
2373b2fbd54SBarry Smith   }
2383b2fbd54SBarry Smith 
2393a40ed3dSBarry Smith   PetscFunctionReturn(0);
2403b2fbd54SBarry Smith }
2413b2fbd54SBarry Smith 
2424a2ae208SSatish Balay #undef __FUNCT__
2434a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ"
244435da068SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
2453b2fbd54SBarry Smith {
2463b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
247606d414cSSatish Balay   int         i,n = a->mbs,ierr;
2483b2fbd54SBarry Smith 
2493a40ed3dSBarry Smith   PetscFunctionBegin;
2503a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
2513b2fbd54SBarry Smith   if (symmetric) {
252606d414cSSatish Balay     ierr = PetscFree(*ia);CHKERRQ(ierr);
253606d414cSSatish Balay     ierr = PetscFree(*ja);CHKERRQ(ierr);
254af5da2bfSBarry Smith   } else if (oshift == 1) {
255435da068SBarry Smith     int nz = a->i[n]-1;
2563b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
2573b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]--;
2583b2fbd54SBarry Smith   }
2593a40ed3dSBarry Smith   PetscFunctionReturn(0);
2603b2fbd54SBarry Smith }
2613b2fbd54SBarry Smith 
2624a2ae208SSatish Balay #undef __FUNCT__
2634a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqBAIJ"
2642d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs)
2652d61bbb3SSatish Balay {
2662d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data;
2672d61bbb3SSatish Balay 
2682d61bbb3SSatish Balay   PetscFunctionBegin;
2692d61bbb3SSatish Balay   *bs = baij->bs;
2702d61bbb3SSatish Balay   PetscFunctionReturn(0);
2712d61bbb3SSatish Balay }
2722d61bbb3SSatish Balay 
2732d61bbb3SSatish Balay 
2744a2ae208SSatish Balay #undef __FUNCT__
2754a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqBAIJ"
276e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
2772d61bbb3SSatish Balay {
2782d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
279e51c0b9cSSatish Balay   int         ierr;
2802d61bbb3SSatish Balay 
281433994e6SBarry Smith   PetscFunctionBegin;
282aa482453SBarry Smith #if defined(PETSC_USE_LOG)
283b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
2842d61bbb3SSatish Balay #endif
285606d414cSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
286606d414cSSatish Balay   if (!a->singlemalloc) {
287606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
288606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
289606d414cSSatish Balay   }
290c38d4ed2SBarry Smith   if (a->row) {
291c38d4ed2SBarry Smith     ierr = ISDestroy(a->row);CHKERRQ(ierr);
292c38d4ed2SBarry Smith   }
293c38d4ed2SBarry Smith   if (a->col) {
294c38d4ed2SBarry Smith     ierr = ISDestroy(a->col);CHKERRQ(ierr);
295c38d4ed2SBarry Smith   }
296606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
297606d414cSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
298606d414cSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
299606d414cSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
300606d414cSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
301e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
302606d414cSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
303563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
304563b5814SBarry Smith   if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);}
305563b5814SBarry Smith #endif
306606d414cSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
3072d61bbb3SSatish Balay   PetscFunctionReturn(0);
3082d61bbb3SSatish Balay }
3092d61bbb3SSatish Balay 
3104a2ae208SSatish Balay #undef __FUNCT__
3114a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqBAIJ"
3122d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
3132d61bbb3SSatish Balay {
3142d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
3152d61bbb3SSatish Balay 
3162d61bbb3SSatish Balay   PetscFunctionBegin;
317aa275fccSKris Buschelman   switch (op) {
318aa275fccSKris Buschelman   case MAT_ROW_ORIENTED:
319aa275fccSKris Buschelman     a->roworiented    = PETSC_TRUE;
320aa275fccSKris Buschelman     break;
321aa275fccSKris Buschelman   case MAT_COLUMN_ORIENTED:
322aa275fccSKris Buschelman     a->roworiented    = PETSC_FALSE;
323aa275fccSKris Buschelman     break;
324aa275fccSKris Buschelman   case MAT_COLUMNS_SORTED:
325aa275fccSKris Buschelman     a->sorted         = PETSC_TRUE;
326aa275fccSKris Buschelman     break;
327aa275fccSKris Buschelman   case MAT_COLUMNS_UNSORTED:
328aa275fccSKris Buschelman     a->sorted         = PETSC_FALSE;
329aa275fccSKris Buschelman     break;
330aa275fccSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
331aa275fccSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
332aa275fccSKris Buschelman     break;
333aa275fccSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
334aa275fccSKris Buschelman     a->nonew          = 1;
335aa275fccSKris Buschelman     break;
336aa275fccSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
337aa275fccSKris Buschelman     a->nonew          = -1;
338aa275fccSKris Buschelman     break;
339aa275fccSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
340aa275fccSKris Buschelman     a->nonew          = -2;
341aa275fccSKris Buschelman     break;
342aa275fccSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
343aa275fccSKris Buschelman     a->nonew          = 0;
344aa275fccSKris Buschelman     break;
345aa275fccSKris Buschelman   case MAT_ROWS_SORTED:
346aa275fccSKris Buschelman   case MAT_ROWS_UNSORTED:
347aa275fccSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
348aa275fccSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
349aa275fccSKris Buschelman   case MAT_USE_HASH_TABLE:
350b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
351aa275fccSKris Buschelman     break;
352aa275fccSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
35329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
354aa275fccSKris Buschelman   default:
35529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
3562d61bbb3SSatish Balay   }
3572d61bbb3SSatish Balay   PetscFunctionReturn(0);
3582d61bbb3SSatish Balay }
3592d61bbb3SSatish Balay 
3604a2ae208SSatish Balay #undef __FUNCT__
3614a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ"
36287828ca2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
3632d61bbb3SSatish Balay {
3642d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
36582502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr;
3663f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
36787828ca2SBarry Smith   PetscScalar  *v_i;
3682d61bbb3SSatish Balay 
3692d61bbb3SSatish Balay   PetscFunctionBegin;
3702d61bbb3SSatish Balay   bs  = a->bs;
3712d61bbb3SSatish Balay   ai  = a->i;
3722d61bbb3SSatish Balay   aj  = a->j;
3732d61bbb3SSatish Balay   aa  = a->a;
3742d61bbb3SSatish Balay   bs2 = a->bs2;
3752d61bbb3SSatish Balay 
376273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
3772d61bbb3SSatish Balay 
3782d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
3792d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
3802d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
3812d61bbb3SSatish Balay   *nz = bs*M;
3822d61bbb3SSatish Balay 
3832d61bbb3SSatish Balay   if (v) {
3842d61bbb3SSatish Balay     *v = 0;
3852d61bbb3SSatish Balay     if (*nz) {
38687828ca2SBarry Smith       ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr);
3872d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
3882d61bbb3SSatish Balay         v_i  = *v + i*bs;
3892d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
3902d61bbb3SSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
3912d61bbb3SSatish Balay       }
3922d61bbb3SSatish Balay     }
3932d61bbb3SSatish Balay   }
3942d61bbb3SSatish Balay 
3952d61bbb3SSatish Balay   if (idx) {
3962d61bbb3SSatish Balay     *idx = 0;
3972d61bbb3SSatish Balay     if (*nz) {
398b0a32e0cSBarry Smith       ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr);
3992d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
4002d61bbb3SSatish Balay         idx_i = *idx + i*bs;
4012d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
4022d61bbb3SSatish Balay         for (j=0; j<bs; j++) {idx_i[j] = itmp++;}
4032d61bbb3SSatish Balay       }
4042d61bbb3SSatish Balay     }
4052d61bbb3SSatish Balay   }
4062d61bbb3SSatish Balay   PetscFunctionReturn(0);
4072d61bbb3SSatish Balay }
4082d61bbb3SSatish Balay 
4094a2ae208SSatish Balay #undef __FUNCT__
4104a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ"
41187828ca2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
4122d61bbb3SSatish Balay {
413606d414cSSatish Balay   int ierr;
414606d414cSSatish Balay 
4152d61bbb3SSatish Balay   PetscFunctionBegin;
416606d414cSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
417606d414cSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
4182d61bbb3SSatish Balay   PetscFunctionReturn(0);
4192d61bbb3SSatish Balay }
4202d61bbb3SSatish Balay 
4214a2ae208SSatish Balay #undef __FUNCT__
4224a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ"
4232d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
4242d61bbb3SSatish Balay {
4252d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
4262d61bbb3SSatish Balay   Mat         C;
4272d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
428273d9f13SBarry Smith   int         *rows,*cols,bs2=a->bs2;
42987828ca2SBarry Smith   PetscScalar *array;
4302d61bbb3SSatish Balay 
4312d61bbb3SSatish Balay   PetscFunctionBegin;
43229bbc08cSBarry Smith   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place");
433b0a32e0cSBarry Smith   ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr);
434549d3d68SSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
4352d61bbb3SSatish Balay 
436375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
43787828ca2SBarry Smith   ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr);
43887828ca2SBarry Smith   for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i];
439375fe846SBarry Smith #else
4403eda8832SBarry Smith   array = a->a;
441375fe846SBarry Smith #endif
442375fe846SBarry Smith 
4432d61bbb3SSatish Balay   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
444273d9f13SBarry Smith   ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
445606d414cSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
446b0a32e0cSBarry Smith   ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr);
4472d61bbb3SSatish Balay   cols = rows + bs;
4482d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4492d61bbb3SSatish Balay     cols[0] = i*bs;
4502d61bbb3SSatish Balay     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
4512d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
4522d61bbb3SSatish Balay     for (j=0; j<len; j++) {
4532d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
4542d61bbb3SSatish Balay       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
4552d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
4562d61bbb3SSatish Balay       array += bs2;
4572d61bbb3SSatish Balay     }
4582d61bbb3SSatish Balay   }
459606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
460375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
461375fe846SBarry Smith   ierr = PetscFree(array);
462375fe846SBarry Smith #endif
4632d61bbb3SSatish Balay 
4642d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4652d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4662d61bbb3SSatish Balay 
467c4992f7dSBarry Smith   if (B) {
4682d61bbb3SSatish Balay     *B = C;
4692d61bbb3SSatish Balay   } else {
470273d9f13SBarry Smith     ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
4712d61bbb3SSatish Balay   }
4722d61bbb3SSatish Balay   PetscFunctionReturn(0);
4732d61bbb3SSatish Balay }
4742d61bbb3SSatish Balay 
4754a2ae208SSatish Balay #undef __FUNCT__
4764a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary"
477b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer)
4782593348eSBarry Smith {
479b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
4809df24120SSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
48187828ca2SBarry Smith   PetscScalar *aa;
482ce6f0cecSBarry Smith   FILE        *file;
4832593348eSBarry Smith 
4843a40ed3dSBarry Smith   PetscFunctionBegin;
485b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
486b0a32e0cSBarry Smith   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
487552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
4883b2fbd54SBarry Smith 
489273d9f13SBarry Smith   col_lens[1] = A->m;
490273d9f13SBarry Smith   col_lens[2] = A->n;
4917e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
4922593348eSBarry Smith 
4932593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
494b6490206SBarry Smith   count = 0;
495b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
496b6490206SBarry Smith     for (j=0; j<bs; j++) {
497b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
498b6490206SBarry Smith     }
4992593348eSBarry Smith   }
500273d9f13SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
501606d414cSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
5022593348eSBarry Smith 
5032593348eSBarry Smith   /* store column indices (zero start index) */
504b0a32e0cSBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
505b6490206SBarry Smith   count = 0;
506b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
507b6490206SBarry Smith     for (j=0; j<bs; j++) {
508b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
509b6490206SBarry Smith         for (l=0; l<bs; l++) {
510b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
5112593348eSBarry Smith         }
5122593348eSBarry Smith       }
513b6490206SBarry Smith     }
514b6490206SBarry Smith   }
5150752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
516606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
5172593348eSBarry Smith 
5182593348eSBarry Smith   /* store nonzero values */
51987828ca2SBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
520b6490206SBarry Smith   count = 0;
521b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
522b6490206SBarry Smith     for (j=0; j<bs; j++) {
523b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
524b6490206SBarry Smith         for (l=0; l<bs; l++) {
5257e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
526b6490206SBarry Smith         }
527b6490206SBarry Smith       }
528b6490206SBarry Smith     }
529b6490206SBarry Smith   }
5300752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
531606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
532ce6f0cecSBarry Smith 
533b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
534ce6f0cecSBarry Smith   if (file) {
535ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
536ce6f0cecSBarry Smith   }
5373a40ed3dSBarry Smith   PetscFunctionReturn(0);
5382593348eSBarry Smith }
5392593348eSBarry Smith 
54004929863SHong Zhang extern int MatMPIBAIJFactorInfo_DSCPACK(Mat,PetscViewer);
54104929863SHong Zhang 
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII"
544b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer)
5452593348eSBarry Smith {
546b6490206SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
547fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
548f3ef73ceSBarry Smith   PetscViewerFormat format;
5492593348eSBarry Smith 
5503a40ed3dSBarry Smith   PetscFunctionBegin;
551b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
552456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
553b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
554fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
555bcd9e38bSBarry Smith     Mat aij;
556bcd9e38bSBarry Smith     ierr = MatConvert(A,MATSEQAIJ,&aij);CHKERRQ(ierr);
557bcd9e38bSBarry Smith     ierr = MatView(aij,viewer);CHKERRQ(ierr);
558bcd9e38bSBarry Smith     ierr = MatDestroy(aij);CHKERRQ(ierr);
55904929863SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
56004929863SHong Zhang #if defined(PETSC_HAVE_DSCPACK) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
56104929863SHong Zhang      ierr = MatMPIBAIJFactorInfo_DSCPACK(A,viewer);CHKERRQ(ierr);
56204929863SHong Zhang #endif
56304929863SHong Zhang      PetscFunctionReturn(0);
564fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
565b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
56644cd7ae7SLois Curfman McInnes     for (i=0; i<a->mbs; i++) {
56744cd7ae7SLois Curfman McInnes       for (j=0; j<bs; j++) {
568b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
56944cd7ae7SLois Curfman McInnes         for (k=a->i[i]; k<a->i[i+1]; k++) {
57044cd7ae7SLois Curfman McInnes           for (l=0; l<bs; l++) {
571aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5720e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
57362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %gi) ",bs*a->j[k]+l,
5740e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
5750e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
57662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %gi) ",bs*a->j[k]+l,
5770e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
5780e6d2581SBarry Smith             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
57962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
5800ef38995SBarry Smith             }
58144cd7ae7SLois Curfman McInnes #else
5820ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
58362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
5840ef38995SBarry Smith             }
58544cd7ae7SLois Curfman McInnes #endif
58644cd7ae7SLois Curfman McInnes           }
58744cd7ae7SLois Curfman McInnes         }
588b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
58944cd7ae7SLois Curfman McInnes       }
59044cd7ae7SLois Curfman McInnes     }
591b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
5920ef38995SBarry Smith   } else {
593b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
594b6490206SBarry Smith     for (i=0; i<a->mbs; i++) {
595b6490206SBarry Smith       for (j=0; j<bs; j++) {
596b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
597b6490206SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
598b6490206SBarry Smith           for (l=0; l<bs; l++) {
599aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6000e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
60162b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
6020e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
6030e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
60462b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
6050e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
6060ef38995SBarry Smith             } else {
60762b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
60888685aaeSLois Curfman McInnes             }
60988685aaeSLois Curfman McInnes #else
61062b941d6SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
61188685aaeSLois Curfman McInnes #endif
6122593348eSBarry Smith           }
6132593348eSBarry Smith         }
614b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6152593348eSBarry Smith       }
6162593348eSBarry Smith     }
617b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
618b6490206SBarry Smith   }
619b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
6212593348eSBarry Smith }
6222593348eSBarry Smith 
6234a2ae208SSatish Balay #undef __FUNCT__
6244a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom"
625b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
6263270192aSSatish Balay {
62777ed5343SBarry Smith   Mat          A = (Mat) Aa;
6283270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
629b65db4caSBarry Smith   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2;
6300e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
6313f1db9ecSBarry Smith   MatScalar    *aa;
632b0a32e0cSBarry Smith   PetscViewer  viewer;
6333270192aSSatish Balay 
6343a40ed3dSBarry Smith   PetscFunctionBegin;
6353270192aSSatish Balay 
636b65db4caSBarry Smith   /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/
63777ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
63877ed5343SBarry Smith 
639b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
64077ed5343SBarry Smith 
6413270192aSSatish Balay   /* loop over matrix elements drawing boxes */
642b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
6433270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
6443270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
645273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
6463270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
6473270192aSSatish Balay       aa = a->a + j*bs2;
6483270192aSSatish Balay       for (k=0; k<bs; k++) {
6493270192aSSatish Balay         for (l=0; l<bs; l++) {
6500e6d2581SBarry Smith           if (PetscRealPart(*aa++) >=  0.) continue;
651b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
6523270192aSSatish Balay         }
6533270192aSSatish Balay       }
6543270192aSSatish Balay     }
6553270192aSSatish Balay   }
656b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
6573270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
6583270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
659273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
6603270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
6613270192aSSatish Balay       aa = a->a + j*bs2;
6623270192aSSatish Balay       for (k=0; k<bs; k++) {
6633270192aSSatish Balay         for (l=0; l<bs; l++) {
6640e6d2581SBarry Smith           if (PetscRealPart(*aa++) != 0.) continue;
665b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
6663270192aSSatish Balay         }
6673270192aSSatish Balay       }
6683270192aSSatish Balay     }
6693270192aSSatish Balay   }
6703270192aSSatish Balay 
671b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
6723270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
6733270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
674273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
6753270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
6763270192aSSatish Balay       aa = a->a + j*bs2;
6773270192aSSatish Balay       for (k=0; k<bs; k++) {
6783270192aSSatish Balay         for (l=0; l<bs; l++) {
6790e6d2581SBarry Smith           if (PetscRealPart(*aa++) <= 0.) continue;
680b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
6813270192aSSatish Balay         }
6823270192aSSatish Balay       }
6833270192aSSatish Balay     }
6843270192aSSatish Balay   }
68577ed5343SBarry Smith   PetscFunctionReturn(0);
68677ed5343SBarry Smith }
6873270192aSSatish Balay 
6884a2ae208SSatish Balay #undef __FUNCT__
6894a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw"
690b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer)
69177ed5343SBarry Smith {
6927dce120fSSatish Balay   int          ierr;
6930e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,w,h;
694b0a32e0cSBarry Smith   PetscDraw    draw;
69577ed5343SBarry Smith   PetscTruth   isnull;
6963270192aSSatish Balay 
69777ed5343SBarry Smith   PetscFunctionBegin;
69877ed5343SBarry Smith 
699b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
700b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
70177ed5343SBarry Smith 
70277ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
703273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
70477ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
705b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
706b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
70777ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
7083a40ed3dSBarry Smith   PetscFunctionReturn(0);
7093270192aSSatish Balay }
7103270192aSSatish Balay 
7114a2ae208SSatish Balay #undef __FUNCT__
7124a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ"
713b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer)
7142593348eSBarry Smith {
71519bcc07fSBarry Smith   int        ierr;
7166831982aSBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
7172593348eSBarry Smith 
7183a40ed3dSBarry Smith   PetscFunctionBegin;
719b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
720b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
721fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
722fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
7230f5bd95cSBarry Smith   if (issocket) {
72429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
7250f5bd95cSBarry Smith   } else if (isascii){
7263a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
7270f5bd95cSBarry Smith   } else if (isbinary) {
7283a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
7290f5bd95cSBarry Smith   } else if (isdraw) {
7303a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
7315cd90555SBarry Smith   } else {
73229bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
7332593348eSBarry Smith   }
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
7352593348eSBarry Smith }
736b6490206SBarry Smith 
737cd0e1443SSatish Balay 
7384a2ae208SSatish Balay #undef __FUNCT__
7394a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ"
74087828ca2SBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
741cd0e1443SSatish Balay {
742cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
7432d61bbb3SSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
7442d61bbb3SSatish Balay   int        *ai = a->i,*ailen = a->ilen;
7452d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
7463f1db9ecSBarry Smith   MatScalar  *ap,*aa = a->a,zero = 0.0;
747cd0e1443SSatish Balay 
7483a40ed3dSBarry Smith   PetscFunctionBegin;
7492d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
750cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
75129bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
752273d9f13SBarry Smith     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
7532d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
7542c3acbe9SBarry Smith     nrow = ailen[brow];
7552d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
75629bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
757273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
7582d61bbb3SSatish Balay       col  = in[l] ;
7592d61bbb3SSatish Balay       bcol = col/bs;
7602d61bbb3SSatish Balay       cidx = col%bs;
7612d61bbb3SSatish Balay       ridx = row%bs;
7622d61bbb3SSatish Balay       high = nrow;
7632d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
7642d61bbb3SSatish Balay       while (high-low > 5) {
765cd0e1443SSatish Balay         t = (low+high)/2;
766cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
767cd0e1443SSatish Balay         else             low  = t;
768cd0e1443SSatish Balay       }
769cd0e1443SSatish Balay       for (i=low; i<high; i++) {
770cd0e1443SSatish Balay         if (rp[i] > bcol) break;
771cd0e1443SSatish Balay         if (rp[i] == bcol) {
7722d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
7732d61bbb3SSatish Balay           goto finished;
774cd0e1443SSatish Balay         }
775cd0e1443SSatish Balay       }
7762d61bbb3SSatish Balay       *v++ = zero;
7772d61bbb3SSatish Balay       finished:;
778cd0e1443SSatish Balay     }
779cd0e1443SSatish Balay   }
7803a40ed3dSBarry Smith   PetscFunctionReturn(0);
781cd0e1443SSatish Balay }
782cd0e1443SSatish Balay 
78395d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
7844a2ae208SSatish Balay #undef __FUNCT__
7854a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
78687828ca2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,PetscScalar *v,InsertMode addv)
78795d5f7c2SBarry Smith {
788563b5814SBarry Smith   Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data;
789563b5814SBarry Smith   int         ierr,i,N = m*n*b->bs2;
790563b5814SBarry Smith   MatScalar   *vsingle;
79195d5f7c2SBarry Smith 
79295d5f7c2SBarry Smith   PetscFunctionBegin;
793563b5814SBarry Smith   if (N > b->setvalueslen) {
794563b5814SBarry Smith     if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);}
795b0a32e0cSBarry Smith     ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr);
796563b5814SBarry Smith     b->setvalueslen  = N;
797563b5814SBarry Smith   }
798563b5814SBarry Smith   vsingle = b->setvaluescopy;
79995d5f7c2SBarry Smith   for (i=0; i<N; i++) {
80095d5f7c2SBarry Smith     vsingle[i] = v[i];
80195d5f7c2SBarry Smith   }
80295d5f7c2SBarry Smith   ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr);
80395d5f7c2SBarry Smith   PetscFunctionReturn(0);
80495d5f7c2SBarry Smith }
80595d5f7c2SBarry Smith #endif
80695d5f7c2SBarry Smith 
8072d61bbb3SSatish Balay 
8084a2ae208SSatish Balay #undef __FUNCT__
8094a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
81095d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is)
81192c4ed94SBarry Smith {
81292c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
8138a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
814273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
815549d3d68SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
816273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
81795d5f7c2SBarry Smith   MatScalar   *value = v,*ap,*aa = a->a,*bap;
81892c4ed94SBarry Smith 
8193a40ed3dSBarry Smith   PetscFunctionBegin;
8200e324ae4SSatish Balay   if (roworiented) {
8210e324ae4SSatish Balay     stepval = (n-1)*bs;
8220e324ae4SSatish Balay   } else {
8230e324ae4SSatish Balay     stepval = (m-1)*bs;
8240e324ae4SSatish Balay   }
82592c4ed94SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
82692c4ed94SBarry Smith     row  = im[k];
8275ef9f2a5SBarry Smith     if (row < 0) continue;
828aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
82929bbc08cSBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
83092c4ed94SBarry Smith #endif
83192c4ed94SBarry Smith     rp   = aj + ai[row];
83292c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
83392c4ed94SBarry Smith     rmax = imax[row];
83492c4ed94SBarry Smith     nrow = ailen[row];
83592c4ed94SBarry Smith     low  = 0;
83692c4ed94SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
8375ef9f2a5SBarry Smith       if (in[l] < 0) continue;
838aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
83929bbc08cSBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
84092c4ed94SBarry Smith #endif
84192c4ed94SBarry Smith       col = in[l];
84292c4ed94SBarry Smith       if (roworiented) {
8430e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
8440e324ae4SSatish Balay       } else {
8450e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
84692c4ed94SBarry Smith       }
84792c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
84892c4ed94SBarry Smith       while (high-low > 7) {
84992c4ed94SBarry Smith         t = (low+high)/2;
85092c4ed94SBarry Smith         if (rp[t] > col) high = t;
85192c4ed94SBarry Smith         else             low  = t;
85292c4ed94SBarry Smith       }
85392c4ed94SBarry Smith       for (i=low; i<high; i++) {
85492c4ed94SBarry Smith         if (rp[i] > col) break;
85592c4ed94SBarry Smith         if (rp[i] == col) {
8568a84c255SSatish Balay           bap  = ap +  bs2*i;
8570e324ae4SSatish Balay           if (roworiented) {
8588a84c255SSatish Balay             if (is == ADD_VALUES) {
859dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
860dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
8618a84c255SSatish Balay                   bap[jj] += *value++;
862dd9472c6SBarry Smith                 }
863dd9472c6SBarry Smith               }
8640e324ae4SSatish Balay             } else {
865dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
866dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
8670e324ae4SSatish Balay                   bap[jj] = *value++;
8688a84c255SSatish Balay                 }
869dd9472c6SBarry Smith               }
870dd9472c6SBarry Smith             }
8710e324ae4SSatish Balay           } else {
8720e324ae4SSatish Balay             if (is == ADD_VALUES) {
873dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
874dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
8750e324ae4SSatish Balay                   *bap++ += *value++;
876dd9472c6SBarry Smith                 }
877dd9472c6SBarry Smith               }
8780e324ae4SSatish Balay             } else {
879dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
880dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
8810e324ae4SSatish Balay                   *bap++  = *value++;
8820e324ae4SSatish Balay                 }
8838a84c255SSatish Balay               }
884dd9472c6SBarry Smith             }
885dd9472c6SBarry Smith           }
886f1241b54SBarry Smith           goto noinsert2;
88792c4ed94SBarry Smith         }
88892c4ed94SBarry Smith       }
88989280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
89029bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
89192c4ed94SBarry Smith       if (nrow >= rmax) {
89292c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
89392c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
8943f1db9ecSBarry Smith         MatScalar *new_a;
89592c4ed94SBarry Smith 
89629bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
89789280ab3SLois Curfman McInnes 
89892c4ed94SBarry Smith         /* malloc new storage space */
8993f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
900b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
90192c4ed94SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
90292c4ed94SBarry Smith         new_i   = new_j + new_nz;
90392c4ed94SBarry Smith 
90492c4ed94SBarry Smith         /* copy over old data into new slots */
90592c4ed94SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
90692c4ed94SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
907549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
90892c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
909549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
910549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
911549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
912549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
91392c4ed94SBarry Smith         /* free up old matrix storage */
914606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
915606d414cSSatish Balay         if (!a->singlemalloc) {
916606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
917606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
918606d414cSSatish Balay         }
91992c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
920c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
92192c4ed94SBarry Smith 
92292c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
92392c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
924b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
92592c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
92692c4ed94SBarry Smith         a->reallocs++;
92792c4ed94SBarry Smith         a->nz++;
92892c4ed94SBarry Smith       }
92992c4ed94SBarry Smith       N = nrow++ - 1;
93092c4ed94SBarry Smith       /* shift up all the later entries in this row */
93192c4ed94SBarry Smith       for (ii=N; ii>=i; ii--) {
93292c4ed94SBarry Smith         rp[ii+1] = rp[ii];
933549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
93492c4ed94SBarry Smith       }
935549d3d68SSatish Balay       if (N >= i) {
936549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
937549d3d68SSatish Balay       }
93892c4ed94SBarry Smith       rp[i] = col;
9398a84c255SSatish Balay       bap   = ap +  bs2*i;
9400e324ae4SSatish Balay       if (roworiented) {
941dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
942dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
9430e324ae4SSatish Balay             bap[jj] = *value++;
944dd9472c6SBarry Smith           }
945dd9472c6SBarry Smith         }
9460e324ae4SSatish Balay       } else {
947dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
948dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++) {
9490e324ae4SSatish Balay             *bap++  = *value++;
9500e324ae4SSatish Balay           }
951dd9472c6SBarry Smith         }
952dd9472c6SBarry Smith       }
953f1241b54SBarry Smith       noinsert2:;
95492c4ed94SBarry Smith       low = i;
95592c4ed94SBarry Smith     }
95692c4ed94SBarry Smith     ailen[row] = nrow;
95792c4ed94SBarry Smith   }
9583a40ed3dSBarry Smith   PetscFunctionReturn(0);
95992c4ed94SBarry Smith }
96092c4ed94SBarry Smith 
9614a2ae208SSatish Balay #undef __FUNCT__
9624a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ"
9638f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
964584200bdSSatish Balay {
965584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
966584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
967273d9f13SBarry Smith   int        m = A->m,*ip,N,*ailen = a->ilen;
968549d3d68SSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
9693f1db9ecSBarry Smith   MatScalar  *aa = a->a,*ap;
970a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
971a56a16c8SHong Zhang   PetscTruth   flag;
972a56a16c8SHong Zhang #endif
973584200bdSSatish Balay 
9743a40ed3dSBarry Smith   PetscFunctionBegin;
9753a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
976584200bdSSatish Balay 
97743ee02c3SBarry Smith   if (m) rmax = ailen[0];
978584200bdSSatish Balay   for (i=1; i<mbs; i++) {
979584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
980584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
981d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
982584200bdSSatish Balay     if (fshift) {
983a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
984584200bdSSatish Balay       N = ailen[i];
985584200bdSSatish Balay       for (j=0; j<N; j++) {
986584200bdSSatish Balay         ip[j-fshift] = ip[j];
987549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
988584200bdSSatish Balay       }
989584200bdSSatish Balay     }
990584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
991584200bdSSatish Balay   }
992584200bdSSatish Balay   if (mbs) {
993584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
994584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
995584200bdSSatish Balay   }
996584200bdSSatish Balay   /* reset ilen and imax for each row */
997584200bdSSatish Balay   for (i=0; i<mbs; i++) {
998584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
999584200bdSSatish Balay   }
1000a7c10996SSatish Balay   a->nz = ai[mbs];
1001584200bdSSatish Balay 
1002584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
1003584200bdSSatish Balay   if (fshift && a->diag) {
1004606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
1005b424e231SHong Zhang     PetscLogObjectMemory(A,-(mbs+1)*sizeof(int));
1006584200bdSSatish Balay     a->diag = 0;
1007584200bdSSatish Balay   }
1008bba1ac68SSatish 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);
1009bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs);
1010b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
1011e2f3b5e9SSatish Balay   a->reallocs          = 0;
10120e6d2581SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
1013cf4441caSHong Zhang 
1014a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
10152c535e4dSHong Zhang   ierr = PetscOptionsHasName(A->prefix,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
1016a56a16c8SHong Zhang   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(A);CHKERRQ(ierr); }
1017a56a16c8SHong Zhang #endif
1018a56a16c8SHong Zhang 
10193a40ed3dSBarry Smith   PetscFunctionReturn(0);
1020584200bdSSatish Balay }
1021584200bdSSatish Balay 
10225a838052SSatish Balay 
1023bea157c4SSatish Balay 
1024bea157c4SSatish Balay /*
1025bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
1026bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
1027bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
1028bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
1029bea157c4SSatish Balay */
10304a2ae208SSatish Balay #undef __FUNCT__
10314a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks"
1032bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
1033d9b7c43dSSatish Balay {
1034bea157c4SSatish Balay   int        i,j,k,row;
1035bea157c4SSatish Balay   PetscTruth flg;
10363a40ed3dSBarry Smith 
1037433994e6SBarry Smith   PetscFunctionBegin;
1038bea157c4SSatish Balay   for (i=0,j=0; i<n; j++) {
1039bea157c4SSatish Balay     row = idx[i];
1040bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
1041bea157c4SSatish Balay       sizes[j] = 1;
1042bea157c4SSatish Balay       i++;
1043e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
1044bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
1045bea157c4SSatish Balay       i++;
1046bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
1047bea157c4SSatish Balay       flg = PETSC_TRUE;
1048bea157c4SSatish Balay       for (k=1; k<bs; k++) {
1049bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
1050bea157c4SSatish Balay           flg = PETSC_FALSE;
1051bea157c4SSatish Balay           break;
1052d9b7c43dSSatish Balay         }
1053bea157c4SSatish Balay       }
1054bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
1055bea157c4SSatish Balay         sizes[j] = bs;
1056bea157c4SSatish Balay         i+= bs;
1057bea157c4SSatish Balay       } else {
1058bea157c4SSatish Balay         sizes[j] = 1;
1059bea157c4SSatish Balay         i++;
1060bea157c4SSatish Balay       }
1061bea157c4SSatish Balay     }
1062bea157c4SSatish Balay   }
1063bea157c4SSatish Balay   *bs_max = j;
10643a40ed3dSBarry Smith   PetscFunctionReturn(0);
1065d9b7c43dSSatish Balay }
1066d9b7c43dSSatish Balay 
10674a2ae208SSatish Balay #undef __FUNCT__
10684a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ"
106987828ca2SBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,PetscScalar *diag)
1070d9b7c43dSSatish Balay {
1071d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
1072b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
1073bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
107487828ca2SBarry Smith   PetscScalar zero = 0.0;
10753f1db9ecSBarry Smith   MatScalar   *aa;
1076d9b7c43dSSatish Balay 
10773a40ed3dSBarry Smith   PetscFunctionBegin;
1078d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
1079b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
1080d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
1081d9b7c43dSSatish Balay 
1082bea157c4SSatish Balay   /* allocate memory for rows,sizes */
1083b0a32e0cSBarry Smith   ierr  = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
1084bea157c4SSatish Balay   sizes = rows + is_n;
1085bea157c4SSatish Balay 
1086563b5814SBarry Smith   /* copy IS values to rows, and sort them */
1087bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
1088bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
1089dffd3267SBarry Smith   if (baij->keepzeroedrows) {
1090dffd3267SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
1091dffd3267SBarry Smith     bs_max = is_n;
1092dffd3267SBarry Smith   } else {
1093bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
1094dffd3267SBarry Smith   }
1095b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
1096bea157c4SSatish Balay 
1097bea157c4SSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
1098bea157c4SSatish Balay     row   = rows[j];
1099273d9f13SBarry Smith     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
1100bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1101bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
1102dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
1103bea157c4SSatish Balay       if (diag) {
1104bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
1105bea157c4SSatish Balay           baij->ilen[row/bs]       = 1;
1106bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
1107bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
1108a07cd24cSSatish Balay         }
1109563b5814SBarry Smith         /* Now insert all the diagonal values for this bs */
1110bea157c4SSatish Balay         for (k=0; k<bs; k++) {
1111bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
1112bea157c4SSatish Balay         }
1113bea157c4SSatish Balay       } else { /* (!diag) */
1114bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
1115bea157c4SSatish Balay       } /* end (!diag) */
1116bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
1117aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
111829bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
1119bea157c4SSatish Balay #endif
1120bea157c4SSatish Balay       for (k=0; k<count; k++) {
1121d9b7c43dSSatish Balay         aa[0] =  zero;
1122d9b7c43dSSatish Balay         aa    += bs;
1123d9b7c43dSSatish Balay       }
1124d9b7c43dSSatish Balay       if (diag) {
1125f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1126d9b7c43dSSatish Balay       }
1127d9b7c43dSSatish Balay     }
1128bea157c4SSatish Balay   }
1129bea157c4SSatish Balay 
1130606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
11319a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11323a40ed3dSBarry Smith   PetscFunctionReturn(0);
1133d9b7c43dSSatish Balay }
11341c351548SSatish Balay 
11354a2ae208SSatish Balay #undef __FUNCT__
11364a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ"
113787828ca2SBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
11382d61bbb3SSatish Balay {
11392d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
11402d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1141273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
11422d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1143549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
1144273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
11453f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
11462d61bbb3SSatish Balay 
11472d61bbb3SSatish Balay   PetscFunctionBegin;
11482d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
11492d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
11505ef9f2a5SBarry Smith     if (row < 0) continue;
1151aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1152273d9f13SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
11532d61bbb3SSatish Balay #endif
11542d61bbb3SSatish Balay     rp   = aj + ai[brow];
11552d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
11562d61bbb3SSatish Balay     rmax = imax[brow];
11572d61bbb3SSatish Balay     nrow = ailen[brow];
11582d61bbb3SSatish Balay     low  = 0;
11592d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
11605ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1161aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1162273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
11632d61bbb3SSatish Balay #endif
11642d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
11652d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
11662d61bbb3SSatish Balay       if (roworiented) {
11675ef9f2a5SBarry Smith         value = v[l + k*n];
11682d61bbb3SSatish Balay       } else {
11692d61bbb3SSatish Balay         value = v[k + l*m];
11702d61bbb3SSatish Balay       }
11712d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
11722d61bbb3SSatish Balay       while (high-low > 7) {
11732d61bbb3SSatish Balay         t = (low+high)/2;
11742d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
11752d61bbb3SSatish Balay         else              low  = t;
11762d61bbb3SSatish Balay       }
11772d61bbb3SSatish Balay       for (i=low; i<high; i++) {
11782d61bbb3SSatish Balay         if (rp[i] > bcol) break;
11792d61bbb3SSatish Balay         if (rp[i] == bcol) {
11802d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
11812d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
11822d61bbb3SSatish Balay           else                  *bap  = value;
11832d61bbb3SSatish Balay           goto noinsert1;
11842d61bbb3SSatish Balay         }
11852d61bbb3SSatish Balay       }
11862d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
118729bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
11882d61bbb3SSatish Balay       if (nrow >= rmax) {
11892d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
11902d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
11913f1db9ecSBarry Smith         MatScalar *new_a;
11922d61bbb3SSatish Balay 
119329bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
11942d61bbb3SSatish Balay 
11952d61bbb3SSatish Balay         /* Malloc new storage space */
11963f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1197b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
11982d61bbb3SSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
11992d61bbb3SSatish Balay         new_i   = new_j + new_nz;
12002d61bbb3SSatish Balay 
12012d61bbb3SSatish Balay         /* copy over old data into new slots */
12022d61bbb3SSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
12032d61bbb3SSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1204549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
12052d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1206549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1207549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1208549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1209549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
12102d61bbb3SSatish Balay         /* free up old matrix storage */
1211606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1212606d414cSSatish Balay         if (!a->singlemalloc) {
1213606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1214606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1215606d414cSSatish Balay         }
12162d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1217c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
12182d61bbb3SSatish Balay 
12192d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
12202d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1221b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
12222d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
12232d61bbb3SSatish Balay         a->reallocs++;
12242d61bbb3SSatish Balay         a->nz++;
12252d61bbb3SSatish Balay       }
12262d61bbb3SSatish Balay       N = nrow++ - 1;
12272d61bbb3SSatish Balay       /* shift up all the later entries in this row */
12282d61bbb3SSatish Balay       for (ii=N; ii>=i; ii--) {
12292d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1230549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
12312d61bbb3SSatish Balay       }
1232549d3d68SSatish Balay       if (N>=i) {
1233549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1234549d3d68SSatish Balay       }
12352d61bbb3SSatish Balay       rp[i]                      = bcol;
12362d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
12372d61bbb3SSatish Balay       noinsert1:;
12382d61bbb3SSatish Balay       low = i;
12392d61bbb3SSatish Balay     }
12402d61bbb3SSatish Balay     ailen[brow] = nrow;
12412d61bbb3SSatish Balay   }
12422d61bbb3SSatish Balay   PetscFunctionReturn(0);
12432d61bbb3SSatish Balay }
12442d61bbb3SSatish Balay 
12452d61bbb3SSatish Balay 
12464a2ae208SSatish Balay #undef __FUNCT__
12474a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ"
12485ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
12492d61bbb3SSatish Balay {
12502d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
12512d61bbb3SSatish Balay   Mat         outA;
12522d61bbb3SSatish Balay   int         ierr;
1253667159a5SBarry Smith   PetscTruth  row_identity,col_identity;
12542d61bbb3SSatish Balay 
12552d61bbb3SSatish Balay   PetscFunctionBegin;
1256d3d32019SBarry Smith   if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
1257667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1258667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1259667159a5SBarry Smith   if (!row_identity || !col_identity) {
126029bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1261667159a5SBarry Smith   }
12622d61bbb3SSatish Balay 
12632d61bbb3SSatish Balay   outA          = inA;
12642d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
12652d61bbb3SSatish Balay 
12662d61bbb3SSatish Balay   if (!a->diag) {
1267c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
12682d61bbb3SSatish Balay   }
1269cf242676SKris Buschelman 
1270c38d4ed2SBarry Smith   a->row        = row;
1271c38d4ed2SBarry Smith   a->col        = col;
1272c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1273c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1274c38d4ed2SBarry Smith 
1275c38d4ed2SBarry Smith   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
12764c49b128SBarry Smith   ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1277b0a32e0cSBarry Smith   PetscLogObjectParent(inA,a->icol);
1278c38d4ed2SBarry Smith 
1279cf242676SKris Buschelman   /*
1280cf242676SKris Buschelman       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1281cf242676SKris Buschelman       for ILU(0) factorization with natural ordering
1282cf242676SKris Buschelman   */
1283cf242676SKris Buschelman   if (a->bs < 8) {
1284cf242676SKris Buschelman     ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr);
1285cf242676SKris Buschelman   } else {
1286c38d4ed2SBarry Smith     if (!a->solve_work) {
128787828ca2SBarry Smith       ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
128887828ca2SBarry Smith       PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar));
1289c38d4ed2SBarry Smith     }
12902d61bbb3SSatish Balay   }
12912d61bbb3SSatish Balay 
1292667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1293667159a5SBarry Smith 
12942d61bbb3SSatish Balay   PetscFunctionReturn(0);
12952d61bbb3SSatish Balay }
12964a2ae208SSatish Balay #undef __FUNCT__
12974a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ"
1298ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1299ba4ca20aSSatish Balay {
1300c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1301ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1302d132466eSBarry Smith   int               ierr;
1303ba4ca20aSSatish Balay 
13043a40ed3dSBarry Smith   PetscFunctionBegin;
1305c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1306d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1307d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
13083a40ed3dSBarry Smith   PetscFunctionReturn(0);
1309ba4ca20aSSatish Balay }
1310d9b7c43dSSatish Balay 
1311fb2e594dSBarry Smith EXTERN_C_BEGIN
13124a2ae208SSatish Balay #undef __FUNCT__
13134a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
131427a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
131527a8da17SBarry Smith {
131627a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
131714db4f74SSatish Balay   int         i,nz,nbs;
131827a8da17SBarry Smith 
131927a8da17SBarry Smith   PetscFunctionBegin;
132014db4f74SSatish Balay   nz  = baij->maxnz/baij->bs2;
132114db4f74SSatish Balay   nbs = baij->nbs;
132227a8da17SBarry Smith   for (i=0; i<nz; i++) {
132327a8da17SBarry Smith     baij->j[i] = indices[i];
132427a8da17SBarry Smith   }
132527a8da17SBarry Smith   baij->nz = nz;
132614db4f74SSatish Balay   for (i=0; i<nbs; i++) {
132727a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
132827a8da17SBarry Smith   }
132927a8da17SBarry Smith 
133027a8da17SBarry Smith   PetscFunctionReturn(0);
133127a8da17SBarry Smith }
1332fb2e594dSBarry Smith EXTERN_C_END
133327a8da17SBarry Smith 
13344a2ae208SSatish Balay #undef __FUNCT__
13354a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices"
133627a8da17SBarry Smith /*@
133727a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
133827a8da17SBarry Smith        in the matrix.
133927a8da17SBarry Smith 
134027a8da17SBarry Smith   Input Parameters:
134127a8da17SBarry Smith +  mat - the SeqBAIJ matrix
134227a8da17SBarry Smith -  indices - the column indices
134327a8da17SBarry Smith 
134415091d37SBarry Smith   Level: advanced
134515091d37SBarry Smith 
134627a8da17SBarry Smith   Notes:
134727a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
134827a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
134927a8da17SBarry Smith   of the MatSetValues() operation.
135027a8da17SBarry Smith 
135127a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
135227a8da17SBarry Smith   MatCreateSeqBAIJ().
135327a8da17SBarry Smith 
135427a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
135527a8da17SBarry Smith 
135627a8da17SBarry Smith @*/
135727a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
135827a8da17SBarry Smith {
135927a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
136027a8da17SBarry Smith 
136127a8da17SBarry Smith   PetscFunctionBegin;
136227a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1363c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
136427a8da17SBarry Smith   if (f) {
136527a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
136627a8da17SBarry Smith   } else {
136729bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
136827a8da17SBarry Smith   }
136927a8da17SBarry Smith   PetscFunctionReturn(0);
137027a8da17SBarry Smith }
137127a8da17SBarry Smith 
13724a2ae208SSatish Balay #undef __FUNCT__
13734a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ"
1374273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v)
1375273d9f13SBarry Smith {
1376273d9f13SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1377273d9f13SBarry Smith   int          ierr,i,j,n,row,bs,*ai,*aj,mbs;
1378273d9f13SBarry Smith   PetscReal    atmp;
137987828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
1380273d9f13SBarry Smith   MatScalar    *aa;
1381273d9f13SBarry Smith   int          ncols,brow,krow,kcol;
1382273d9f13SBarry Smith 
1383273d9f13SBarry Smith   PetscFunctionBegin;
1384273d9f13SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1385273d9f13SBarry Smith   bs   = a->bs;
1386273d9f13SBarry Smith   aa   = a->a;
1387273d9f13SBarry Smith   ai   = a->i;
1388273d9f13SBarry Smith   aj   = a->j;
1389273d9f13SBarry Smith   mbs = a->mbs;
1390273d9f13SBarry Smith 
1391273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1392273d9f13SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1393273d9f13SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1394273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1395273d9f13SBarry Smith   for (i=0; i<mbs; i++) {
1396273d9f13SBarry Smith     ncols = ai[1] - ai[0]; ai++;
1397273d9f13SBarry Smith     brow  = bs*i;
1398273d9f13SBarry Smith     for (j=0; j<ncols; j++){
1399273d9f13SBarry Smith       /* bcol = bs*(*aj); */
1400273d9f13SBarry Smith       for (kcol=0; kcol<bs; kcol++){
1401273d9f13SBarry Smith         for (krow=0; krow<bs; krow++){
1402273d9f13SBarry Smith           atmp = PetscAbsScalar(*aa); aa++;
1403273d9f13SBarry Smith           row = brow + krow;    /* row index */
1404273d9f13SBarry Smith           /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */
1405273d9f13SBarry Smith           if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp;
1406273d9f13SBarry Smith         }
1407273d9f13SBarry Smith       }
1408273d9f13SBarry Smith       aj++;
1409273d9f13SBarry Smith     }
1410273d9f13SBarry Smith   }
1411273d9f13SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1412273d9f13SBarry Smith   PetscFunctionReturn(0);
1413273d9f13SBarry Smith }
1414273d9f13SBarry Smith 
14154a2ae208SSatish Balay #undef __FUNCT__
14164a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ"
1417273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A)
1418273d9f13SBarry Smith {
1419273d9f13SBarry Smith   int        ierr;
1420273d9f13SBarry Smith 
1421273d9f13SBarry Smith   PetscFunctionBegin;
1422273d9f13SBarry Smith   ierr =  MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1423273d9f13SBarry Smith   PetscFunctionReturn(0);
1424273d9f13SBarry Smith }
1425273d9f13SBarry Smith 
14264a2ae208SSatish Balay #undef __FUNCT__
14274a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ"
142887828ca2SBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar **array)
1429f2a5309cSSatish Balay {
1430f2a5309cSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1431f2a5309cSSatish Balay   PetscFunctionBegin;
1432f2a5309cSSatish Balay   *array = a->a;
1433f2a5309cSSatish Balay   PetscFunctionReturn(0);
1434f2a5309cSSatish Balay }
1435f2a5309cSSatish Balay 
14364a2ae208SSatish Balay #undef __FUNCT__
14374a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ"
143887828ca2SBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar **array)
1439f2a5309cSSatish Balay {
1440f2a5309cSSatish Balay   PetscFunctionBegin;
1441f2a5309cSSatish Balay   PetscFunctionReturn(0);
1442f2a5309cSSatish Balay }
1443f2a5309cSSatish Balay 
144442ee4b1aSHong Zhang #include "petscblaslapack.h"
144542ee4b1aSHong Zhang #undef __FUNCT__
144642ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqBAIJ"
144742ee4b1aSHong Zhang int MatAXPY_SeqBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
144842ee4b1aSHong Zhang {
144942ee4b1aSHong Zhang   int          ierr,one=1;
145042ee4b1aSHong Zhang   Mat_SeqBAIJ  *x  = (Mat_SeqBAIJ *)X->data,*y = (Mat_SeqBAIJ *)Y->data;
145142ee4b1aSHong Zhang 
145242ee4b1aSHong Zhang   PetscFunctionBegin;
145342ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
145442ee4b1aSHong Zhang     BLaxpy_(&x->nz,a,x->a,&one,y->a,&one);
145542ee4b1aSHong Zhang   } else {
145642ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
145742ee4b1aSHong Zhang   }
145842ee4b1aSHong Zhang   PetscFunctionReturn(0);
145942ee4b1aSHong Zhang }
146042ee4b1aSHong Zhang 
14612593348eSBarry Smith /* -------------------------------------------------------------------*/
1462cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1463cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1464cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1465cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1466cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
14677c922b88SBarry Smith        MatMultTranspose_SeqBAIJ,
14687c922b88SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1469cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1470cc2dc46cSBarry Smith        0,
1471cc2dc46cSBarry Smith        0,
1472cc2dc46cSBarry Smith        0,
1473cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1474cc2dc46cSBarry Smith        0,
1475b6490206SBarry Smith        0,
1476f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1477cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1478cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1479cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1480cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1481cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1482b6490206SBarry Smith        0,
1483cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1484cc2dc46cSBarry Smith        0,
1485cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1486cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1487cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1488cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1489cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1490cc2dc46cSBarry Smith        0,
1491cc2dc46cSBarry Smith        0,
1492273d9f13SBarry Smith        MatSetUpPreallocation_SeqBAIJ,
1493cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1494cc2dc46cSBarry Smith        0,
1495f2a5309cSSatish Balay        MatGetArray_SeqBAIJ,
1496f2a5309cSSatish Balay        MatRestoreArray_SeqBAIJ,
14972e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1498cc2dc46cSBarry Smith        0,
1499cc2dc46cSBarry Smith        0,
1500cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1501cc2dc46cSBarry Smith        0,
150242ee4b1aSHong Zhang        MatAXPY_SeqBAIJ,
1503cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1504cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1505cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1506cc2dc46cSBarry Smith        0,
1507cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1508cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1509cc2dc46cSBarry Smith        0,
1510cc2dc46cSBarry Smith        0,
1511cc2dc46cSBarry Smith        0,
1512cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
15133b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
151492c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1515cc2dc46cSBarry Smith        0,
1516cc2dc46cSBarry Smith        0,
1517cc2dc46cSBarry Smith        0,
1518cc2dc46cSBarry Smith        0,
1519cc2dc46cSBarry Smith        0,
1520cc2dc46cSBarry Smith        0,
1521d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1522cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1523b9b97703SBarry Smith        MatDestroy_SeqBAIJ,
1524b9b97703SBarry Smith        MatView_SeqBAIJ,
15253a64cc32SBarry Smith        MatGetPetscMaps_Petsc,
1526273d9f13SBarry Smith        0,
1527273d9f13SBarry Smith        0,
1528273d9f13SBarry Smith        0,
1529273d9f13SBarry Smith        0,
1530273d9f13SBarry Smith        0,
1531273d9f13SBarry Smith        0,
1532273d9f13SBarry Smith        MatGetRowMax_SeqBAIJ,
1533273d9f13SBarry Smith        MatConvert_Basic};
15342593348eSBarry Smith 
15353e90b805SBarry Smith EXTERN_C_BEGIN
15364a2ae208SSatish Balay #undef __FUNCT__
15374a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ"
15383e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
15393e90b805SBarry Smith {
15403e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1541273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1542d132466eSBarry Smith   int         ierr;
15433e90b805SBarry Smith 
15443e90b805SBarry Smith   PetscFunctionBegin;
15453e90b805SBarry Smith   if (aij->nonew != 1) {
154629bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15473e90b805SBarry Smith   }
15483e90b805SBarry Smith 
15493e90b805SBarry Smith   /* allocate space for values if not already there */
15503e90b805SBarry Smith   if (!aij->saved_values) {
155187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
15523e90b805SBarry Smith   }
15533e90b805SBarry Smith 
15543e90b805SBarry Smith   /* copy values over */
155587828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
15563e90b805SBarry Smith   PetscFunctionReturn(0);
15573e90b805SBarry Smith }
15583e90b805SBarry Smith EXTERN_C_END
15593e90b805SBarry Smith 
15603e90b805SBarry Smith EXTERN_C_BEGIN
15614a2ae208SSatish Balay #undef __FUNCT__
15624a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ"
156332e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
15643e90b805SBarry Smith {
15653e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1566273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
15673e90b805SBarry Smith 
15683e90b805SBarry Smith   PetscFunctionBegin;
15693e90b805SBarry Smith   if (aij->nonew != 1) {
157029bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
15713e90b805SBarry Smith   }
15723e90b805SBarry Smith   if (!aij->saved_values) {
157329bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
15743e90b805SBarry Smith   }
15753e90b805SBarry Smith 
15763e90b805SBarry Smith   /* copy values over */
157787828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
15783e90b805SBarry Smith   PetscFunctionReturn(0);
15793e90b805SBarry Smith }
15803e90b805SBarry Smith EXTERN_C_END
15813e90b805SBarry Smith 
1582273d9f13SBarry Smith EXTERN_C_BEGIN
1583273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*);
1584273d9f13SBarry Smith EXTERN_C_END
1585273d9f13SBarry Smith 
1586273d9f13SBarry Smith EXTERN_C_BEGIN
15874a2ae208SSatish Balay #undef __FUNCT__
15884a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ"
1589273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B)
15902593348eSBarry Smith {
1591273d9f13SBarry Smith   int         ierr,size;
1592b6490206SBarry Smith   Mat_SeqBAIJ *b;
15933b2fbd54SBarry Smith 
15943a40ed3dSBarry Smith   PetscFunctionBegin;
1595273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
159629bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1597b6490206SBarry Smith 
1598273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1599273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1600b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr);
1601b0a32e0cSBarry Smith   B->data = (void*)b;
1602549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1603549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
16042593348eSBarry Smith   B->factor           = 0;
16052593348eSBarry Smith   B->lupivotthreshold = 1.0;
160690f02eecSBarry Smith   B->mapping          = 0;
16072593348eSBarry Smith   b->row              = 0;
16082593348eSBarry Smith   b->col              = 0;
1609e51c0b9cSSatish Balay   b->icol             = 0;
16102593348eSBarry Smith   b->reallocs         = 0;
16113e90b805SBarry Smith   b->saved_values     = 0;
1612cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
1613563b5814SBarry Smith   b->setvalueslen     = 0;
1614563b5814SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1615563b5814SBarry Smith #endif
16162593348eSBarry Smith 
16173a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
16183a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1619a5ae1ecdSBarry Smith 
1620c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1621c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
16222593348eSBarry Smith   b->nonew            = 0;
16232593348eSBarry Smith   b->diag             = 0;
16242593348eSBarry Smith   b->solve_work       = 0;
1625de6a44a3SBarry Smith   b->mult_work        = 0;
16262a1b7f2aSHong Zhang   B->spptr            = 0;
16270e6d2581SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1628883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
16294e220ebcSLois Curfman McInnes 
1630f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
16313e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1632bc4b532fSSatish Balay                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1633f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
16343e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1635bc4b532fSSatish Balay                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1636f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
163727a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1638bc4b532fSSatish Balay                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1639a6175056SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_seqaij_C",
1640273d9f13SBarry Smith                                      "MatConvert_SeqBAIJ_SeqAIJ",
1641273d9f13SBarry Smith                                       MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr);
16423a40ed3dSBarry Smith   PetscFunctionReturn(0);
16432593348eSBarry Smith }
1644273d9f13SBarry Smith EXTERN_C_END
16452593348eSBarry Smith 
16464a2ae208SSatish Balay #undef __FUNCT__
16474a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ"
16482e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
16492593348eSBarry Smith {
16502593348eSBarry Smith   Mat         C;
1651b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
165227a8da17SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1653de6a44a3SBarry Smith 
16543a40ed3dSBarry Smith   PetscFunctionBegin;
165529bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
16562593348eSBarry Smith 
16572593348eSBarry Smith   *B = 0;
1658273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1659273d9f13SBarry Smith   ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
1660273d9f13SBarry Smith   c    = (Mat_SeqBAIJ*)C->data;
166144cd7ae7SLois Curfman McInnes 
1662b6490206SBarry Smith   c->bs         = a->bs;
16639df24120SSatish Balay   c->bs2        = a->bs2;
1664b6490206SBarry Smith   c->mbs        = a->mbs;
1665b6490206SBarry Smith   c->nbs        = a->nbs;
166635d8aa7fSBarry Smith   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
16672593348eSBarry Smith 
1668b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1669b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
1670b6490206SBarry Smith   for (i=0; i<mbs; i++) {
16712593348eSBarry Smith     c->imax[i] = a->imax[i];
16722593348eSBarry Smith     c->ilen[i] = a->ilen[i];
16732593348eSBarry Smith   }
16742593348eSBarry Smith 
16752593348eSBarry Smith   /* allocate the matrix space */
1676c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
16773f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1678b0a32e0cSBarry Smith   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
16797e67e3f9SSatish Balay   c->j = (int*)(c->a + nz*bs2);
1680de6a44a3SBarry Smith   c->i = c->j + nz;
1681549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1682b6490206SBarry Smith   if (mbs > 0) {
1683549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
16842e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1685549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
16862e8a6d31SBarry Smith     } else {
1687549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
16882593348eSBarry Smith     }
16892593348eSBarry Smith   }
16902593348eSBarry Smith 
1691b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
16922593348eSBarry Smith   c->sorted      = a->sorted;
16932593348eSBarry Smith   c->roworiented = a->roworiented;
16942593348eSBarry Smith   c->nonew       = a->nonew;
16952593348eSBarry Smith 
16962593348eSBarry Smith   if (a->diag) {
1697b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1698b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1699b6490206SBarry Smith     for (i=0; i<mbs; i++) {
17002593348eSBarry Smith       c->diag[i] = a->diag[i];
17012593348eSBarry Smith     }
170298305bb5SBarry Smith   } else c->diag        = 0;
17032593348eSBarry Smith   c->nz                 = a->nz;
17042593348eSBarry Smith   c->maxnz              = a->maxnz;
17052593348eSBarry Smith   c->solve_work         = 0;
17062a1b7f2aSHong Zhang   C->spptr              = 0;     /* Dangerous -I'm throwing away a->spptr */
17077fc0212eSBarry Smith   c->mult_work          = 0;
1708273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
1709273d9f13SBarry Smith   C->assembled          = PETSC_TRUE;
17102593348eSBarry Smith   *B = C;
1711b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
17123a40ed3dSBarry Smith   PetscFunctionReturn(0);
17132593348eSBarry Smith }
17142593348eSBarry Smith 
1715273d9f13SBarry Smith EXTERN_C_BEGIN
17164a2ae208SSatish Balay #undef __FUNCT__
17174a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ"
1718b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A)
17192593348eSBarry Smith {
1720b6490206SBarry Smith   Mat_SeqBAIJ  *a;
17212593348eSBarry Smith   Mat          B;
1722f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1723b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
172435aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1725a7c10996SSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
172687828ca2SBarry Smith   PetscScalar  *aa;
172719bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1728*ecc4ab6dSBarry Smith #if defined(PETSC_HAVE_DSCPACK)
1729*ecc4ab6dSBarry Smith   PetscTruth   flag;
1730*ecc4ab6dSBarry Smith #endif
17312593348eSBarry Smith 
17323a40ed3dSBarry Smith   PetscFunctionBegin;
1733b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1734de6a44a3SBarry Smith   bs2  = bs*bs;
1735b6490206SBarry Smith 
1736d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
173729bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1738b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
17390752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1740552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
17412593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
17422593348eSBarry Smith 
1743d64ed03dSBarry Smith   if (header[3] < 0) {
174429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ");
1745d64ed03dSBarry Smith   }
1746d64ed03dSBarry Smith 
174729bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
174835aab85fSBarry Smith 
174935aab85fSBarry Smith   /*
175035aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
175135aab85fSBarry Smith     divisible by the blocksize
175235aab85fSBarry Smith   */
1753b6490206SBarry Smith   mbs        = M/bs;
175435aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
175535aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
175635aab85fSBarry Smith   else                  mbs++;
175735aab85fSBarry Smith   if (extra_rows) {
1758b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
175935aab85fSBarry Smith   }
1760b6490206SBarry Smith 
17612593348eSBarry Smith   /* read in row lengths */
1762b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
17630752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
176435aab85fSBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
17652593348eSBarry Smith 
1766b6490206SBarry Smith   /* read in column indices */
1767b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
17680752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
176935aab85fSBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1770b6490206SBarry Smith 
1771b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1772b0a32e0cSBarry Smith   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
1773549d3d68SSatish Balay   ierr     = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1774b0a32e0cSBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
1775549d3d68SSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
177635aab85fSBarry Smith   masked   = mask + mbs;
1777b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1778b6490206SBarry Smith   for (i=0; i<mbs; i++) {
177935aab85fSBarry Smith     nmask = 0;
1780b6490206SBarry Smith     for (j=0; j<bs; j++) {
1781b6490206SBarry Smith       kmax = rowlengths[rowcount];
1782b6490206SBarry Smith       for (k=0; k<kmax; k++) {
178335aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
178435aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1785b6490206SBarry Smith       }
1786b6490206SBarry Smith       rowcount++;
1787b6490206SBarry Smith     }
178835aab85fSBarry Smith     browlengths[i] += nmask;
178935aab85fSBarry Smith     /* zero out the mask elements we set */
179035aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1791b6490206SBarry Smith   }
1792b6490206SBarry Smith 
17932593348eSBarry Smith   /* create our matrix */
17943eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
17952593348eSBarry Smith   B = *A;
1796b6490206SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
17972593348eSBarry Smith 
17982593348eSBarry Smith   /* set matrix "i" values */
1799de6a44a3SBarry Smith   a->i[0] = 0;
1800b6490206SBarry Smith   for (i=1; i<= mbs; i++) {
1801b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1802b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
18032593348eSBarry Smith   }
1804b6490206SBarry Smith   a->nz         = 0;
1805b6490206SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
18062593348eSBarry Smith 
1807b6490206SBarry Smith   /* read in nonzero values */
180887828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
18090752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
181035aab85fSBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1811b6490206SBarry Smith 
1812b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1813b6490206SBarry Smith   nzcount = 0; jcount = 0;
1814b6490206SBarry Smith   for (i=0; i<mbs; i++) {
1815b6490206SBarry Smith     nzcountb = nzcount;
181635aab85fSBarry Smith     nmask    = 0;
1817b6490206SBarry Smith     for (j=0; j<bs; j++) {
1818b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1819b6490206SBarry Smith       for (k=0; k<kmax; k++) {
182035aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
182135aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1822b6490206SBarry Smith       }
1823b6490206SBarry Smith     }
1824de6a44a3SBarry Smith     /* sort the masked values */
1825433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1826de6a44a3SBarry Smith 
1827b6490206SBarry Smith     /* set "j" values into matrix */
1828b6490206SBarry Smith     maskcount = 1;
182935aab85fSBarry Smith     for (j=0; j<nmask; j++) {
183035aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1831de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1832b6490206SBarry Smith     }
1833b6490206SBarry Smith     /* set "a" values into matrix */
1834de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1835b6490206SBarry Smith     for (j=0; j<bs; j++) {
1836b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1837b6490206SBarry Smith       for (k=0; k<kmax; k++) {
1838de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1839de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1840de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1841de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1842375fe846SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1843b6490206SBarry Smith       }
1844b6490206SBarry Smith     }
184535aab85fSBarry Smith     /* zero out the mask elements we set */
184635aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1847b6490206SBarry Smith   }
184829bbc08cSBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
1849b6490206SBarry Smith 
1850606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1851606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1852606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
1853606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
1854606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
1855b6490206SBarry Smith 
1856b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1857de6a44a3SBarry Smith 
1858*ecc4ab6dSBarry Smith #if defined(PETSC_HAVE_DSCPACK)
1859*ecc4ab6dSBarry Smith   ierr = PetscOptionsHasName(B->prefix,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
1860*ecc4ab6dSBarry Smith   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(B);CHKERRQ(ierr); }
1861*ecc4ab6dSBarry Smith #endif
1862*ecc4ab6dSBarry Smith 
18639c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
18643a40ed3dSBarry Smith   PetscFunctionReturn(0);
18652593348eSBarry Smith }
1866273d9f13SBarry Smith EXTERN_C_END
18672593348eSBarry Smith 
18684a2ae208SSatish Balay #undef __FUNCT__
18694a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ"
1870273d9f13SBarry Smith /*@C
1871273d9f13SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1872273d9f13SBarry Smith    compressed row) format.  For good matrix assembly performance the
1873273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1874273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1875273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
18762593348eSBarry Smith 
1877273d9f13SBarry Smith    Collective on MPI_Comm
1878273d9f13SBarry Smith 
1879273d9f13SBarry Smith    Input Parameters:
1880273d9f13SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
1881273d9f13SBarry Smith .  bs - size of block
1882273d9f13SBarry Smith .  m - number of rows
1883273d9f13SBarry Smith .  n - number of columns
188435d8aa7fSBarry Smith .  nz - number of nonzero blocks  per block row (same for all rows)
188535d8aa7fSBarry Smith -  nnz - array containing the number of nonzero blocks in the various block rows
1886273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1887273d9f13SBarry Smith 
1888273d9f13SBarry Smith    Output Parameter:
1889273d9f13SBarry Smith .  A - the matrix
1890273d9f13SBarry Smith 
1891273d9f13SBarry Smith    Options Database Keys:
1892273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1893273d9f13SBarry Smith                      block calculations (much slower)
1894273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1895273d9f13SBarry Smith 
1896273d9f13SBarry Smith    Level: intermediate
1897273d9f13SBarry Smith 
1898273d9f13SBarry Smith    Notes:
189935d8aa7fSBarry Smith    A nonzero block is any block that as 1 or more nonzeros in it
190035d8aa7fSBarry Smith 
1901273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1902273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1903273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1904273d9f13SBarry Smith 
1905273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1906273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1907273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1908273d9f13SBarry Smith    matrices.
1909273d9f13SBarry Smith 
1910273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1911273d9f13SBarry Smith @*/
1912273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1913273d9f13SBarry Smith {
1914273d9f13SBarry Smith   int ierr;
1915273d9f13SBarry Smith 
1916273d9f13SBarry Smith   PetscFunctionBegin;
1917273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1918273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr);
1919273d9f13SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1920273d9f13SBarry Smith   PetscFunctionReturn(0);
1921273d9f13SBarry Smith }
1922273d9f13SBarry Smith 
19234a2ae208SSatish Balay #undef __FUNCT__
19244a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation"
1925273d9f13SBarry Smith /*@C
1926273d9f13SBarry Smith    MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros
1927273d9f13SBarry Smith    per row in the matrix. For good matrix assembly performance the
1928273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1929273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1930273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1931273d9f13SBarry Smith 
1932273d9f13SBarry Smith    Collective on MPI_Comm
1933273d9f13SBarry Smith 
1934273d9f13SBarry Smith    Input Parameters:
1935273d9f13SBarry Smith +  A - the matrix
1936273d9f13SBarry Smith .  bs - size of block
1937273d9f13SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1938273d9f13SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
1939273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1940273d9f13SBarry Smith 
1941273d9f13SBarry Smith    Options Database Keys:
1942273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1943273d9f13SBarry Smith                      block calculations (much slower)
1944273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1945273d9f13SBarry Smith 
1946273d9f13SBarry Smith    Level: intermediate
1947273d9f13SBarry Smith 
1948273d9f13SBarry Smith    Notes:
1949273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1950273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1951273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1952273d9f13SBarry Smith 
1953273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1954273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1955273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1956273d9f13SBarry Smith    matrices.
1957273d9f13SBarry Smith 
1958273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1959273d9f13SBarry Smith @*/
1960273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
1961273d9f13SBarry Smith {
1962273d9f13SBarry Smith   Mat_SeqBAIJ *b;
1963273d9f13SBarry Smith   int         i,len,ierr,mbs,nbs,bs2,newbs = bs;
1964273d9f13SBarry Smith   PetscTruth  flg;
1965273d9f13SBarry Smith 
1966273d9f13SBarry Smith   PetscFunctionBegin;
1967273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr);
1968273d9f13SBarry Smith   if (!flg) PetscFunctionReturn(0);
1969273d9f13SBarry Smith 
1970273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1971b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1972273d9f13SBarry Smith   if (nnz && newbs != bs) {
1973273d9f13SBarry Smith     SETERRQ(1,"Cannot change blocksize from command line if setting nnz");
1974273d9f13SBarry Smith   }
1975273d9f13SBarry Smith   bs   = newbs;
1976273d9f13SBarry Smith 
1977273d9f13SBarry Smith   mbs  = B->m/bs;
1978273d9f13SBarry Smith   nbs  = B->n/bs;
1979273d9f13SBarry Smith   bs2  = bs*bs;
1980273d9f13SBarry Smith 
1981273d9f13SBarry Smith   if (mbs*bs!=B->m || nbs*bs!=B->n) {
198235d8aa7fSBarry Smith     SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs);
1983273d9f13SBarry Smith   }
1984273d9f13SBarry Smith 
1985435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1986435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
1987273d9f13SBarry Smith   if (nnz) {
1988273d9f13SBarry Smith     for (i=0; i<mbs; i++) {
1989273d9f13SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
19903a7fca6bSBarry Smith       if (nnz[i] > nbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %d value %d rowlength %d",i,nnz[i],nbs);
1991273d9f13SBarry Smith     }
1992273d9f13SBarry Smith   }
1993273d9f13SBarry Smith 
1994273d9f13SBarry Smith   b       = (Mat_SeqBAIJ*)B->data;
1995b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
199654138f6bSKris Buschelman   B->ops->solve               = MatSolve_SeqBAIJ_Update;
199754138f6bSKris Buschelman   B->ops->solvetranspose      = MatSolveTranspose_SeqBAIJ_Update;
1998273d9f13SBarry Smith   if (!flg) {
1999273d9f13SBarry Smith     switch (bs) {
2000273d9f13SBarry Smith     case 1:
2001273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
2002273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
2003273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
2004273d9f13SBarry Smith       break;
2005273d9f13SBarry Smith     case 2:
2006273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
2007273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
2008273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
2009273d9f13SBarry Smith       break;
2010273d9f13SBarry Smith     case 3:
2011273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
2012273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
2013273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
2014273d9f13SBarry Smith       break;
2015273d9f13SBarry Smith     case 4:
2016273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
2017273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
2018273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
2019273d9f13SBarry Smith       break;
2020273d9f13SBarry Smith     case 5:
2021273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
2022273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
2023273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
2024273d9f13SBarry Smith       break;
2025273d9f13SBarry Smith     case 6:
2026273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
2027273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
2028273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
2029273d9f13SBarry Smith       break;
2030273d9f13SBarry Smith     case 7:
203154138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
2032273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
2033273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
2034273d9f13SBarry Smith       break;
2035273d9f13SBarry Smith     default:
203654138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
2037273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_N;
2038273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_N;
2039273d9f13SBarry Smith       break;
2040273d9f13SBarry Smith     }
2041273d9f13SBarry Smith   }
2042273d9f13SBarry Smith   b->bs      = bs;
2043273d9f13SBarry Smith   b->mbs     = mbs;
2044273d9f13SBarry Smith   b->nbs     = nbs;
2045b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
2046273d9f13SBarry Smith   if (!nnz) {
2047435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
2048273d9f13SBarry Smith     else if (nz <= 0)        nz = 1;
2049273d9f13SBarry Smith     for (i=0; i<mbs; i++) b->imax[i] = nz;
2050273d9f13SBarry Smith     nz = nz*mbs;
2051273d9f13SBarry Smith   } else {
2052273d9f13SBarry Smith     nz = 0;
2053273d9f13SBarry Smith     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
2054273d9f13SBarry Smith   }
2055273d9f13SBarry Smith 
2056273d9f13SBarry Smith   /* allocate the matrix space */
2057273d9f13SBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
2058b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&b->a);CHKERRQ(ierr);
2059273d9f13SBarry Smith   ierr  = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
2060273d9f13SBarry Smith   b->j  = (int*)(b->a + nz*bs2);
2061273d9f13SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
2062273d9f13SBarry Smith   b->i  = b->j + nz;
2063273d9f13SBarry Smith   b->singlemalloc = PETSC_TRUE;
2064273d9f13SBarry Smith 
2065273d9f13SBarry Smith   b->i[0] = 0;
2066273d9f13SBarry Smith   for (i=1; i<mbs+1; i++) {
2067273d9f13SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
2068273d9f13SBarry Smith   }
2069273d9f13SBarry Smith 
2070273d9f13SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
207116d6aa21SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
2072b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
2073273d9f13SBarry Smith   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
2074273d9f13SBarry Smith 
2075273d9f13SBarry Smith   b->bs               = bs;
2076273d9f13SBarry Smith   b->bs2              = bs2;
2077273d9f13SBarry Smith   b->mbs              = mbs;
2078273d9f13SBarry Smith   b->nz               = 0;
2079273d9f13SBarry Smith   b->maxnz            = nz*bs2;
2080273d9f13SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
2081273d9f13SBarry Smith   PetscFunctionReturn(0);
2082273d9f13SBarry Smith }
20832593348eSBarry Smith 
2084