xref: /petsc/src/mat/impls/baij/seq/baij.c (revision 04929863ded087c2368f8050cef69e2db0544a05)
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 
1295d5f7c2SBarry Smith /*  UGLY, ugly, ugly
1387828ca2SBarry Smith    When MatScalar == PetscScalar the function MatSetValuesBlocked_SeqBAIJ_MatScalar() does
143477af42SKris Buschelman    not exist. Otherwise ..._MatScalar() takes matrix dlements in single precision and
1595d5f7c2SBarry Smith    inserts them into the single precision data structure. The function MatSetValuesBlocked_SeqBAIJ()
1695d5f7c2SBarry Smith    converts the entries into single precision and then calls ..._MatScalar() to put them
1795d5f7c2SBarry Smith    into the single precision data structures.
1895d5f7c2SBarry Smith */
1995d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
20ca44d042SBarry Smith EXTERN int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat,int,int*,int,int*,MatScalar*,InsertMode);
2195d5f7c2SBarry Smith #else
2295d5f7c2SBarry Smith #define MatSetValuesBlocked_SeqBAIJ_MatScalar MatSetValuesBlocked_SeqBAIJ
2395d5f7c2SBarry Smith #endif
24*04929863SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
25*04929863SHong Zhang EXTERN int MatUseDSCPACK_MPIBAIJ(Mat);
26*04929863SHong Zhang #endif
2795d5f7c2SBarry Smith 
282d61bbb3SSatish Balay #define CHUNKSIZE  10
29de6a44a3SBarry Smith 
30be5855fcSBarry Smith /*
31be5855fcSBarry Smith      Checks for missing diagonals
32be5855fcSBarry Smith */
334a2ae208SSatish Balay #undef __FUNCT__
344a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqBAIJ"
35c4992f7dSBarry Smith int MatMissingDiagonal_SeqBAIJ(Mat A)
36be5855fcSBarry Smith {
37be5855fcSBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
38883fce79SBarry Smith   int         *diag,*jj = a->j,i,ierr;
39be5855fcSBarry Smith 
40be5855fcSBarry Smith   PetscFunctionBegin;
41c4992f7dSBarry Smith   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
42883fce79SBarry Smith   diag = a->diag;
430e8e8aceSBarry Smith   for (i=0; i<a->mbs; i++) {
44be5855fcSBarry Smith     if (jj[diag[i]] != i) {
4529bbc08cSBarry Smith       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
46be5855fcSBarry Smith     }
47be5855fcSBarry Smith   }
48be5855fcSBarry Smith   PetscFunctionReturn(0);
49be5855fcSBarry Smith }
50be5855fcSBarry Smith 
514a2ae208SSatish Balay #undef __FUNCT__
524a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqBAIJ"
53c4992f7dSBarry Smith int MatMarkDiagonal_SeqBAIJ(Mat A)
54de6a44a3SBarry Smith {
55de6a44a3SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
5682502324SSatish Balay   int         i,j,*diag,m = a->mbs,ierr;
57de6a44a3SBarry Smith 
583a40ed3dSBarry Smith   PetscFunctionBegin;
59883fce79SBarry Smith   if (a->diag) PetscFunctionReturn(0);
60883fce79SBarry Smith 
61b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
62b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
637fc0212eSBarry Smith   for (i=0; i<m; i++) {
64ecc615b2SBarry Smith     diag[i] = a->i[i+1];
65de6a44a3SBarry Smith     for (j=a->i[i]; j<a->i[i+1]; j++) {
66de6a44a3SBarry Smith       if (a->j[j] == i) {
67de6a44a3SBarry Smith         diag[i] = j;
68de6a44a3SBarry Smith         break;
69de6a44a3SBarry Smith       }
70de6a44a3SBarry Smith     }
71de6a44a3SBarry Smith   }
72de6a44a3SBarry Smith   a->diag = diag;
733a40ed3dSBarry Smith   PetscFunctionReturn(0);
74de6a44a3SBarry Smith }
752593348eSBarry Smith 
762593348eSBarry Smith 
77ca44d042SBarry Smith EXTERN int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
783b2fbd54SBarry Smith 
794a2ae208SSatish Balay #undef __FUNCT__
804a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqBAIJ"
810e6d2581SBarry Smith static int MatGetRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
823b2fbd54SBarry Smith {
833b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
843b2fbd54SBarry Smith   int         ierr,n = a->mbs,i;
853b2fbd54SBarry Smith 
863a40ed3dSBarry Smith   PetscFunctionBegin;
873b2fbd54SBarry Smith   *nn = n;
883a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
893b2fbd54SBarry Smith   if (symmetric) {
903b2fbd54SBarry Smith     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
913b2fbd54SBarry Smith   } else if (oshift == 1) {
923b2fbd54SBarry Smith     /* temporarily add 1 to i and j indices */
93435da068SBarry Smith     int nz = a->i[n];
943b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
953b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]++;
963b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
973b2fbd54SBarry Smith   } else {
983b2fbd54SBarry Smith     *ia = a->i; *ja = a->j;
993b2fbd54SBarry Smith   }
1003b2fbd54SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1023b2fbd54SBarry Smith }
1033b2fbd54SBarry Smith 
1044a2ae208SSatish Balay #undef __FUNCT__
1054a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqBAIJ"
106435da068SBarry Smith static int MatRestoreRowIJ_SeqBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
1073b2fbd54SBarry Smith {
1083b2fbd54SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
109606d414cSSatish Balay   int         i,n = a->mbs,ierr;
1103b2fbd54SBarry Smith 
1113a40ed3dSBarry Smith   PetscFunctionBegin;
1123a40ed3dSBarry Smith   if (!ia) PetscFunctionReturn(0);
1133b2fbd54SBarry Smith   if (symmetric) {
114606d414cSSatish Balay     ierr = PetscFree(*ia);CHKERRQ(ierr);
115606d414cSSatish Balay     ierr = PetscFree(*ja);CHKERRQ(ierr);
116af5da2bfSBarry Smith   } else if (oshift == 1) {
117435da068SBarry Smith     int nz = a->i[n]-1;
1183b2fbd54SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
1193b2fbd54SBarry Smith     for (i=0; i<n+1; i++) a->i[i]--;
1203b2fbd54SBarry Smith   }
1213a40ed3dSBarry Smith   PetscFunctionReturn(0);
1223b2fbd54SBarry Smith }
1233b2fbd54SBarry Smith 
1244a2ae208SSatish Balay #undef __FUNCT__
1254a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqBAIJ"
1262d61bbb3SSatish Balay int MatGetBlockSize_SeqBAIJ(Mat mat,int *bs)
1272d61bbb3SSatish Balay {
1282d61bbb3SSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ*)mat->data;
1292d61bbb3SSatish Balay 
1302d61bbb3SSatish Balay   PetscFunctionBegin;
1312d61bbb3SSatish Balay   *bs = baij->bs;
1322d61bbb3SSatish Balay   PetscFunctionReturn(0);
1332d61bbb3SSatish Balay }
1342d61bbb3SSatish Balay 
1352d61bbb3SSatish Balay 
1364a2ae208SSatish Balay #undef __FUNCT__
1374a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqBAIJ"
138e1311b90SBarry Smith int MatDestroy_SeqBAIJ(Mat A)
1392d61bbb3SSatish Balay {
1402d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
141e51c0b9cSSatish Balay   int         ierr;
1422d61bbb3SSatish Balay 
143433994e6SBarry Smith   PetscFunctionBegin;
144aa482453SBarry Smith #if defined(PETSC_USE_LOG)
145b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, Cols=%d, NZ=%d",A->m,A->n,a->nz);
1462d61bbb3SSatish Balay #endif
147606d414cSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
148606d414cSSatish Balay   if (!a->singlemalloc) {
149606d414cSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
150606d414cSSatish Balay     ierr = PetscFree(a->j);CHKERRQ(ierr);
151606d414cSSatish Balay   }
152c38d4ed2SBarry Smith   if (a->row) {
153c38d4ed2SBarry Smith     ierr = ISDestroy(a->row);CHKERRQ(ierr);
154c38d4ed2SBarry Smith   }
155c38d4ed2SBarry Smith   if (a->col) {
156c38d4ed2SBarry Smith     ierr = ISDestroy(a->col);CHKERRQ(ierr);
157c38d4ed2SBarry Smith   }
158606d414cSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
159606d414cSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
160606d414cSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
161606d414cSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
162606d414cSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
163e51c0b9cSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
164606d414cSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
165563b5814SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
166563b5814SBarry Smith   if (a->setvaluescopy) {ierr = PetscFree(a->setvaluescopy);CHKERRQ(ierr);}
167563b5814SBarry Smith #endif
168606d414cSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
1692d61bbb3SSatish Balay   PetscFunctionReturn(0);
1702d61bbb3SSatish Balay }
1712d61bbb3SSatish Balay 
1724a2ae208SSatish Balay #undef __FUNCT__
1734a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqBAIJ"
1742d61bbb3SSatish Balay int MatSetOption_SeqBAIJ(Mat A,MatOption op)
1752d61bbb3SSatish Balay {
1762d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1772d61bbb3SSatish Balay 
1782d61bbb3SSatish Balay   PetscFunctionBegin;
179aa275fccSKris Buschelman   switch (op) {
180aa275fccSKris Buschelman   case MAT_ROW_ORIENTED:
181aa275fccSKris Buschelman     a->roworiented    = PETSC_TRUE;
182aa275fccSKris Buschelman     break;
183aa275fccSKris Buschelman   case MAT_COLUMN_ORIENTED:
184aa275fccSKris Buschelman     a->roworiented    = PETSC_FALSE;
185aa275fccSKris Buschelman     break;
186aa275fccSKris Buschelman   case MAT_COLUMNS_SORTED:
187aa275fccSKris Buschelman     a->sorted         = PETSC_TRUE;
188aa275fccSKris Buschelman     break;
189aa275fccSKris Buschelman   case MAT_COLUMNS_UNSORTED:
190aa275fccSKris Buschelman     a->sorted         = PETSC_FALSE;
191aa275fccSKris Buschelman     break;
192aa275fccSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
193aa275fccSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
194aa275fccSKris Buschelman     break;
195aa275fccSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
196aa275fccSKris Buschelman     a->nonew          = 1;
197aa275fccSKris Buschelman     break;
198aa275fccSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
199aa275fccSKris Buschelman     a->nonew          = -1;
200aa275fccSKris Buschelman     break;
201aa275fccSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
202aa275fccSKris Buschelman     a->nonew          = -2;
203aa275fccSKris Buschelman     break;
204aa275fccSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
205aa275fccSKris Buschelman     a->nonew          = 0;
206aa275fccSKris Buschelman     break;
207aa275fccSKris Buschelman   case MAT_ROWS_SORTED:
208aa275fccSKris Buschelman   case MAT_ROWS_UNSORTED:
209aa275fccSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
210aa275fccSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
211aa275fccSKris Buschelman   case MAT_USE_HASH_TABLE:
212b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
213aa275fccSKris Buschelman     break;
214aa275fccSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
21529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
216bd648c37SKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
217bd648c37SKris Buschelman     if (a->bs==4) {
21894ee7fc8SKris Buschelman       PetscTruth sse_enabled_local,sse_enabled_global;
219bd648c37SKris Buschelman       int        ierr;
22094ee7fc8SKris Buschelman 
22194ee7fc8SKris Buschelman       sse_enabled_local  = PETSC_FALSE;
22294ee7fc8SKris Buschelman       sse_enabled_global = PETSC_FALSE;
22394ee7fc8SKris Buschelman 
22494ee7fc8SKris Buschelman       ierr = PetscSSEIsEnabled(A->comm,&sse_enabled_local,&sse_enabled_global);CHKERRQ(ierr);
225e64df4b9SKris Buschelman #if defined(PETSC_HAVE_SSE)
22694ee7fc8SKris Buschelman       if (sse_enabled_local) {
22754138f6bSKris Buschelman           a->single_precision_solves = PETSC_TRUE;
22854138f6bSKris Buschelman           A->ops->solve              = MatSolve_SeqBAIJ_Update;
22954138f6bSKris Buschelman           A->ops->solvetranspose     = MatSolveTranspose_SeqBAIJ_Update;
230cf242676SKris Buschelman           PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES set\n");
23154138f6bSKris Buschelman           break;
23294ee7fc8SKris Buschelman       } else {
23394ee7fc8SKris Buschelman         PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
2348661488fSKris Buschelman       }
235e64df4b9SKris Buschelman #else
236bd648c37SKris Buschelman       PetscLogInfo(A,"MatSetOption_SeqBAIJ:Option MAT_USE_SINGLE_PRECISION_SOLVES ignored\n");
237e64df4b9SKris Buschelman #endif
238bd648c37SKris Buschelman     }
239bd648c37SKris Buschelman     break;
240aa275fccSKris Buschelman   default:
24129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
2422d61bbb3SSatish Balay   }
2432d61bbb3SSatish Balay   PetscFunctionReturn(0);
2442d61bbb3SSatish Balay }
2452d61bbb3SSatish Balay 
2464a2ae208SSatish Balay #undef __FUNCT__
2474a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqBAIJ"
24887828ca2SBarry Smith int MatGetRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
2492d61bbb3SSatish Balay {
2502d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
25182502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*idx_i,bs2,ierr;
2523f1db9ecSBarry Smith   MatScalar    *aa,*aa_i;
25387828ca2SBarry Smith   PetscScalar  *v_i;
2542d61bbb3SSatish Balay 
2552d61bbb3SSatish Balay   PetscFunctionBegin;
2562d61bbb3SSatish Balay   bs  = a->bs;
2572d61bbb3SSatish Balay   ai  = a->i;
2582d61bbb3SSatish Balay   aj  = a->j;
2592d61bbb3SSatish Balay   aa  = a->a;
2602d61bbb3SSatish Balay   bs2 = a->bs2;
2612d61bbb3SSatish Balay 
262273d9f13SBarry Smith   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
2632d61bbb3SSatish Balay 
2642d61bbb3SSatish Balay   bn  = row/bs;   /* Block number */
2652d61bbb3SSatish Balay   bp  = row % bs; /* Block Position */
2662d61bbb3SSatish Balay   M   = ai[bn+1] - ai[bn];
2672d61bbb3SSatish Balay   *nz = bs*M;
2682d61bbb3SSatish Balay 
2692d61bbb3SSatish Balay   if (v) {
2702d61bbb3SSatish Balay     *v = 0;
2712d61bbb3SSatish Balay     if (*nz) {
27287828ca2SBarry Smith       ierr = PetscMalloc((*nz)*sizeof(PetscScalar),v);CHKERRQ(ierr);
2732d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
2742d61bbb3SSatish Balay         v_i  = *v + i*bs;
2752d61bbb3SSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
2762d61bbb3SSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
2772d61bbb3SSatish Balay       }
2782d61bbb3SSatish Balay     }
2792d61bbb3SSatish Balay   }
2802d61bbb3SSatish Balay 
2812d61bbb3SSatish Balay   if (idx) {
2822d61bbb3SSatish Balay     *idx = 0;
2832d61bbb3SSatish Balay     if (*nz) {
284b0a32e0cSBarry Smith       ierr = PetscMalloc((*nz)*sizeof(int),idx);CHKERRQ(ierr);
2852d61bbb3SSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
2862d61bbb3SSatish Balay         idx_i = *idx + i*bs;
2872d61bbb3SSatish Balay         itmp  = bs*aj[ai[bn] + i];
2882d61bbb3SSatish Balay         for (j=0; j<bs; j++) {idx_i[j] = itmp++;}
2892d61bbb3SSatish Balay       }
2902d61bbb3SSatish Balay     }
2912d61bbb3SSatish Balay   }
2922d61bbb3SSatish Balay   PetscFunctionReturn(0);
2932d61bbb3SSatish Balay }
2942d61bbb3SSatish Balay 
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqBAIJ"
29787828ca2SBarry Smith int MatRestoreRow_SeqBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
2982d61bbb3SSatish Balay {
299606d414cSSatish Balay   int ierr;
300606d414cSSatish Balay 
3012d61bbb3SSatish Balay   PetscFunctionBegin;
302606d414cSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
303606d414cSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
3042d61bbb3SSatish Balay   PetscFunctionReturn(0);
3052d61bbb3SSatish Balay }
3062d61bbb3SSatish Balay 
3074a2ae208SSatish Balay #undef __FUNCT__
3084a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqBAIJ"
3092d61bbb3SSatish Balay int MatTranspose_SeqBAIJ(Mat A,Mat *B)
3102d61bbb3SSatish Balay {
3112d61bbb3SSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
3122d61bbb3SSatish Balay   Mat         C;
3132d61bbb3SSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
314273d9f13SBarry Smith   int         *rows,*cols,bs2=a->bs2;
31587828ca2SBarry Smith   PetscScalar *array;
3162d61bbb3SSatish Balay 
3172d61bbb3SSatish Balay   PetscFunctionBegin;
31829bbc08cSBarry Smith   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Square matrix only for in-place");
319b0a32e0cSBarry Smith   ierr = PetscMalloc((1+nbs)*sizeof(int),&col);CHKERRQ(ierr);
320549d3d68SSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
3212d61bbb3SSatish Balay 
322375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
32387828ca2SBarry Smith   ierr = PetscMalloc(a->bs2*a->nz*sizeof(PetscScalar),&array);CHKERRQ(ierr);
32487828ca2SBarry Smith   for (i=0; i<a->bs2*a->nz; i++) array[i] = (PetscScalar)a->a[i];
325375fe846SBarry Smith #else
3263eda8832SBarry Smith   array = a->a;
327375fe846SBarry Smith #endif
328375fe846SBarry Smith 
3292d61bbb3SSatish Balay   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
330273d9f13SBarry Smith   ierr = MatCreateSeqBAIJ(A->comm,bs,A->n,A->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
331606d414cSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
332b0a32e0cSBarry Smith   ierr = PetscMalloc(2*bs*sizeof(int),&rows);CHKERRQ(ierr);
3332d61bbb3SSatish Balay   cols = rows + bs;
3342d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3352d61bbb3SSatish Balay     cols[0] = i*bs;
3362d61bbb3SSatish Balay     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
3372d61bbb3SSatish Balay     len = ai[i+1] - ai[i];
3382d61bbb3SSatish Balay     for (j=0; j<len; j++) {
3392d61bbb3SSatish Balay       rows[0] = (*aj++)*bs;
3402d61bbb3SSatish Balay       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
3412d61bbb3SSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
3422d61bbb3SSatish Balay       array += bs2;
3432d61bbb3SSatish Balay     }
3442d61bbb3SSatish Balay   }
345606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
346375fe846SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
347375fe846SBarry Smith   ierr = PetscFree(array);
348375fe846SBarry Smith #endif
3492d61bbb3SSatish Balay 
3502d61bbb3SSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3512d61bbb3SSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3522d61bbb3SSatish Balay 
353c4992f7dSBarry Smith   if (B) {
3542d61bbb3SSatish Balay     *B = C;
3552d61bbb3SSatish Balay   } else {
356273d9f13SBarry Smith     ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
3572d61bbb3SSatish Balay   }
3582d61bbb3SSatish Balay   PetscFunctionReturn(0);
3592d61bbb3SSatish Balay }
3602d61bbb3SSatish Balay 
3614a2ae208SSatish Balay #undef __FUNCT__
3624a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Binary"
363b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Binary(Mat A,PetscViewer viewer)
3642593348eSBarry Smith {
365b6490206SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
3669df24120SSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
36787828ca2SBarry Smith   PetscScalar *aa;
368ce6f0cecSBarry Smith   FILE        *file;
3692593348eSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
371b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
372b0a32e0cSBarry Smith   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
373552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
3743b2fbd54SBarry Smith 
375273d9f13SBarry Smith   col_lens[1] = A->m;
376273d9f13SBarry Smith   col_lens[2] = A->n;
3777e67e3f9SSatish Balay   col_lens[3] = a->nz*bs2;
3782593348eSBarry Smith 
3792593348eSBarry Smith   /* store lengths of each row and write (including header) to file */
380b6490206SBarry Smith   count = 0;
381b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
382b6490206SBarry Smith     for (j=0; j<bs; j++) {
383b6490206SBarry Smith       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
384b6490206SBarry Smith     }
3852593348eSBarry Smith   }
386273d9f13SBarry Smith   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
387606d414cSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
3882593348eSBarry Smith 
3892593348eSBarry Smith   /* store column indices (zero start index) */
390b0a32e0cSBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
391b6490206SBarry Smith   count = 0;
392b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
393b6490206SBarry Smith     for (j=0; j<bs; j++) {
394b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
395b6490206SBarry Smith         for (l=0; l<bs; l++) {
396b6490206SBarry Smith           jj[count++] = bs*a->j[k] + l;
3972593348eSBarry Smith         }
3982593348eSBarry Smith       }
399b6490206SBarry Smith     }
400b6490206SBarry Smith   }
4010752156aSBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
402606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
4032593348eSBarry Smith 
4042593348eSBarry Smith   /* store nonzero values */
40587828ca2SBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
406b6490206SBarry Smith   count = 0;
407b6490206SBarry Smith   for (i=0; i<a->mbs; i++) {
408b6490206SBarry Smith     for (j=0; j<bs; j++) {
409b6490206SBarry Smith       for (k=a->i[i]; k<a->i[i+1]; k++) {
410b6490206SBarry Smith         for (l=0; l<bs; l++) {
4117e67e3f9SSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
412b6490206SBarry Smith         }
413b6490206SBarry Smith       }
414b6490206SBarry Smith     }
415b6490206SBarry Smith   }
4160752156aSBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
417606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
418ce6f0cecSBarry Smith 
419b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
420ce6f0cecSBarry Smith   if (file) {
421ce6f0cecSBarry Smith     fprintf(file,"-matload_block_size %d\n",a->bs);
422ce6f0cecSBarry Smith   }
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
4242593348eSBarry Smith }
4252593348eSBarry Smith 
426*04929863SHong Zhang extern int MatMPIBAIJFactorInfo_DSCPACK(Mat,PetscViewer);
427*04929863SHong Zhang 
4284a2ae208SSatish Balay #undef __FUNCT__
4294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_ASCII"
430b0a32e0cSBarry Smith static int MatView_SeqBAIJ_ASCII(Mat A,PetscViewer viewer)
4312593348eSBarry Smith {
432b6490206SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
433fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
434f3ef73ceSBarry Smith   PetscViewerFormat format;
4352593348eSBarry Smith 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
437b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
438fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
439b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
440fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
44129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
442*04929863SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
443*04929863SHong Zhang #if defined(PETSC_HAVE_DSCPACK) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
444*04929863SHong Zhang      ierr = MatMPIBAIJFactorInfo_DSCPACK(A,viewer);CHKERRQ(ierr);
445*04929863SHong Zhang #endif
446*04929863SHong Zhang      PetscFunctionReturn(0);
447fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
448b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
44944cd7ae7SLois Curfman McInnes     for (i=0; i<a->mbs; i++) {
45044cd7ae7SLois Curfman McInnes       for (j=0; j<bs; j++) {
451b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
45244cd7ae7SLois Curfman McInnes         for (k=a->i[i]; k<a->i[i+1]; k++) {
45344cd7ae7SLois Curfman McInnes           for (l=0; l<bs; l++) {
454aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4550e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
456b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %gi",bs*a->j[k]+l,
4570e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4580e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
459b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %gi",bs*a->j[k]+l,
4600e6d2581SBarry Smith                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4610e6d2581SBarry Smith             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
462b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4630ef38995SBarry Smith             }
46444cd7ae7SLois Curfman McInnes #else
4650ef38995SBarry Smith             if (a->a[bs2*k + l*bs + j] != 0.0) {
466b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
4670ef38995SBarry Smith             }
46844cd7ae7SLois Curfman McInnes #endif
46944cd7ae7SLois Curfman McInnes           }
47044cd7ae7SLois Curfman McInnes         }
471b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
47244cd7ae7SLois Curfman McInnes       }
47344cd7ae7SLois Curfman McInnes     }
474b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
4750ef38995SBarry Smith   } else {
476b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
477b6490206SBarry Smith     for (i=0; i<a->mbs; i++) {
478b6490206SBarry Smith       for (j=0; j<bs; j++) {
479b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
480b6490206SBarry Smith         for (k=a->i[i]; k<a->i[i+1]; k++) {
481b6490206SBarry Smith           for (l=0; l<bs; l++) {
482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4830e6d2581SBarry Smith             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
484b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
4850e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4860e6d2581SBarry Smith             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
487b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
4880e6d2581SBarry Smith                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
4890ef38995SBarry Smith             } else {
490b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
49188685aaeSLois Curfman McInnes             }
49288685aaeSLois Curfman McInnes #else
493b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
49488685aaeSLois Curfman McInnes #endif
4952593348eSBarry Smith           }
4962593348eSBarry Smith         }
497b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
4982593348eSBarry Smith       }
4992593348eSBarry Smith     }
500b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
501b6490206SBarry Smith   }
502b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
5033a40ed3dSBarry Smith   PetscFunctionReturn(0);
5042593348eSBarry Smith }
5052593348eSBarry Smith 
5064a2ae208SSatish Balay #undef __FUNCT__
5074a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw_Zoom"
508b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
5093270192aSSatish Balay {
51077ed5343SBarry Smith   Mat          A = (Mat) Aa;
5113270192aSSatish Balay   Mat_SeqBAIJ  *a=(Mat_SeqBAIJ*)A->data;
512b65db4caSBarry Smith   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2;
5130e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
5143f1db9ecSBarry Smith   MatScalar    *aa;
515b0a32e0cSBarry Smith   PetscViewer  viewer;
5163270192aSSatish Balay 
5173a40ed3dSBarry Smith   PetscFunctionBegin;
5183270192aSSatish Balay 
519b65db4caSBarry Smith   /* still need to add support for contour plot of nonzeros; see MatView_SeqAIJ_Draw_Zoom()*/
52077ed5343SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
52177ed5343SBarry Smith 
522b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
52377ed5343SBarry Smith 
5243270192aSSatish Balay   /* loop over matrix elements drawing boxes */
525b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
5263270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5273270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
528273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5293270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5303270192aSSatish Balay       aa = a->a + j*bs2;
5313270192aSSatish Balay       for (k=0; k<bs; k++) {
5323270192aSSatish Balay         for (l=0; l<bs; l++) {
5330e6d2581SBarry Smith           if (PetscRealPart(*aa++) >=  0.) continue;
534b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5353270192aSSatish Balay         }
5363270192aSSatish Balay       }
5373270192aSSatish Balay     }
5383270192aSSatish Balay   }
539b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
5403270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5413270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
542273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5433270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5443270192aSSatish Balay       aa = a->a + j*bs2;
5453270192aSSatish Balay       for (k=0; k<bs; k++) {
5463270192aSSatish Balay         for (l=0; l<bs; l++) {
5470e6d2581SBarry Smith           if (PetscRealPart(*aa++) != 0.) continue;
548b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5493270192aSSatish Balay         }
5503270192aSSatish Balay       }
5513270192aSSatish Balay     }
5523270192aSSatish Balay   }
5533270192aSSatish Balay 
554b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
5553270192aSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
5563270192aSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
557273d9f13SBarry Smith       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
5583270192aSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
5593270192aSSatish Balay       aa = a->a + j*bs2;
5603270192aSSatish Balay       for (k=0; k<bs; k++) {
5613270192aSSatish Balay         for (l=0; l<bs; l++) {
5620e6d2581SBarry Smith           if (PetscRealPart(*aa++) <= 0.) continue;
563b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
5643270192aSSatish Balay         }
5653270192aSSatish Balay       }
5663270192aSSatish Balay     }
5673270192aSSatish Balay   }
56877ed5343SBarry Smith   PetscFunctionReturn(0);
56977ed5343SBarry Smith }
5703270192aSSatish Balay 
5714a2ae208SSatish Balay #undef __FUNCT__
5724a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ_Draw"
573b0a32e0cSBarry Smith static int MatView_SeqBAIJ_Draw(Mat A,PetscViewer viewer)
57477ed5343SBarry Smith {
5757dce120fSSatish Balay   int          ierr;
5760e6d2581SBarry Smith   PetscReal    xl,yl,xr,yr,w,h;
577b0a32e0cSBarry Smith   PetscDraw    draw;
57877ed5343SBarry Smith   PetscTruth   isnull;
5793270192aSSatish Balay 
58077ed5343SBarry Smith   PetscFunctionBegin;
58177ed5343SBarry Smith 
582b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
583b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
58477ed5343SBarry Smith 
58577ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
586273d9f13SBarry Smith   xr  = A->n; yr = A->m; h = yr/10.0; w = xr/10.0;
58777ed5343SBarry Smith   xr += w;    yr += h;  xl = -w;     yl = -h;
588b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
589b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
59077ed5343SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
5913a40ed3dSBarry Smith   PetscFunctionReturn(0);
5923270192aSSatish Balay }
5933270192aSSatish Balay 
5944a2ae208SSatish Balay #undef __FUNCT__
5954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqBAIJ"
596b0a32e0cSBarry Smith int MatView_SeqBAIJ(Mat A,PetscViewer viewer)
5972593348eSBarry Smith {
59819bcc07fSBarry Smith   int        ierr;
5996831982aSBarry Smith   PetscTruth issocket,isascii,isbinary,isdraw;
6002593348eSBarry Smith 
6013a40ed3dSBarry Smith   PetscFunctionBegin;
602b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
603b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
604fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
605fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
6060f5bd95cSBarry Smith   if (issocket) {
60729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
6080f5bd95cSBarry Smith   } else if (isascii){
6093a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
6100f5bd95cSBarry Smith   } else if (isbinary) {
6113a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Binary(A,viewer);CHKERRQ(ierr);
6120f5bd95cSBarry Smith   } else if (isdraw) {
6133a40ed3dSBarry Smith     ierr = MatView_SeqBAIJ_Draw(A,viewer);CHKERRQ(ierr);
6145cd90555SBarry Smith   } else {
61529bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
6162593348eSBarry Smith   }
6173a40ed3dSBarry Smith   PetscFunctionReturn(0);
6182593348eSBarry Smith }
619b6490206SBarry Smith 
620cd0e1443SSatish Balay 
6214a2ae208SSatish Balay #undef __FUNCT__
6224a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqBAIJ"
62387828ca2SBarry Smith int MatGetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
624cd0e1443SSatish Balay {
625cd0e1443SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6262d61bbb3SSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
6272d61bbb3SSatish Balay   int        *ai = a->i,*ailen = a->ilen;
6282d61bbb3SSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
6293f1db9ecSBarry Smith   MatScalar  *ap,*aa = a->a,zero = 0.0;
630cd0e1443SSatish Balay 
6313a40ed3dSBarry Smith   PetscFunctionBegin;
6322d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
633cd0e1443SSatish Balay     row  = im[k]; brow = row/bs;
63429bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
635273d9f13SBarry Smith     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
6362d61bbb3SSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
6372c3acbe9SBarry Smith     nrow = ailen[brow];
6382d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
63929bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
640273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
6412d61bbb3SSatish Balay       col  = in[l] ;
6422d61bbb3SSatish Balay       bcol = col/bs;
6432d61bbb3SSatish Balay       cidx = col%bs;
6442d61bbb3SSatish Balay       ridx = row%bs;
6452d61bbb3SSatish Balay       high = nrow;
6462d61bbb3SSatish Balay       low  = 0; /* assume unsorted */
6472d61bbb3SSatish Balay       while (high-low > 5) {
648cd0e1443SSatish Balay         t = (low+high)/2;
649cd0e1443SSatish Balay         if (rp[t] > bcol) high = t;
650cd0e1443SSatish Balay         else             low  = t;
651cd0e1443SSatish Balay       }
652cd0e1443SSatish Balay       for (i=low; i<high; i++) {
653cd0e1443SSatish Balay         if (rp[i] > bcol) break;
654cd0e1443SSatish Balay         if (rp[i] == bcol) {
6552d61bbb3SSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
6562d61bbb3SSatish Balay           goto finished;
657cd0e1443SSatish Balay         }
658cd0e1443SSatish Balay       }
6592d61bbb3SSatish Balay       *v++ = zero;
6602d61bbb3SSatish Balay       finished:;
661cd0e1443SSatish Balay     }
662cd0e1443SSatish Balay   }
6633a40ed3dSBarry Smith   PetscFunctionReturn(0);
664cd0e1443SSatish Balay }
665cd0e1443SSatish Balay 
66695d5f7c2SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
6674a2ae208SSatish Balay #undef __FUNCT__
6684a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
66987828ca2SBarry Smith int MatSetValuesBlocked_SeqBAIJ(Mat mat,int m,int *im,int n,int *in,PetscScalar *v,InsertMode addv)
67095d5f7c2SBarry Smith {
671563b5814SBarry Smith   Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)mat->data;
672563b5814SBarry Smith   int         ierr,i,N = m*n*b->bs2;
673563b5814SBarry Smith   MatScalar   *vsingle;
67495d5f7c2SBarry Smith 
67595d5f7c2SBarry Smith   PetscFunctionBegin;
676563b5814SBarry Smith   if (N > b->setvalueslen) {
677563b5814SBarry Smith     if (b->setvaluescopy) {ierr = PetscFree(b->setvaluescopy);CHKERRQ(ierr);}
678b0a32e0cSBarry Smith     ierr = PetscMalloc(N*sizeof(MatScalar),&b->setvaluescopy);CHKERRQ(ierr);
679563b5814SBarry Smith     b->setvalueslen  = N;
680563b5814SBarry Smith   }
681563b5814SBarry Smith   vsingle = b->setvaluescopy;
68295d5f7c2SBarry Smith   for (i=0; i<N; i++) {
68395d5f7c2SBarry Smith     vsingle[i] = v[i];
68495d5f7c2SBarry Smith   }
68595d5f7c2SBarry Smith   ierr = MatSetValuesBlocked_SeqBAIJ_MatScalar(mat,m,im,n,in,vsingle,addv);CHKERRQ(ierr);
68695d5f7c2SBarry Smith   PetscFunctionReturn(0);
68795d5f7c2SBarry Smith }
68895d5f7c2SBarry Smith #endif
68995d5f7c2SBarry Smith 
6902d61bbb3SSatish Balay 
6914a2ae208SSatish Balay #undef __FUNCT__
6924a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqBAIJ"
69395d5f7c2SBarry Smith int MatSetValuesBlocked_SeqBAIJ_MatScalar(Mat A,int m,int *im,int n,int *in,MatScalar *v,InsertMode is)
69492c4ed94SBarry Smith {
69592c4ed94SBarry Smith   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6968a84c255SSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
697273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
698549d3d68SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
699273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
70095d5f7c2SBarry Smith   MatScalar   *value = v,*ap,*aa = a->a,*bap;
70192c4ed94SBarry Smith 
7023a40ed3dSBarry Smith   PetscFunctionBegin;
7030e324ae4SSatish Balay   if (roworiented) {
7040e324ae4SSatish Balay     stepval = (n-1)*bs;
7050e324ae4SSatish Balay   } else {
7060e324ae4SSatish Balay     stepval = (m-1)*bs;
7070e324ae4SSatish Balay   }
70892c4ed94SBarry Smith   for (k=0; k<m; k++) { /* loop over added rows */
70992c4ed94SBarry Smith     row  = im[k];
7105ef9f2a5SBarry Smith     if (row < 0) continue;
711aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
71229bbc08cSBarry Smith     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
71392c4ed94SBarry Smith #endif
71492c4ed94SBarry Smith     rp   = aj + ai[row];
71592c4ed94SBarry Smith     ap   = aa + bs2*ai[row];
71692c4ed94SBarry Smith     rmax = imax[row];
71792c4ed94SBarry Smith     nrow = ailen[row];
71892c4ed94SBarry Smith     low  = 0;
71992c4ed94SBarry Smith     for (l=0; l<n; l++) { /* loop over added columns */
7205ef9f2a5SBarry Smith       if (in[l] < 0) continue;
721aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
72229bbc08cSBarry Smith       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
72392c4ed94SBarry Smith #endif
72492c4ed94SBarry Smith       col = in[l];
72592c4ed94SBarry Smith       if (roworiented) {
7260e324ae4SSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
7270e324ae4SSatish Balay       } else {
7280e324ae4SSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
72992c4ed94SBarry Smith       }
73092c4ed94SBarry Smith       if (!sorted) low = 0; high = nrow;
73192c4ed94SBarry Smith       while (high-low > 7) {
73292c4ed94SBarry Smith         t = (low+high)/2;
73392c4ed94SBarry Smith         if (rp[t] > col) high = t;
73492c4ed94SBarry Smith         else             low  = t;
73592c4ed94SBarry Smith       }
73692c4ed94SBarry Smith       for (i=low; i<high; i++) {
73792c4ed94SBarry Smith         if (rp[i] > col) break;
73892c4ed94SBarry Smith         if (rp[i] == col) {
7398a84c255SSatish Balay           bap  = ap +  bs2*i;
7400e324ae4SSatish Balay           if (roworiented) {
7418a84c255SSatish Balay             if (is == ADD_VALUES) {
742dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
743dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
7448a84c255SSatish Balay                   bap[jj] += *value++;
745dd9472c6SBarry Smith                 }
746dd9472c6SBarry Smith               }
7470e324ae4SSatish Balay             } else {
748dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
749dd9472c6SBarry Smith                 for (jj=ii; jj<bs2; jj+=bs) {
7500e324ae4SSatish Balay                   bap[jj] = *value++;
7518a84c255SSatish Balay                 }
752dd9472c6SBarry Smith               }
753dd9472c6SBarry Smith             }
7540e324ae4SSatish Balay           } else {
7550e324ae4SSatish Balay             if (is == ADD_VALUES) {
756dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
757dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
7580e324ae4SSatish Balay                   *bap++ += *value++;
759dd9472c6SBarry Smith                 }
760dd9472c6SBarry Smith               }
7610e324ae4SSatish Balay             } else {
762dd9472c6SBarry Smith               for (ii=0; ii<bs; ii++,value+=stepval) {
763dd9472c6SBarry Smith                 for (jj=0; jj<bs; jj++) {
7640e324ae4SSatish Balay                   *bap++  = *value++;
7650e324ae4SSatish Balay                 }
7668a84c255SSatish Balay               }
767dd9472c6SBarry Smith             }
768dd9472c6SBarry Smith           }
769f1241b54SBarry Smith           goto noinsert2;
77092c4ed94SBarry Smith         }
77192c4ed94SBarry Smith       }
77289280ab3SLois Curfman McInnes       if (nonew == 1) goto noinsert2;
77329bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
77492c4ed94SBarry Smith       if (nrow >= rmax) {
77592c4ed94SBarry Smith         /* there is no extra room in row, therefore enlarge */
77692c4ed94SBarry Smith         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
7773f1db9ecSBarry Smith         MatScalar *new_a;
77892c4ed94SBarry Smith 
77929bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
78089280ab3SLois Curfman McInnes 
78192c4ed94SBarry Smith         /* malloc new storage space */
7823f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
783b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
78492c4ed94SBarry Smith         new_j   = (int*)(new_a + bs2*new_nz);
78592c4ed94SBarry Smith         new_i   = new_j + new_nz;
78692c4ed94SBarry Smith 
78792c4ed94SBarry Smith         /* copy over old data into new slots */
78892c4ed94SBarry Smith         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
78992c4ed94SBarry Smith         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
790549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
79192c4ed94SBarry Smith         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
792549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
793549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
794549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
795549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
79692c4ed94SBarry Smith         /* free up old matrix storage */
797606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
798606d414cSSatish Balay         if (!a->singlemalloc) {
799606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
800606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
801606d414cSSatish Balay         }
80292c4ed94SBarry Smith         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
803c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
80492c4ed94SBarry Smith 
80592c4ed94SBarry Smith         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
80692c4ed94SBarry Smith         rmax = imax[row] = imax[row] + CHUNKSIZE;
807b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
80892c4ed94SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
80992c4ed94SBarry Smith         a->reallocs++;
81092c4ed94SBarry Smith         a->nz++;
81192c4ed94SBarry Smith       }
81292c4ed94SBarry Smith       N = nrow++ - 1;
81392c4ed94SBarry Smith       /* shift up all the later entries in this row */
81492c4ed94SBarry Smith       for (ii=N; ii>=i; ii--) {
81592c4ed94SBarry Smith         rp[ii+1] = rp[ii];
816549d3d68SSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
81792c4ed94SBarry Smith       }
818549d3d68SSatish Balay       if (N >= i) {
819549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
820549d3d68SSatish Balay       }
82192c4ed94SBarry Smith       rp[i] = col;
8228a84c255SSatish Balay       bap   = ap +  bs2*i;
8230e324ae4SSatish Balay       if (roworiented) {
824dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
825dd9472c6SBarry Smith           for (jj=ii; jj<bs2; jj+=bs) {
8260e324ae4SSatish Balay             bap[jj] = *value++;
827dd9472c6SBarry Smith           }
828dd9472c6SBarry Smith         }
8290e324ae4SSatish Balay       } else {
830dd9472c6SBarry Smith         for (ii=0; ii<bs; ii++,value+=stepval) {
831dd9472c6SBarry Smith           for (jj=0; jj<bs; jj++) {
8320e324ae4SSatish Balay             *bap++  = *value++;
8330e324ae4SSatish Balay           }
834dd9472c6SBarry Smith         }
835dd9472c6SBarry Smith       }
836f1241b54SBarry Smith       noinsert2:;
83792c4ed94SBarry Smith       low = i;
83892c4ed94SBarry Smith     }
83992c4ed94SBarry Smith     ailen[row] = nrow;
84092c4ed94SBarry Smith   }
8413a40ed3dSBarry Smith   PetscFunctionReturn(0);
84292c4ed94SBarry Smith }
84392c4ed94SBarry Smith 
844*04929863SHong Zhang /* EXTERN int MatUseDSCPACK_MPIBAIJ(Mat); */
8454a2ae208SSatish Balay #undef __FUNCT__
8464a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqBAIJ"
8478f6be9afSLois Curfman McInnes int MatAssemblyEnd_SeqBAIJ(Mat A,MatAssemblyType mode)
848584200bdSSatish Balay {
849584200bdSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
850584200bdSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
851273d9f13SBarry Smith   int        m = A->m,*ip,N,*ailen = a->ilen;
852549d3d68SSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
8533f1db9ecSBarry Smith   MatScalar  *aa = a->a,*ap;
854a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
855a56a16c8SHong Zhang   PetscTruth   flag;
856a56a16c8SHong Zhang #endif
857584200bdSSatish Balay 
8583a40ed3dSBarry Smith   PetscFunctionBegin;
8593a40ed3dSBarry Smith   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
860584200bdSSatish Balay 
86143ee02c3SBarry Smith   if (m) rmax = ailen[0];
862584200bdSSatish Balay   for (i=1; i<mbs; i++) {
863584200bdSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
864584200bdSSatish Balay     fshift += imax[i-1] - ailen[i-1];
865d402145bSBarry Smith     rmax   = PetscMax(rmax,ailen[i]);
866584200bdSSatish Balay     if (fshift) {
867a7c10996SSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
868584200bdSSatish Balay       N = ailen[i];
869584200bdSSatish Balay       for (j=0; j<N; j++) {
870584200bdSSatish Balay         ip[j-fshift] = ip[j];
871549d3d68SSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
872584200bdSSatish Balay       }
873584200bdSSatish Balay     }
874584200bdSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
875584200bdSSatish Balay   }
876584200bdSSatish Balay   if (mbs) {
877584200bdSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
878584200bdSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
879584200bdSSatish Balay   }
880584200bdSSatish Balay   /* reset ilen and imax for each row */
881584200bdSSatish Balay   for (i=0; i<mbs; i++) {
882584200bdSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
883584200bdSSatish Balay   }
884a7c10996SSatish Balay   a->nz = ai[mbs];
885584200bdSSatish Balay 
886584200bdSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
887584200bdSSatish Balay   if (fshift && a->diag) {
888606d414cSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
889b424e231SHong Zhang     PetscLogObjectMemory(A,-(mbs+1)*sizeof(int));
890584200bdSSatish Balay     a->diag = 0;
891584200bdSSatish Balay   }
892bba1ac68SSatish 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);
893bba1ac68SSatish Balay   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Number of mallocs during MatSetValues is %d\n",a->reallocs);
894b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
895e2f3b5e9SSatish Balay   a->reallocs          = 0;
8960e6d2581SBarry Smith   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
897a56a16c8SHong Zhang #if defined(PETSC_HAVE_DSCPACK)
898a56a16c8SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_baij_dscpack",&flag);CHKERRQ(ierr);
899a56a16c8SHong Zhang   if (flag) { ierr = MatUseDSCPACK_MPIBAIJ(A);CHKERRQ(ierr); }
900a56a16c8SHong Zhang #endif
901a56a16c8SHong Zhang 
9024e220ebcSLois Curfman McInnes 
9033a40ed3dSBarry Smith   PetscFunctionReturn(0);
904584200bdSSatish Balay }
905584200bdSSatish Balay 
9065a838052SSatish Balay 
907bea157c4SSatish Balay 
908bea157c4SSatish Balay /*
909bea157c4SSatish Balay    This function returns an array of flags which indicate the locations of contiguous
910bea157c4SSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
911bea157c4SSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
912bea157c4SSatish Balay    Assume: sizes should be long enough to hold all the values.
913bea157c4SSatish Balay */
9144a2ae208SSatish Balay #undef __FUNCT__
9154a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ_Check_Blocks"
916bea157c4SSatish Balay static int MatZeroRows_SeqBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
917d9b7c43dSSatish Balay {
918bea157c4SSatish Balay   int        i,j,k,row;
919bea157c4SSatish Balay   PetscTruth flg;
9203a40ed3dSBarry Smith 
921433994e6SBarry Smith   PetscFunctionBegin;
922bea157c4SSatish Balay   for (i=0,j=0; i<n; j++) {
923bea157c4SSatish Balay     row = idx[i];
924bea157c4SSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
925bea157c4SSatish Balay       sizes[j] = 1;
926bea157c4SSatish Balay       i++;
927e4fda26cSSatish Balay     } else if (i+bs > n) { /* complete block doesn't exist (at idx end) */
928bea157c4SSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
929bea157c4SSatish Balay       i++;
930bea157c4SSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
931bea157c4SSatish Balay       flg = PETSC_TRUE;
932bea157c4SSatish Balay       for (k=1; k<bs; k++) {
933bea157c4SSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
934bea157c4SSatish Balay           flg = PETSC_FALSE;
935bea157c4SSatish Balay           break;
936d9b7c43dSSatish Balay         }
937bea157c4SSatish Balay       }
938bea157c4SSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
939bea157c4SSatish Balay         sizes[j] = bs;
940bea157c4SSatish Balay         i+= bs;
941bea157c4SSatish Balay       } else {
942bea157c4SSatish Balay         sizes[j] = 1;
943bea157c4SSatish Balay         i++;
944bea157c4SSatish Balay       }
945bea157c4SSatish Balay     }
946bea157c4SSatish Balay   }
947bea157c4SSatish Balay   *bs_max = j;
9483a40ed3dSBarry Smith   PetscFunctionReturn(0);
949d9b7c43dSSatish Balay }
950d9b7c43dSSatish Balay 
9514a2ae208SSatish Balay #undef __FUNCT__
9524a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqBAIJ"
95387828ca2SBarry Smith int MatZeroRows_SeqBAIJ(Mat A,IS is,PetscScalar *diag)
954d9b7c43dSSatish Balay {
955d9b7c43dSSatish Balay   Mat_SeqBAIJ *baij=(Mat_SeqBAIJ*)A->data;
956b6815fffSBarry Smith   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
957bea157c4SSatish Balay   int         bs=baij->bs,bs2=baij->bs2,*sizes,row,bs_max;
95887828ca2SBarry Smith   PetscScalar zero = 0.0;
9593f1db9ecSBarry Smith   MatScalar   *aa;
960d9b7c43dSSatish Balay 
9613a40ed3dSBarry Smith   PetscFunctionBegin;
962d9b7c43dSSatish Balay   /* Make a copy of the IS and  sort it */
963b9b97703SBarry Smith   ierr = ISGetLocalSize(is,&is_n);CHKERRQ(ierr);
964d9b7c43dSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
965d9b7c43dSSatish Balay 
966bea157c4SSatish Balay   /* allocate memory for rows,sizes */
967b0a32e0cSBarry Smith   ierr  = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
968bea157c4SSatish Balay   sizes = rows + is_n;
969bea157c4SSatish Balay 
970563b5814SBarry Smith   /* copy IS values to rows, and sort them */
971bea157c4SSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
972bea157c4SSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
973dffd3267SBarry Smith   if (baij->keepzeroedrows) {
974dffd3267SBarry Smith     for (i=0; i<is_n; i++) { sizes[i] = 1; }
975dffd3267SBarry Smith     bs_max = is_n;
976dffd3267SBarry Smith   } else {
977bea157c4SSatish Balay     ierr = MatZeroRows_SeqBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
978dffd3267SBarry Smith   }
979b97a56abSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
980bea157c4SSatish Balay 
981bea157c4SSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
982bea157c4SSatish Balay     row   = rows[j];
983273d9f13SBarry Smith     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
984bea157c4SSatish Balay     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
985bea157c4SSatish Balay     aa    = baij->a + baij->i[row/bs]*bs2 + (row%bs);
986dffd3267SBarry Smith     if (sizes[i] == bs && !baij->keepzeroedrows) {
987bea157c4SSatish Balay       if (diag) {
988bea157c4SSatish Balay         if (baij->ilen[row/bs] > 0) {
989bea157c4SSatish Balay           baij->ilen[row/bs]       = 1;
990bea157c4SSatish Balay           baij->j[baij->i[row/bs]] = row/bs;
991bea157c4SSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
992a07cd24cSSatish Balay         }
993563b5814SBarry Smith         /* Now insert all the diagonal values for this bs */
994bea157c4SSatish Balay         for (k=0; k<bs; k++) {
995bea157c4SSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
996bea157c4SSatish Balay         }
997bea157c4SSatish Balay       } else { /* (!diag) */
998bea157c4SSatish Balay         baij->ilen[row/bs] = 0;
999bea157c4SSatish Balay       } /* end (!diag) */
1000bea157c4SSatish Balay     } else { /* (sizes[i] != bs) */
1001aa482453SBarry Smith #if defined (PETSC_USE_DEBUG)
100229bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
1003bea157c4SSatish Balay #endif
1004bea157c4SSatish Balay       for (k=0; k<count; k++) {
1005d9b7c43dSSatish Balay         aa[0] =  zero;
1006d9b7c43dSSatish Balay         aa    += bs;
1007d9b7c43dSSatish Balay       }
1008d9b7c43dSSatish Balay       if (diag) {
1009f830108cSBarry Smith         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
1010d9b7c43dSSatish Balay       }
1011d9b7c43dSSatish Balay     }
1012bea157c4SSatish Balay   }
1013bea157c4SSatish Balay 
1014606d414cSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
10159a8dea36SBarry Smith   ierr = MatAssemblyEnd_SeqBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
10163a40ed3dSBarry Smith   PetscFunctionReturn(0);
1017d9b7c43dSSatish Balay }
10181c351548SSatish Balay 
10194a2ae208SSatish Balay #undef __FUNCT__
10204a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqBAIJ"
102187828ca2SBarry Smith int MatSetValues_SeqBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
10222d61bbb3SSatish Balay {
10232d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
10242d61bbb3SSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
1025273d9f13SBarry Smith   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
10262d61bbb3SSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
1027549d3d68SSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
1028273d9f13SBarry Smith   PetscTruth  roworiented=a->roworiented;
10293f1db9ecSBarry Smith   MatScalar   *ap,value,*aa=a->a,*bap;
10302d61bbb3SSatish Balay 
10312d61bbb3SSatish Balay   PetscFunctionBegin;
10322d61bbb3SSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
10332d61bbb3SSatish Balay     row  = im[k]; brow = row/bs;
10345ef9f2a5SBarry Smith     if (row < 0) continue;
1035aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1036273d9f13SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
10372d61bbb3SSatish Balay #endif
10382d61bbb3SSatish Balay     rp   = aj + ai[brow];
10392d61bbb3SSatish Balay     ap   = aa + bs2*ai[brow];
10402d61bbb3SSatish Balay     rmax = imax[brow];
10412d61bbb3SSatish Balay     nrow = ailen[brow];
10422d61bbb3SSatish Balay     low  = 0;
10432d61bbb3SSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
10445ef9f2a5SBarry Smith       if (in[l] < 0) continue;
1045aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g)
1046273d9f13SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n);
10472d61bbb3SSatish Balay #endif
10482d61bbb3SSatish Balay       col = in[l]; bcol = col/bs;
10492d61bbb3SSatish Balay       ridx = row % bs; cidx = col % bs;
10502d61bbb3SSatish Balay       if (roworiented) {
10515ef9f2a5SBarry Smith         value = v[l + k*n];
10522d61bbb3SSatish Balay       } else {
10532d61bbb3SSatish Balay         value = v[k + l*m];
10542d61bbb3SSatish Balay       }
10552d61bbb3SSatish Balay       if (!sorted) low = 0; high = nrow;
10562d61bbb3SSatish Balay       while (high-low > 7) {
10572d61bbb3SSatish Balay         t = (low+high)/2;
10582d61bbb3SSatish Balay         if (rp[t] > bcol) high = t;
10592d61bbb3SSatish Balay         else              low  = t;
10602d61bbb3SSatish Balay       }
10612d61bbb3SSatish Balay       for (i=low; i<high; i++) {
10622d61bbb3SSatish Balay         if (rp[i] > bcol) break;
10632d61bbb3SSatish Balay         if (rp[i] == bcol) {
10642d61bbb3SSatish Balay           bap  = ap +  bs2*i + bs*cidx + ridx;
10652d61bbb3SSatish Balay           if (is == ADD_VALUES) *bap += value;
10662d61bbb3SSatish Balay           else                  *bap  = value;
10672d61bbb3SSatish Balay           goto noinsert1;
10682d61bbb3SSatish Balay         }
10692d61bbb3SSatish Balay       }
10702d61bbb3SSatish Balay       if (nonew == 1) goto noinsert1;
107129bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
10722d61bbb3SSatish Balay       if (nrow >= rmax) {
10732d61bbb3SSatish Balay         /* there is no extra room in row, therefore enlarge */
10742d61bbb3SSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
10753f1db9ecSBarry Smith         MatScalar *new_a;
10762d61bbb3SSatish Balay 
107729bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
10782d61bbb3SSatish Balay 
10792d61bbb3SSatish Balay         /* Malloc new storage space */
10803f1db9ecSBarry Smith         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
1081b0a32e0cSBarry Smith 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
10822d61bbb3SSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
10832d61bbb3SSatish Balay         new_i   = new_j + new_nz;
10842d61bbb3SSatish Balay 
10852d61bbb3SSatish Balay         /* copy over old data into new slots */
10862d61bbb3SSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
10872d61bbb3SSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
1088549d3d68SSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
10892d61bbb3SSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
1090549d3d68SSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
1091549d3d68SSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1092549d3d68SSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
1093549d3d68SSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
10942d61bbb3SSatish Balay         /* free up old matrix storage */
1095606d414cSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
1096606d414cSSatish Balay         if (!a->singlemalloc) {
1097606d414cSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
1098606d414cSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
1099606d414cSSatish Balay         }
11002d61bbb3SSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
1101c4992f7dSBarry Smith         a->singlemalloc = PETSC_TRUE;
11022d61bbb3SSatish Balay 
11032d61bbb3SSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
11042d61bbb3SSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1105b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
11062d61bbb3SSatish Balay         a->maxnz += bs2*CHUNKSIZE;
11072d61bbb3SSatish Balay         a->reallocs++;
11082d61bbb3SSatish Balay         a->nz++;
11092d61bbb3SSatish Balay       }
11102d61bbb3SSatish Balay       N = nrow++ - 1;
11112d61bbb3SSatish Balay       /* shift up all the later entries in this row */
11122d61bbb3SSatish Balay       for (ii=N; ii>=i; ii--) {
11132d61bbb3SSatish Balay         rp[ii+1] = rp[ii];
1114549d3d68SSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
11152d61bbb3SSatish Balay       }
1116549d3d68SSatish Balay       if (N>=i) {
1117549d3d68SSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1118549d3d68SSatish Balay       }
11192d61bbb3SSatish Balay       rp[i]                      = bcol;
11202d61bbb3SSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
11212d61bbb3SSatish Balay       noinsert1:;
11222d61bbb3SSatish Balay       low = i;
11232d61bbb3SSatish Balay     }
11242d61bbb3SSatish Balay     ailen[brow] = nrow;
11252d61bbb3SSatish Balay   }
11262d61bbb3SSatish Balay   PetscFunctionReturn(0);
11272d61bbb3SSatish Balay }
11282d61bbb3SSatish Balay 
11292d61bbb3SSatish Balay 
11304a2ae208SSatish Balay #undef __FUNCT__
11314a2ae208SSatish Balay #define __FUNCT__ "MatILUFactor_SeqBAIJ"
11325ef9f2a5SBarry Smith int MatILUFactor_SeqBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
11332d61bbb3SSatish Balay {
11342d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
11352d61bbb3SSatish Balay   Mat         outA;
11362d61bbb3SSatish Balay   int         ierr;
1137667159a5SBarry Smith   PetscTruth  row_identity,col_identity;
11382d61bbb3SSatish Balay 
11392d61bbb3SSatish Balay   PetscFunctionBegin;
114029bbc08cSBarry Smith   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
1141667159a5SBarry Smith   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1142667159a5SBarry Smith   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
1143667159a5SBarry Smith   if (!row_identity || !col_identity) {
114429bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ILU");
1145667159a5SBarry Smith   }
11462d61bbb3SSatish Balay 
11472d61bbb3SSatish Balay   outA          = inA;
11482d61bbb3SSatish Balay   inA->factor   = FACTOR_LU;
11492d61bbb3SSatish Balay 
11502d61bbb3SSatish Balay   if (!a->diag) {
1151c4992f7dSBarry Smith     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
11522d61bbb3SSatish Balay   }
1153cf242676SKris Buschelman 
1154c38d4ed2SBarry Smith   a->row        = row;
1155c38d4ed2SBarry Smith   a->col        = col;
1156c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1157c38d4ed2SBarry Smith   ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
1158c38d4ed2SBarry Smith 
1159c38d4ed2SBarry Smith   /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
11604c49b128SBarry Smith   ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr);
1161b0a32e0cSBarry Smith   PetscLogObjectParent(inA,a->icol);
1162c38d4ed2SBarry Smith 
1163cf242676SKris Buschelman   /*
1164cf242676SKris Buschelman       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
1165cf242676SKris Buschelman       for ILU(0) factorization with natural ordering
1166cf242676SKris Buschelman   */
1167cf242676SKris Buschelman   if (a->bs < 8) {
1168cf242676SKris Buschelman     ierr = MatSeqBAIJ_UpdateFactorNumeric_NaturalOrdering(inA);CHKERRQ(ierr);
1169cf242676SKris Buschelman   } else {
1170c38d4ed2SBarry Smith     if (!a->solve_work) {
117187828ca2SBarry Smith       ierr = PetscMalloc((inA->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
117287828ca2SBarry Smith       PetscLogObjectMemory(inA,(inA->m+a->bs)*sizeof(PetscScalar));
1173c38d4ed2SBarry Smith     }
11742d61bbb3SSatish Balay   }
11752d61bbb3SSatish Balay 
1176667159a5SBarry Smith   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
1177667159a5SBarry Smith 
11782d61bbb3SSatish Balay   PetscFunctionReturn(0);
11792d61bbb3SSatish Balay }
11804a2ae208SSatish Balay #undef __FUNCT__
11814a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqBAIJ"
1182ba4ca20aSSatish Balay int MatPrintHelp_SeqBAIJ(Mat A)
1183ba4ca20aSSatish Balay {
1184c38d4ed2SBarry Smith   static PetscTruth called = PETSC_FALSE;
1185ba4ca20aSSatish Balay   MPI_Comm          comm = A->comm;
1186d132466eSBarry Smith   int               ierr;
1187ba4ca20aSSatish Balay 
11883a40ed3dSBarry Smith   PetscFunctionBegin;
1189c38d4ed2SBarry Smith   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
1190d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQBAIJ and MATMPIBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
1191d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
11923a40ed3dSBarry Smith   PetscFunctionReturn(0);
1193ba4ca20aSSatish Balay }
1194d9b7c43dSSatish Balay 
1195fb2e594dSBarry Smith EXTERN_C_BEGIN
11964a2ae208SSatish Balay #undef __FUNCT__
11974a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices_SeqBAIJ"
119827a8da17SBarry Smith int MatSeqBAIJSetColumnIndices_SeqBAIJ(Mat mat,int *indices)
119927a8da17SBarry Smith {
120027a8da17SBarry Smith   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
120114db4f74SSatish Balay   int         i,nz,nbs;
120227a8da17SBarry Smith 
120327a8da17SBarry Smith   PetscFunctionBegin;
120414db4f74SSatish Balay   nz  = baij->maxnz/baij->bs2;
120514db4f74SSatish Balay   nbs = baij->nbs;
120627a8da17SBarry Smith   for (i=0; i<nz; i++) {
120727a8da17SBarry Smith     baij->j[i] = indices[i];
120827a8da17SBarry Smith   }
120927a8da17SBarry Smith   baij->nz = nz;
121014db4f74SSatish Balay   for (i=0; i<nbs; i++) {
121127a8da17SBarry Smith     baij->ilen[i] = baij->imax[i];
121227a8da17SBarry Smith   }
121327a8da17SBarry Smith 
121427a8da17SBarry Smith   PetscFunctionReturn(0);
121527a8da17SBarry Smith }
1216fb2e594dSBarry Smith EXTERN_C_END
121727a8da17SBarry Smith 
12184a2ae208SSatish Balay #undef __FUNCT__
12194a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetColumnIndices"
122027a8da17SBarry Smith /*@
122127a8da17SBarry Smith     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
122227a8da17SBarry Smith        in the matrix.
122327a8da17SBarry Smith 
122427a8da17SBarry Smith   Input Parameters:
122527a8da17SBarry Smith +  mat - the SeqBAIJ matrix
122627a8da17SBarry Smith -  indices - the column indices
122727a8da17SBarry Smith 
122815091d37SBarry Smith   Level: advanced
122915091d37SBarry Smith 
123027a8da17SBarry Smith   Notes:
123127a8da17SBarry Smith     This can be called if you have precomputed the nonzero structure of the
123227a8da17SBarry Smith   matrix and want to provide it to the matrix object to improve the performance
123327a8da17SBarry Smith   of the MatSetValues() operation.
123427a8da17SBarry Smith 
123527a8da17SBarry Smith     You MUST have set the correct numbers of nonzeros per row in the call to
123627a8da17SBarry Smith   MatCreateSeqBAIJ().
123727a8da17SBarry Smith 
123827a8da17SBarry Smith     MUST be called before any calls to MatSetValues();
123927a8da17SBarry Smith 
124027a8da17SBarry Smith @*/
124127a8da17SBarry Smith int MatSeqBAIJSetColumnIndices(Mat mat,int *indices)
124227a8da17SBarry Smith {
124327a8da17SBarry Smith   int ierr,(*f)(Mat,int *);
124427a8da17SBarry Smith 
124527a8da17SBarry Smith   PetscFunctionBegin;
124627a8da17SBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1247c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
124827a8da17SBarry Smith   if (f) {
124927a8da17SBarry Smith     ierr = (*f)(mat,indices);CHKERRQ(ierr);
125027a8da17SBarry Smith   } else {
125129bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
125227a8da17SBarry Smith   }
125327a8da17SBarry Smith   PetscFunctionReturn(0);
125427a8da17SBarry Smith }
125527a8da17SBarry Smith 
12564a2ae208SSatish Balay #undef __FUNCT__
12574a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqBAIJ"
1258273d9f13SBarry Smith int MatGetRowMax_SeqBAIJ(Mat A,Vec v)
1259273d9f13SBarry Smith {
1260273d9f13SBarry Smith   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1261273d9f13SBarry Smith   int          ierr,i,j,n,row,bs,*ai,*aj,mbs;
1262273d9f13SBarry Smith   PetscReal    atmp;
126387828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
1264273d9f13SBarry Smith   MatScalar    *aa;
1265273d9f13SBarry Smith   int          ncols,brow,krow,kcol;
1266273d9f13SBarry Smith 
1267273d9f13SBarry Smith   PetscFunctionBegin;
1268273d9f13SBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1269273d9f13SBarry Smith   bs   = a->bs;
1270273d9f13SBarry Smith   aa   = a->a;
1271273d9f13SBarry Smith   ai   = a->i;
1272273d9f13SBarry Smith   aj   = a->j;
1273273d9f13SBarry Smith   mbs = a->mbs;
1274273d9f13SBarry Smith 
1275273d9f13SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1276273d9f13SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1277273d9f13SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1278273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1279273d9f13SBarry Smith   for (i=0; i<mbs; i++) {
1280273d9f13SBarry Smith     ncols = ai[1] - ai[0]; ai++;
1281273d9f13SBarry Smith     brow  = bs*i;
1282273d9f13SBarry Smith     for (j=0; j<ncols; j++){
1283273d9f13SBarry Smith       /* bcol = bs*(*aj); */
1284273d9f13SBarry Smith       for (kcol=0; kcol<bs; kcol++){
1285273d9f13SBarry Smith         for (krow=0; krow<bs; krow++){
1286273d9f13SBarry Smith           atmp = PetscAbsScalar(*aa); aa++;
1287273d9f13SBarry Smith           row = brow + krow;    /* row index */
1288273d9f13SBarry Smith           /* printf("val[%d,%d]: %g\n",row,bcol+kcol,atmp); */
1289273d9f13SBarry Smith           if (PetscAbsScalar(x[row]) < atmp) x[row] = atmp;
1290273d9f13SBarry Smith         }
1291273d9f13SBarry Smith       }
1292273d9f13SBarry Smith       aj++;
1293273d9f13SBarry Smith     }
1294273d9f13SBarry Smith   }
1295273d9f13SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1296273d9f13SBarry Smith   PetscFunctionReturn(0);
1297273d9f13SBarry Smith }
1298273d9f13SBarry Smith 
12994a2ae208SSatish Balay #undef __FUNCT__
13004a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqBAIJ"
1301273d9f13SBarry Smith int MatSetUpPreallocation_SeqBAIJ(Mat A)
1302273d9f13SBarry Smith {
1303273d9f13SBarry Smith   int        ierr;
1304273d9f13SBarry Smith 
1305273d9f13SBarry Smith   PetscFunctionBegin;
1306273d9f13SBarry Smith   ierr =  MatSeqBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1307273d9f13SBarry Smith   PetscFunctionReturn(0);
1308273d9f13SBarry Smith }
1309273d9f13SBarry Smith 
13104a2ae208SSatish Balay #undef __FUNCT__
13114a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqBAIJ"
131287828ca2SBarry Smith int MatGetArray_SeqBAIJ(Mat A,PetscScalar **array)
1313f2a5309cSSatish Balay {
1314f2a5309cSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1315f2a5309cSSatish Balay   PetscFunctionBegin;
1316f2a5309cSSatish Balay   *array = a->a;
1317f2a5309cSSatish Balay   PetscFunctionReturn(0);
1318f2a5309cSSatish Balay }
1319f2a5309cSSatish Balay 
13204a2ae208SSatish Balay #undef __FUNCT__
13214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqBAIJ"
132287828ca2SBarry Smith int MatRestoreArray_SeqBAIJ(Mat A,PetscScalar **array)
1323f2a5309cSSatish Balay {
1324f2a5309cSSatish Balay   PetscFunctionBegin;
1325f2a5309cSSatish Balay   PetscFunctionReturn(0);
1326f2a5309cSSatish Balay }
1327f2a5309cSSatish Balay 
13282593348eSBarry Smith /* -------------------------------------------------------------------*/
1329cc2dc46cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqBAIJ,
1330cc2dc46cSBarry Smith        MatGetRow_SeqBAIJ,
1331cc2dc46cSBarry Smith        MatRestoreRow_SeqBAIJ,
1332cc2dc46cSBarry Smith        MatMult_SeqBAIJ_N,
1333cc2dc46cSBarry Smith        MatMultAdd_SeqBAIJ_N,
13347c922b88SBarry Smith        MatMultTranspose_SeqBAIJ,
13357c922b88SBarry Smith        MatMultTransposeAdd_SeqBAIJ,
1336cc2dc46cSBarry Smith        MatSolve_SeqBAIJ_N,
1337cc2dc46cSBarry Smith        0,
1338cc2dc46cSBarry Smith        0,
1339cc2dc46cSBarry Smith        0,
1340cc2dc46cSBarry Smith        MatLUFactor_SeqBAIJ,
1341cc2dc46cSBarry Smith        0,
1342b6490206SBarry Smith        0,
1343f2501298SSatish Balay        MatTranspose_SeqBAIJ,
1344cc2dc46cSBarry Smith        MatGetInfo_SeqBAIJ,
1345cc2dc46cSBarry Smith        MatEqual_SeqBAIJ,
1346cc2dc46cSBarry Smith        MatGetDiagonal_SeqBAIJ,
1347cc2dc46cSBarry Smith        MatDiagonalScale_SeqBAIJ,
1348cc2dc46cSBarry Smith        MatNorm_SeqBAIJ,
1349b6490206SBarry Smith        0,
1350cc2dc46cSBarry Smith        MatAssemblyEnd_SeqBAIJ,
1351cc2dc46cSBarry Smith        0,
1352cc2dc46cSBarry Smith        MatSetOption_SeqBAIJ,
1353cc2dc46cSBarry Smith        MatZeroEntries_SeqBAIJ,
1354cc2dc46cSBarry Smith        MatZeroRows_SeqBAIJ,
1355cc2dc46cSBarry Smith        MatLUFactorSymbolic_SeqBAIJ,
1356cc2dc46cSBarry Smith        MatLUFactorNumeric_SeqBAIJ_N,
1357cc2dc46cSBarry Smith        0,
1358cc2dc46cSBarry Smith        0,
1359273d9f13SBarry Smith        MatSetUpPreallocation_SeqBAIJ,
1360cc2dc46cSBarry Smith        MatILUFactorSymbolic_SeqBAIJ,
1361cc2dc46cSBarry Smith        0,
1362f2a5309cSSatish Balay        MatGetArray_SeqBAIJ,
1363f2a5309cSSatish Balay        MatRestoreArray_SeqBAIJ,
13642e8a6d31SBarry Smith        MatDuplicate_SeqBAIJ,
1365cc2dc46cSBarry Smith        0,
1366cc2dc46cSBarry Smith        0,
1367cc2dc46cSBarry Smith        MatILUFactor_SeqBAIJ,
1368cc2dc46cSBarry Smith        0,
1369cc2dc46cSBarry Smith        0,
1370cc2dc46cSBarry Smith        MatGetSubMatrices_SeqBAIJ,
1371cc2dc46cSBarry Smith        MatIncreaseOverlap_SeqBAIJ,
1372cc2dc46cSBarry Smith        MatGetValues_SeqBAIJ,
1373cc2dc46cSBarry Smith        0,
1374cc2dc46cSBarry Smith        MatPrintHelp_SeqBAIJ,
1375cc2dc46cSBarry Smith        MatScale_SeqBAIJ,
1376cc2dc46cSBarry Smith        0,
1377cc2dc46cSBarry Smith        0,
1378cc2dc46cSBarry Smith        0,
1379cc2dc46cSBarry Smith        MatGetBlockSize_SeqBAIJ,
13803b2fbd54SBarry Smith        MatGetRowIJ_SeqBAIJ,
138192c4ed94SBarry Smith        MatRestoreRowIJ_SeqBAIJ,
1382cc2dc46cSBarry Smith        0,
1383cc2dc46cSBarry Smith        0,
1384cc2dc46cSBarry Smith        0,
1385cc2dc46cSBarry Smith        0,
1386cc2dc46cSBarry Smith        0,
1387cc2dc46cSBarry Smith        0,
1388d3825aa8SBarry Smith        MatSetValuesBlocked_SeqBAIJ,
1389cc2dc46cSBarry Smith        MatGetSubMatrix_SeqBAIJ,
1390b9b97703SBarry Smith        MatDestroy_SeqBAIJ,
1391b9b97703SBarry Smith        MatView_SeqBAIJ,
13923a64cc32SBarry Smith        MatGetPetscMaps_Petsc,
1393273d9f13SBarry Smith        0,
1394273d9f13SBarry Smith        0,
1395273d9f13SBarry Smith        0,
1396273d9f13SBarry Smith        0,
1397273d9f13SBarry Smith        0,
1398273d9f13SBarry Smith        0,
1399273d9f13SBarry Smith        MatGetRowMax_SeqBAIJ,
1400273d9f13SBarry Smith        MatConvert_Basic};
14012593348eSBarry Smith 
14023e90b805SBarry Smith EXTERN_C_BEGIN
14034a2ae208SSatish Balay #undef __FUNCT__
14044a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqBAIJ"
14053e90b805SBarry Smith int MatStoreValues_SeqBAIJ(Mat mat)
14063e90b805SBarry Smith {
14073e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1408273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1409d132466eSBarry Smith   int         ierr;
14103e90b805SBarry Smith 
14113e90b805SBarry Smith   PetscFunctionBegin;
14123e90b805SBarry Smith   if (aij->nonew != 1) {
141329bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
14143e90b805SBarry Smith   }
14153e90b805SBarry Smith 
14163e90b805SBarry Smith   /* allocate space for values if not already there */
14173e90b805SBarry Smith   if (!aij->saved_values) {
141887828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
14193e90b805SBarry Smith   }
14203e90b805SBarry Smith 
14213e90b805SBarry Smith   /* copy values over */
142287828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
14233e90b805SBarry Smith   PetscFunctionReturn(0);
14243e90b805SBarry Smith }
14253e90b805SBarry Smith EXTERN_C_END
14263e90b805SBarry Smith 
14273e90b805SBarry Smith EXTERN_C_BEGIN
14284a2ae208SSatish Balay #undef __FUNCT__
14294a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqBAIJ"
143032e87ba3SBarry Smith int MatRetrieveValues_SeqBAIJ(Mat mat)
14313e90b805SBarry Smith {
14323e90b805SBarry Smith   Mat_SeqBAIJ *aij = (Mat_SeqBAIJ *)mat->data;
1433273d9f13SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
14343e90b805SBarry Smith 
14353e90b805SBarry Smith   PetscFunctionBegin;
14363e90b805SBarry Smith   if (aij->nonew != 1) {
143729bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
14383e90b805SBarry Smith   }
14393e90b805SBarry Smith   if (!aij->saved_values) {
144029bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
14413e90b805SBarry Smith   }
14423e90b805SBarry Smith 
14433e90b805SBarry Smith   /* copy values over */
144487828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
14453e90b805SBarry Smith   PetscFunctionReturn(0);
14463e90b805SBarry Smith }
14473e90b805SBarry Smith EXTERN_C_END
14483e90b805SBarry Smith 
1449273d9f13SBarry Smith EXTERN_C_BEGIN
1450273d9f13SBarry Smith extern int MatConvert_SeqBAIJ_SeqAIJ(Mat,MatType,Mat*);
1451273d9f13SBarry Smith EXTERN_C_END
1452273d9f13SBarry Smith 
1453273d9f13SBarry Smith EXTERN_C_BEGIN
14544a2ae208SSatish Balay #undef __FUNCT__
14554a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqBAIJ"
1456273d9f13SBarry Smith int MatCreate_SeqBAIJ(Mat B)
14572593348eSBarry Smith {
1458273d9f13SBarry Smith   int         ierr,size;
1459b6490206SBarry Smith   Mat_SeqBAIJ *b;
14603b2fbd54SBarry Smith 
14613a40ed3dSBarry Smith   PetscFunctionBegin;
1462273d9f13SBarry Smith   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
146329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1464b6490206SBarry Smith 
1465273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1466273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1467b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqBAIJ,&b);CHKERRQ(ierr);
1468b0a32e0cSBarry Smith   B->data = (void*)b;
1469549d3d68SSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqBAIJ));CHKERRQ(ierr);
1470549d3d68SSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
14712593348eSBarry Smith   B->factor           = 0;
14722593348eSBarry Smith   B->lupivotthreshold = 1.0;
147390f02eecSBarry Smith   B->mapping          = 0;
14742593348eSBarry Smith   b->row              = 0;
14752593348eSBarry Smith   b->col              = 0;
1476e51c0b9cSSatish Balay   b->icol             = 0;
14772593348eSBarry Smith   b->reallocs         = 0;
14783e90b805SBarry Smith   b->saved_values     = 0;
1479cfce73b9SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
1480563b5814SBarry Smith   b->setvalueslen     = 0;
1481563b5814SBarry Smith   b->setvaluescopy    = PETSC_NULL;
1482563b5814SBarry Smith #endif
14838661488fSKris Buschelman   b->single_precision_solves = PETSC_FALSE;
14842593348eSBarry Smith 
14853a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr);
14863a64cc32SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr);
1487a5ae1ecdSBarry Smith 
1488c4992f7dSBarry Smith   b->sorted           = PETSC_FALSE;
1489c4992f7dSBarry Smith   b->roworiented      = PETSC_TRUE;
14902593348eSBarry Smith   b->nonew            = 0;
14912593348eSBarry Smith   b->diag             = 0;
14922593348eSBarry Smith   b->solve_work       = 0;
1493de6a44a3SBarry Smith   b->mult_work        = 0;
14942593348eSBarry Smith   b->spptr            = 0;
14950e6d2581SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1496883fce79SBarry Smith   b->keepzeroedrows   = PETSC_FALSE;
14974e220ebcSLois Curfman McInnes 
1498f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
14993e90b805SBarry Smith                                      "MatStoreValues_SeqBAIJ",
1500bc4b532fSSatish Balay                                       MatStoreValues_SeqBAIJ);CHKERRQ(ierr);
1501f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
15023e90b805SBarry Smith                                      "MatRetrieveValues_SeqBAIJ",
1503bc4b532fSSatish Balay                                       MatRetrieveValues_SeqBAIJ);CHKERRQ(ierr);
1504f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqBAIJSetColumnIndices_C",
150527a8da17SBarry Smith                                      "MatSeqBAIJSetColumnIndices_SeqBAIJ",
1506bc4b532fSSatish Balay                                       MatSeqBAIJSetColumnIndices_SeqBAIJ);CHKERRQ(ierr);
1507273d9f13SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_SeqBAIJ_SeqAIJ_C",
1508273d9f13SBarry Smith                                      "MatConvert_SeqBAIJ_SeqAIJ",
1509273d9f13SBarry Smith                                       MatConvert_SeqBAIJ_SeqAIJ);CHKERRQ(ierr);
15103a40ed3dSBarry Smith   PetscFunctionReturn(0);
15112593348eSBarry Smith }
1512273d9f13SBarry Smith EXTERN_C_END
15132593348eSBarry Smith 
15144a2ae208SSatish Balay #undef __FUNCT__
15154a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqBAIJ"
15162e8a6d31SBarry Smith int MatDuplicate_SeqBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
15172593348eSBarry Smith {
15182593348eSBarry Smith   Mat         C;
1519b6490206SBarry Smith   Mat_SeqBAIJ *c,*a = (Mat_SeqBAIJ*)A->data;
152027a8da17SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
1521de6a44a3SBarry Smith 
15223a40ed3dSBarry Smith   PetscFunctionBegin;
152329bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
15242593348eSBarry Smith 
15252593348eSBarry Smith   *B = 0;
1526273d9f13SBarry Smith   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1527273d9f13SBarry Smith   ierr = MatSetType(C,MATSEQBAIJ);CHKERRQ(ierr);
1528273d9f13SBarry Smith   c    = (Mat_SeqBAIJ*)C->data;
152944cd7ae7SLois Curfman McInnes 
1530b6490206SBarry Smith   c->bs         = a->bs;
15319df24120SSatish Balay   c->bs2        = a->bs2;
1532b6490206SBarry Smith   c->mbs        = a->mbs;
1533b6490206SBarry Smith   c->nbs        = a->nbs;
153435d8aa7fSBarry Smith   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
15352593348eSBarry Smith 
1536b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1537b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
1538b6490206SBarry Smith   for (i=0; i<mbs; i++) {
15392593348eSBarry Smith     c->imax[i] = a->imax[i];
15402593348eSBarry Smith     c->ilen[i] = a->ilen[i];
15412593348eSBarry Smith   }
15422593348eSBarry Smith 
15432593348eSBarry Smith   /* allocate the matrix space */
1544c4992f7dSBarry Smith   c->singlemalloc = PETSC_TRUE;
15453f1db9ecSBarry Smith   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
1546b0a32e0cSBarry Smith   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
15477e67e3f9SSatish Balay   c->j = (int*)(c->a + nz*bs2);
1548de6a44a3SBarry Smith   c->i = c->j + nz;
1549549d3d68SSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
1550b6490206SBarry Smith   if (mbs > 0) {
1551549d3d68SSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
15522e8a6d31SBarry Smith     if (cpvalues == MAT_COPY_VALUES) {
1553549d3d68SSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
15542e8a6d31SBarry Smith     } else {
1555549d3d68SSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
15562593348eSBarry Smith     }
15572593348eSBarry Smith   }
15582593348eSBarry Smith 
1559b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
15602593348eSBarry Smith   c->sorted      = a->sorted;
15612593348eSBarry Smith   c->roworiented = a->roworiented;
15622593348eSBarry Smith   c->nonew       = a->nonew;
15632593348eSBarry Smith 
15642593348eSBarry Smith   if (a->diag) {
1565b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1566b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
1567b6490206SBarry Smith     for (i=0; i<mbs; i++) {
15682593348eSBarry Smith       c->diag[i] = a->diag[i];
15692593348eSBarry Smith     }
157098305bb5SBarry Smith   } else c->diag        = 0;
15712593348eSBarry Smith   c->nz                 = a->nz;
15722593348eSBarry Smith   c->maxnz              = a->maxnz;
15732593348eSBarry Smith   c->solve_work         = 0;
15742593348eSBarry Smith   c->spptr              = 0;     /* Dangerous -I'm throwing away a->spptr */
15757fc0212eSBarry Smith   c->mult_work          = 0;
1576273d9f13SBarry Smith   C->preallocated       = PETSC_TRUE;
1577273d9f13SBarry Smith   C->assembled          = PETSC_TRUE;
15782593348eSBarry Smith   *B = C;
1579b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
15803a40ed3dSBarry Smith   PetscFunctionReturn(0);
15812593348eSBarry Smith }
15822593348eSBarry Smith 
1583273d9f13SBarry Smith EXTERN_C_BEGIN
15844a2ae208SSatish Balay #undef __FUNCT__
15854a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqBAIJ"
1586b0a32e0cSBarry Smith int MatLoad_SeqBAIJ(PetscViewer viewer,MatType type,Mat *A)
15872593348eSBarry Smith {
1588b6490206SBarry Smith   Mat_SeqBAIJ  *a;
15892593348eSBarry Smith   Mat          B;
1590f1af5d2fSBarry Smith   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1591b6490206SBarry Smith   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
159235aab85fSBarry Smith   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
1593a7c10996SSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
159487828ca2SBarry Smith   PetscScalar  *aa;
159519bcc07fSBarry Smith   MPI_Comm     comm = ((PetscObject)viewer)->comm;
15962593348eSBarry Smith 
15973a40ed3dSBarry Smith   PetscFunctionBegin;
1598b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1599de6a44a3SBarry Smith   bs2  = bs*bs;
1600b6490206SBarry Smith 
1601d132466eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
160229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1603b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
16040752156aSBarry Smith   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1605552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
16062593348eSBarry Smith   M = header[1]; N = header[2]; nz = header[3];
16072593348eSBarry Smith 
1608d64ed03dSBarry Smith   if (header[3] < 0) {
160929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqBAIJ");
1610d64ed03dSBarry Smith   }
1611d64ed03dSBarry Smith 
161229bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
161335aab85fSBarry Smith 
161435aab85fSBarry Smith   /*
161535aab85fSBarry Smith      This code adds extra rows to make sure the number of rows is
161635aab85fSBarry Smith     divisible by the blocksize
161735aab85fSBarry Smith   */
1618b6490206SBarry Smith   mbs        = M/bs;
161935aab85fSBarry Smith   extra_rows = bs - M + bs*(mbs);
162035aab85fSBarry Smith   if (extra_rows == bs) extra_rows = 0;
162135aab85fSBarry Smith   else                  mbs++;
162235aab85fSBarry Smith   if (extra_rows) {
1623b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqBAIJ:Padding loaded matrix to match blocksize\n");
162435aab85fSBarry Smith   }
1625b6490206SBarry Smith 
16262593348eSBarry Smith   /* read in row lengths */
1627b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
16280752156aSBarry Smith   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
162935aab85fSBarry Smith   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
16302593348eSBarry Smith 
1631b6490206SBarry Smith   /* read in column indices */
1632b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
16330752156aSBarry Smith   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
163435aab85fSBarry Smith   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
1635b6490206SBarry Smith 
1636b6490206SBarry Smith   /* loop over row lengths determining block row lengths */
1637b0a32e0cSBarry Smith   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
1638549d3d68SSatish Balay   ierr     = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
1639b0a32e0cSBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
1640549d3d68SSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
164135aab85fSBarry Smith   masked   = mask + mbs;
1642b6490206SBarry Smith   rowcount = 0; nzcount = 0;
1643b6490206SBarry Smith   for (i=0; i<mbs; i++) {
164435aab85fSBarry Smith     nmask = 0;
1645b6490206SBarry Smith     for (j=0; j<bs; j++) {
1646b6490206SBarry Smith       kmax = rowlengths[rowcount];
1647b6490206SBarry Smith       for (k=0; k<kmax; k++) {
164835aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
164935aab85fSBarry Smith         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
1650b6490206SBarry Smith       }
1651b6490206SBarry Smith       rowcount++;
1652b6490206SBarry Smith     }
165335aab85fSBarry Smith     browlengths[i] += nmask;
165435aab85fSBarry Smith     /* zero out the mask elements we set */
165535aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1656b6490206SBarry Smith   }
1657b6490206SBarry Smith 
16582593348eSBarry Smith   /* create our matrix */
16593eee16b0SBarry Smith   ierr = MatCreateSeqBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
16602593348eSBarry Smith   B = *A;
1661b6490206SBarry Smith   a = (Mat_SeqBAIJ*)B->data;
16622593348eSBarry Smith 
16632593348eSBarry Smith   /* set matrix "i" values */
1664de6a44a3SBarry Smith   a->i[0] = 0;
1665b6490206SBarry Smith   for (i=1; i<= mbs; i++) {
1666b6490206SBarry Smith     a->i[i]      = a->i[i-1] + browlengths[i-1];
1667b6490206SBarry Smith     a->ilen[i-1] = browlengths[i-1];
16682593348eSBarry Smith   }
1669b6490206SBarry Smith   a->nz         = 0;
1670b6490206SBarry Smith   for (i=0; i<mbs; i++) a->nz += browlengths[i];
16712593348eSBarry Smith 
1672b6490206SBarry Smith   /* read in nonzero values */
167387828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
16740752156aSBarry Smith   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
167535aab85fSBarry Smith   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
1676b6490206SBarry Smith 
1677b6490206SBarry Smith   /* set "a" and "j" values into matrix */
1678b6490206SBarry Smith   nzcount = 0; jcount = 0;
1679b6490206SBarry Smith   for (i=0; i<mbs; i++) {
1680b6490206SBarry Smith     nzcountb = nzcount;
168135aab85fSBarry Smith     nmask    = 0;
1682b6490206SBarry Smith     for (j=0; j<bs; j++) {
1683b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1684b6490206SBarry Smith       for (k=0; k<kmax; k++) {
168535aab85fSBarry Smith         tmp = jj[nzcount++]/bs;
168635aab85fSBarry Smith 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
1687b6490206SBarry Smith       }
1688b6490206SBarry Smith     }
1689de6a44a3SBarry Smith     /* sort the masked values */
1690433994e6SBarry Smith     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
1691de6a44a3SBarry Smith 
1692b6490206SBarry Smith     /* set "j" values into matrix */
1693b6490206SBarry Smith     maskcount = 1;
169435aab85fSBarry Smith     for (j=0; j<nmask; j++) {
169535aab85fSBarry Smith       a->j[jcount++]  = masked[j];
1696de6a44a3SBarry Smith       mask[masked[j]] = maskcount++;
1697b6490206SBarry Smith     }
1698b6490206SBarry Smith     /* set "a" values into matrix */
1699de6a44a3SBarry Smith     ishift = bs2*a->i[i];
1700b6490206SBarry Smith     for (j=0; j<bs; j++) {
1701b6490206SBarry Smith       kmax = rowlengths[i*bs+j];
1702b6490206SBarry Smith       for (k=0; k<kmax; k++) {
1703de6a44a3SBarry Smith         tmp       = jj[nzcountb]/bs ;
1704de6a44a3SBarry Smith         block     = mask[tmp] - 1;
1705de6a44a3SBarry Smith         point     = jj[nzcountb] - bs*tmp;
1706de6a44a3SBarry Smith         idx       = ishift + bs2*block + j + bs*point;
1707375fe846SBarry Smith         a->a[idx] = (MatScalar)aa[nzcountb++];
1708b6490206SBarry Smith       }
1709b6490206SBarry Smith     }
171035aab85fSBarry Smith     /* zero out the mask elements we set */
171135aab85fSBarry Smith     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
1712b6490206SBarry Smith   }
171329bbc08cSBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
1714b6490206SBarry Smith 
1715606d414cSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1716606d414cSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1717606d414cSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
1718606d414cSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
1719606d414cSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
1720b6490206SBarry Smith 
1721b6490206SBarry Smith   B->assembled = PETSC_TRUE;
1722de6a44a3SBarry Smith 
17239c01be13SBarry Smith   ierr = MatView_Private(B);CHKERRQ(ierr);
17243a40ed3dSBarry Smith   PetscFunctionReturn(0);
17252593348eSBarry Smith }
1726273d9f13SBarry Smith EXTERN_C_END
17272593348eSBarry Smith 
17284a2ae208SSatish Balay #undef __FUNCT__
17294a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqBAIJ"
1730273d9f13SBarry Smith /*@C
1731273d9f13SBarry Smith    MatCreateSeqBAIJ - Creates a sparse matrix in block AIJ (block
1732273d9f13SBarry Smith    compressed row) format.  For good matrix assembly performance the
1733273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1734273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1735273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
17362593348eSBarry Smith 
1737273d9f13SBarry Smith    Collective on MPI_Comm
1738273d9f13SBarry Smith 
1739273d9f13SBarry Smith    Input Parameters:
1740273d9f13SBarry Smith +  comm - MPI communicator, set to PETSC_COMM_SELF
1741273d9f13SBarry Smith .  bs - size of block
1742273d9f13SBarry Smith .  m - number of rows
1743273d9f13SBarry Smith .  n - number of columns
174435d8aa7fSBarry Smith .  nz - number of nonzero blocks  per block row (same for all rows)
174535d8aa7fSBarry Smith -  nnz - array containing the number of nonzero blocks in the various block rows
1746273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1747273d9f13SBarry Smith 
1748273d9f13SBarry Smith    Output Parameter:
1749273d9f13SBarry Smith .  A - the matrix
1750273d9f13SBarry Smith 
1751273d9f13SBarry Smith    Options Database Keys:
1752273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1753273d9f13SBarry Smith                      block calculations (much slower)
1754273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1755273d9f13SBarry Smith 
1756273d9f13SBarry Smith    Level: intermediate
1757273d9f13SBarry Smith 
1758273d9f13SBarry Smith    Notes:
175935d8aa7fSBarry Smith    A nonzero block is any block that as 1 or more nonzeros in it
176035d8aa7fSBarry Smith 
1761273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1762273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1763273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1764273d9f13SBarry Smith 
1765273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1766273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1767273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1768273d9f13SBarry Smith    matrices.
1769273d9f13SBarry Smith 
1770273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1771273d9f13SBarry Smith @*/
1772273d9f13SBarry Smith int MatCreateSeqBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1773273d9f13SBarry Smith {
1774273d9f13SBarry Smith   int ierr;
1775273d9f13SBarry Smith 
1776273d9f13SBarry Smith   PetscFunctionBegin;
1777273d9f13SBarry Smith   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1778273d9f13SBarry Smith   ierr = MatSetType(*A,MATSEQBAIJ);CHKERRQ(ierr);
1779273d9f13SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
1780273d9f13SBarry Smith   PetscFunctionReturn(0);
1781273d9f13SBarry Smith }
1782273d9f13SBarry Smith 
17834a2ae208SSatish Balay #undef __FUNCT__
17844a2ae208SSatish Balay #define __FUNCT__ "MatSeqBAIJSetPreallocation"
1785273d9f13SBarry Smith /*@C
1786273d9f13SBarry Smith    MatSeqBAIJSetPreallocation - Sets the block size and expected nonzeros
1787273d9f13SBarry Smith    per row in the matrix. For good matrix assembly performance the
1788273d9f13SBarry Smith    user should preallocate the matrix storage by setting the parameter nz
1789273d9f13SBarry Smith    (or the array nnz).  By setting these parameters accurately, performance
1790273d9f13SBarry Smith    during matrix assembly can be increased by more than a factor of 50.
1791273d9f13SBarry Smith 
1792273d9f13SBarry Smith    Collective on MPI_Comm
1793273d9f13SBarry Smith 
1794273d9f13SBarry Smith    Input Parameters:
1795273d9f13SBarry Smith +  A - the matrix
1796273d9f13SBarry Smith .  bs - size of block
1797273d9f13SBarry Smith .  nz - number of block nonzeros per block row (same for all rows)
1798273d9f13SBarry Smith -  nnz - array containing the number of block nonzeros in the various block rows
1799273d9f13SBarry Smith          (possibly different for each block row) or PETSC_NULL
1800273d9f13SBarry Smith 
1801273d9f13SBarry Smith    Options Database Keys:
1802273d9f13SBarry Smith .   -mat_no_unroll - uses code that does not unroll the loops in the
1803273d9f13SBarry Smith                      block calculations (much slower)
1804273d9f13SBarry Smith .    -mat_block_size - size of the blocks to use
1805273d9f13SBarry Smith 
1806273d9f13SBarry Smith    Level: intermediate
1807273d9f13SBarry Smith 
1808273d9f13SBarry Smith    Notes:
1809273d9f13SBarry Smith    The block AIJ format is fully compatible with standard Fortran 77
1810273d9f13SBarry Smith    storage.  That is, the stored row and column indices can begin at
1811273d9f13SBarry Smith    either one (as in Fortran) or zero.  See the users' manual for details.
1812273d9f13SBarry Smith 
1813273d9f13SBarry Smith    Specify the preallocated storage with either nz or nnz (not both).
1814273d9f13SBarry Smith    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1815273d9f13SBarry Smith    allocation.  For additional details, see the users manual chapter on
1816273d9f13SBarry Smith    matrices.
1817273d9f13SBarry Smith 
1818273d9f13SBarry Smith .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1819273d9f13SBarry Smith @*/
1820273d9f13SBarry Smith int MatSeqBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
1821273d9f13SBarry Smith {
1822273d9f13SBarry Smith   Mat_SeqBAIJ *b;
1823273d9f13SBarry Smith   int         i,len,ierr,mbs,nbs,bs2,newbs = bs;
1824273d9f13SBarry Smith   PetscTruth  flg;
1825273d9f13SBarry Smith 
1826273d9f13SBarry Smith   PetscFunctionBegin;
1827273d9f13SBarry Smith   ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flg);CHKERRQ(ierr);
1828273d9f13SBarry Smith   if (!flg) PetscFunctionReturn(0);
1829273d9f13SBarry Smith 
1830273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1831b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&newbs,PETSC_NULL);CHKERRQ(ierr);
1832273d9f13SBarry Smith   if (nnz && newbs != bs) {
1833273d9f13SBarry Smith     SETERRQ(1,"Cannot change blocksize from command line if setting nnz");
1834273d9f13SBarry Smith   }
1835273d9f13SBarry Smith   bs   = newbs;
1836273d9f13SBarry Smith 
1837273d9f13SBarry Smith   mbs  = B->m/bs;
1838273d9f13SBarry Smith   nbs  = B->n/bs;
1839273d9f13SBarry Smith   bs2  = bs*bs;
1840273d9f13SBarry Smith 
1841273d9f13SBarry Smith   if (mbs*bs!=B->m || nbs*bs!=B->n) {
184235d8aa7fSBarry Smith     SETERRQ3(PETSC_ERR_ARG_SIZ,"Number rows %d, cols %d must be divisible by blocksize %d",B->m,B->n,bs);
1843273d9f13SBarry Smith   }
1844273d9f13SBarry Smith 
1845435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1846435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
1847273d9f13SBarry Smith   if (nnz) {
1848273d9f13SBarry Smith     for (i=0; i<mbs; i++) {
1849273d9f13SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
18503a7fca6bSBarry 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);
1851273d9f13SBarry Smith     }
1852273d9f13SBarry Smith   }
1853273d9f13SBarry Smith 
1854273d9f13SBarry Smith   b       = (Mat_SeqBAIJ*)B->data;
1855b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
185654138f6bSKris Buschelman   B->ops->solve               = MatSolve_SeqBAIJ_Update;
185754138f6bSKris Buschelman   B->ops->solvetranspose      = MatSolveTranspose_SeqBAIJ_Update;
1858273d9f13SBarry Smith   if (!flg) {
1859273d9f13SBarry Smith     switch (bs) {
1860273d9f13SBarry Smith     case 1:
1861273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_1;
1862273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_1;
1863273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_1;
1864273d9f13SBarry Smith       break;
1865273d9f13SBarry Smith     case 2:
1866273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2;
1867273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_2;
1868273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_2;
1869273d9f13SBarry Smith       break;
1870273d9f13SBarry Smith     case 3:
1871273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3;
1872273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_3;
1873273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_3;
1874273d9f13SBarry Smith       break;
1875273d9f13SBarry Smith     case 4:
1876273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4;
1877273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_4;
1878273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_4;
1879273d9f13SBarry Smith       break;
1880273d9f13SBarry Smith     case 5:
1881273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5;
1882273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_5;
1883273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_5;
1884273d9f13SBarry Smith       break;
1885273d9f13SBarry Smith     case 6:
1886273d9f13SBarry Smith       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6;
1887273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_6;
1888273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_6;
1889273d9f13SBarry Smith       break;
1890273d9f13SBarry Smith     case 7:
189154138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7;
1892273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_7;
1893273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_7;
1894273d9f13SBarry Smith       break;
1895273d9f13SBarry Smith     default:
189654138f6bSKris Buschelman       B->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_N;
1897273d9f13SBarry Smith       B->ops->mult            = MatMult_SeqBAIJ_N;
1898273d9f13SBarry Smith       B->ops->multadd         = MatMultAdd_SeqBAIJ_N;
1899273d9f13SBarry Smith       break;
1900273d9f13SBarry Smith     }
1901273d9f13SBarry Smith   }
1902273d9f13SBarry Smith   b->bs      = bs;
1903273d9f13SBarry Smith   b->mbs     = mbs;
1904273d9f13SBarry Smith   b->nbs     = nbs;
1905b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
1906273d9f13SBarry Smith   if (!nnz) {
1907435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1908273d9f13SBarry Smith     else if (nz <= 0)        nz = 1;
1909273d9f13SBarry Smith     for (i=0; i<mbs; i++) b->imax[i] = nz;
1910273d9f13SBarry Smith     nz = nz*mbs;
1911273d9f13SBarry Smith   } else {
1912273d9f13SBarry Smith     nz = 0;
1913273d9f13SBarry Smith     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1914273d9f13SBarry Smith   }
1915273d9f13SBarry Smith 
1916273d9f13SBarry Smith   /* allocate the matrix space */
1917273d9f13SBarry Smith   len   = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
1918b0a32e0cSBarry Smith   ierr  = PetscMalloc(len,&b->a);CHKERRQ(ierr);
1919273d9f13SBarry Smith   ierr  = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1920273d9f13SBarry Smith   b->j  = (int*)(b->a + nz*bs2);
1921273d9f13SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
1922273d9f13SBarry Smith   b->i  = b->j + nz;
1923273d9f13SBarry Smith   b->singlemalloc = PETSC_TRUE;
1924273d9f13SBarry Smith 
1925273d9f13SBarry Smith   b->i[0] = 0;
1926273d9f13SBarry Smith   for (i=1; i<mbs+1; i++) {
1927273d9f13SBarry Smith     b->i[i] = b->i[i-1] + b->imax[i-1];
1928273d9f13SBarry Smith   }
1929273d9f13SBarry Smith 
1930273d9f13SBarry Smith   /* b->ilen will count nonzeros in each block row so far. */
193116d6aa21SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1932b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqBAIJ));
1933273d9f13SBarry Smith   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
1934273d9f13SBarry Smith 
1935273d9f13SBarry Smith   b->bs               = bs;
1936273d9f13SBarry Smith   b->bs2              = bs2;
1937273d9f13SBarry Smith   b->mbs              = mbs;
1938273d9f13SBarry Smith   b->nz               = 0;
1939273d9f13SBarry Smith   b->maxnz            = nz*bs2;
1940273d9f13SBarry Smith   B->info.nz_unneeded = (PetscReal)b->maxnz;
1941273d9f13SBarry Smith   PetscFunctionReturn(0);
1942273d9f13SBarry Smith }
19432593348eSBarry Smith 
1944