xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 4ce68768e41e382c1fc596e135fd299def09fddc)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay /*
3a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
449b5e25fSSatish Balay   matrix storage format.
549b5e25fSSatish Balay */
69e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h"         /*I "petscmat.h" I*/
749b5e25fSSatish Balay #include "src/inline/spops.h"
8aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
949b5e25fSSatish Balay 
1049b5e25fSSatish Balay #define CHUNKSIZE  10
1149b5e25fSSatish Balay 
1249b5e25fSSatish Balay /*
1349b5e25fSSatish Balay      Checks for missing diagonals
1449b5e25fSSatish Balay */
154a2ae208SSatish Balay #undef __FUNCT__
164a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
17dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A)
1849b5e25fSSatish Balay {
19045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
206849ba73SBarry Smith   PetscErrorCode ierr;
216849ba73SBarry Smith   int         *diag,*jj = a->j,i;
2249b5e25fSSatish Balay 
2349b5e25fSSatish Balay   PetscFunctionBegin;
24045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2549b5e25fSSatish Balay   diag = a->diag;
2649b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
27e11b0e14SHong Zhang     if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i);
2849b5e25fSSatish Balay   }
2949b5e25fSSatish Balay   PetscFunctionReturn(0);
3049b5e25fSSatish Balay }
3149b5e25fSSatish Balay 
324a2ae208SSatish Balay #undef __FUNCT__
334a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
34dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
3549b5e25fSSatish Balay {
36045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
376849ba73SBarry Smith   PetscErrorCode ierr;
386849ba73SBarry Smith   int          i,mbs = a->mbs;
3949b5e25fSSatish Balay 
4049b5e25fSSatish Balay   PetscFunctionBegin;
4149b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4249b5e25fSSatish Balay 
43b424e231SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr);
44b424e231SHong Zhang   PetscLogObjectMemory(A,(mbs+1)*sizeof(int));
45b424e231SHong Zhang   for (i=0; i<mbs; i++) a->diag[i] = a->i[i];
4649b5e25fSSatish Balay   PetscFunctionReturn(0);
4749b5e25fSSatish Balay }
4849b5e25fSSatish Balay 
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
516849ba73SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done)
5249b5e25fSSatish Balay {
53a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
549e37e433SSatish Balay   int         n = a->mbs,i;
5549b5e25fSSatish Balay 
5649b5e25fSSatish Balay   PetscFunctionBegin;
57d3e5a4abSHong Zhang   *nn = n;
58a1373b80SHong Zhang   if (!ia) PetscFunctionReturn(0);
59a1373b80SHong Zhang 
60a6ece127SHong Zhang   if (oshift == 1) {
61a6ece127SHong Zhang     /* temporarily add 1 to i and j indices */
626c6c5352SBarry Smith     int nz = a->i[n];
636c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
64a1373b80SHong Zhang     for (i=0; i<n+1; i++) a->i[i]++;
65a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
66a1373b80SHong Zhang   } else {
67a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
68a6ece127SHong Zhang   }
6949b5e25fSSatish Balay   PetscFunctionReturn(0);
7049b5e25fSSatish Balay }
7149b5e25fSSatish Balay 
724a2ae208SSatish Balay #undef __FUNCT__
734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
746849ba73SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done)
7549b5e25fSSatish Balay {
76b7aaefc3SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
77a15ac5e4SHong Zhang   int         i,n = a->mbs;
78a6ece127SHong Zhang 
7949b5e25fSSatish Balay   PetscFunctionBegin;
8049b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
81a6ece127SHong Zhang 
82a6ece127SHong Zhang   if (oshift == 1) {
836c6c5352SBarry Smith     int nz = a->i[n]-1;
846c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
85a6ece127SHong Zhang     for (i=0; i<n+1; i++) a->i[i]--;
86a6ece127SHong Zhang   }
87a6ece127SHong Zhang   PetscFunctionReturn(0);
8849b5e25fSSatish Balay }
8949b5e25fSSatish Balay 
904a2ae208SSatish Balay #undef __FUNCT__
914a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ"
92dfbe8321SBarry Smith PetscErrorCode MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
9349b5e25fSSatish Balay {
9449b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
9549b5e25fSSatish Balay 
9649b5e25fSSatish Balay   PetscFunctionBegin;
9749b5e25fSSatish Balay   *bs = sbaij->bs;
9849b5e25fSSatish Balay   PetscFunctionReturn(0);
9949b5e25fSSatish Balay }
10049b5e25fSSatish Balay 
1014a2ae208SSatish Balay #undef __FUNCT__
1024a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
103dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
10449b5e25fSSatish Balay {
10549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
106dfbe8321SBarry Smith   PetscErrorCode ierr;
10749b5e25fSSatish Balay 
10849b5e25fSSatish Balay   PetscFunctionBegin;
10949b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
1106c6c5352SBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, NZ=%d",A->m,a->nz);
11149b5e25fSSatish Balay #endif
11249b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
11349b5e25fSSatish Balay   if (!a->singlemalloc) {
11449b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
11563c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
11649b5e25fSSatish Balay   }
11749b5e25fSSatish Balay   if (a->row) {
11849b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
11949b5e25fSSatish Balay   }
12049b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
12149b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
12249b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
12349b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
124d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
125d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
126d59c15a7SBarry Smith   if (a->mult_work)   {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
12749b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
128c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
1291a3463dfSHong Zhang 
1301a3463dfSHong Zhang   if (a->inew){
1311a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1321a3463dfSHong Zhang     a->inew = 0;
1331a3463dfSHong Zhang   }
13449b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13549b5e25fSSatish Balay   PetscFunctionReturn(0);
13649b5e25fSSatish Balay }
13749b5e25fSSatish Balay 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
140dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14149b5e25fSSatish Balay {
142045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14349b5e25fSSatish Balay 
14449b5e25fSSatish Balay   PetscFunctionBegin;
1454d9d31abSKris Buschelman   switch (op) {
1464d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1474d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1484d9d31abSKris Buschelman     break;
1494d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1504d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1514d9d31abSKris Buschelman     break;
1524d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1534d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1564d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1594d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1624d9d31abSKris Buschelman     a->nonew          = 1;
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1654d9d31abSKris Buschelman     a->nonew          = -1;
1664d9d31abSKris Buschelman     break;
1674d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1684d9d31abSKris Buschelman     a->nonew          = -2;
1694d9d31abSKris Buschelman     break;
1704d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1714d9d31abSKris Buschelman     a->nonew          = 0;
1724d9d31abSKris Buschelman     break;
1734d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1744d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1754d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1764d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1774d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
178b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1794d9d31abSKris Buschelman     break;
1804d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
18129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1829a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
1839a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1849a4540c5SBarry Smith   case MAT_HERMITIAN:
1859a4540c5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
18677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
18777e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
1889a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
1899a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1909a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
19177e54ba9SKris Buschelman     break;
1924d9d31abSKris Buschelman   default:
19329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
19449b5e25fSSatish Balay   }
19549b5e25fSSatish Balay   PetscFunctionReturn(0);
19649b5e25fSSatish Balay }
19749b5e25fSSatish Balay 
1984a2ae208SSatish Balay #undef __FUNCT__
1994a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
200dfbe8321SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
20149b5e25fSSatish Balay {
20249b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
2036849ba73SBarry Smith   PetscErrorCode ierr;
2046849ba73SBarry Smith   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
20549b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
20687828ca2SBarry Smith   PetscScalar  *v_i;
20749b5e25fSSatish Balay 
20849b5e25fSSatish Balay   PetscFunctionBegin;
20949b5e25fSSatish Balay   bs  = a->bs;
21049b5e25fSSatish Balay   ai  = a->i;
21149b5e25fSSatish Balay   aj  = a->j;
21249b5e25fSSatish Balay   aa  = a->a;
21349b5e25fSSatish Balay   bs2 = a->bs2;
21449b5e25fSSatish Balay 
215a45adfd6SMatthew Knepley   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %d out of range", row);
21649b5e25fSSatish Balay 
21749b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21849b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21949b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
22049b5e25fSSatish Balay   *ncols = bs*M;
22149b5e25fSSatish Balay 
22249b5e25fSSatish Balay   if (v) {
22349b5e25fSSatish Balay     *v = 0;
22449b5e25fSSatish Balay     if (*ncols) {
22587828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
22649b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22749b5e25fSSatish Balay         v_i  = *v + i*bs;
22849b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22949b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
23049b5e25fSSatish Balay       }
23149b5e25fSSatish Balay     }
23249b5e25fSSatish Balay   }
23349b5e25fSSatish Balay 
23449b5e25fSSatish Balay   if (cols) {
23549b5e25fSSatish Balay     *cols = 0;
23649b5e25fSSatish Balay     if (*ncols) {
23782502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
23849b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23949b5e25fSSatish Balay         cols_i = *cols + i*bs;
24049b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24149b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24249b5e25fSSatish Balay       }
24349b5e25fSSatish Balay     }
24449b5e25fSSatish Balay   }
24549b5e25fSSatish Balay 
24649b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2475ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2485ddb2528SHong Zhang #ifdef column_search
24949b5e25fSSatish Balay   v_i    = *v    + M*bs;
25049b5e25fSSatish Balay   cols_i = *cols + M*bs;
25149b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25249b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25349b5e25fSSatish Balay     for (j=0; j<M; j++){
25449b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25549b5e25fSSatish Balay       if (itmp == bn){
25649b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25749b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25849b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25949b5e25fSSatish Balay           *v_i++    = aa_i[k];
26049b5e25fSSatish Balay         }
26149b5e25fSSatish Balay         *ncols += bs;
26249b5e25fSSatish Balay         break;
26349b5e25fSSatish Balay       }
26449b5e25fSSatish Balay     }
26549b5e25fSSatish Balay   }
2665ddb2528SHong Zhang #endif
26749b5e25fSSatish Balay 
26849b5e25fSSatish Balay   PetscFunctionReturn(0);
26949b5e25fSSatish Balay }
27049b5e25fSSatish Balay 
2714a2ae208SSatish Balay #undef __FUNCT__
2724a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
273dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
27449b5e25fSSatish Balay {
275dfbe8321SBarry Smith   PetscErrorCode ierr;
27649b5e25fSSatish Balay 
27749b5e25fSSatish Balay   PetscFunctionBegin;
27849b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27949b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
28049b5e25fSSatish Balay   PetscFunctionReturn(0);
28149b5e25fSSatish Balay }
28249b5e25fSSatish Balay 
2834a2ae208SSatish Balay #undef __FUNCT__
2844a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
285dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B)
28649b5e25fSSatish Balay {
287dfbe8321SBarry Smith   PetscErrorCode ierr;
28849b5e25fSSatish Balay   PetscFunctionBegin;
289999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2908115998fSBarry Smith   PetscFunctionReturn(0);
29149b5e25fSSatish Balay }
29249b5e25fSSatish Balay 
2934a2ae208SSatish Balay #undef __FUNCT__
2944a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
2956849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
29649b5e25fSSatish Balay {
29749b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
298dfbe8321SBarry Smith   PetscErrorCode ierr;
299dfbe8321SBarry Smith   int i,j,bs = a->bs,k,l,bs2=a->bs2;
300fb9695e5SSatish Balay   char              *name;
301f3ef73ceSBarry Smith   PetscViewerFormat format;
30249b5e25fSSatish Balay 
30349b5e25fSSatish Balay   PetscFunctionBegin;
30480fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
305b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
306456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
307b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
308fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
30929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
310fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
311b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
31249b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
31349b5e25fSSatish Balay       for (j=0; j<bs; j++) {
314b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
31549b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
31649b5e25fSSatish Balay           for (l=0; l<bs; l++) {
31749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
31849b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
31962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
32049b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
32149b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
32262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
32349b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
32449b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
32562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
32649b5e25fSSatish Balay             }
32749b5e25fSSatish Balay #else
32849b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
32962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
33049b5e25fSSatish Balay             }
33149b5e25fSSatish Balay #endif
33249b5e25fSSatish Balay           }
33349b5e25fSSatish Balay         }
334b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
33549b5e25fSSatish Balay       }
33649b5e25fSSatish Balay     }
337b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
33849b5e25fSSatish Balay   } else {
339b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
34049b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
34149b5e25fSSatish Balay       for (j=0; j<bs; j++) {
342b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
34349b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
34449b5e25fSSatish Balay           for (l=0; l<bs; l++) {
34549b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
34649b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
34762b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
34849b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34949b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
35062b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
35149b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35249b5e25fSSatish Balay             } else {
35362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35449b5e25fSSatish Balay             }
35549b5e25fSSatish Balay #else
356b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
35749b5e25fSSatish Balay #endif
35849b5e25fSSatish Balay           }
35949b5e25fSSatish Balay         }
360b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
36149b5e25fSSatish Balay       }
36249b5e25fSSatish Balay     }
363b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
36449b5e25fSSatish Balay   }
365b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
36649b5e25fSSatish Balay   PetscFunctionReturn(0);
36749b5e25fSSatish Balay }
36849b5e25fSSatish Balay 
3694a2ae208SSatish Balay #undef __FUNCT__
3704a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
3716849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
37249b5e25fSSatish Balay {
37349b5e25fSSatish Balay   Mat          A = (Mat) Aa;
37449b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
3756849ba73SBarry Smith   PetscErrorCode ierr;
3766849ba73SBarry Smith   int          row,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
37749b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
37849b5e25fSSatish Balay   MatScalar    *aa;
37949b5e25fSSatish Balay   MPI_Comm     comm;
380b0a32e0cSBarry Smith   PetscViewer  viewer;
38149b5e25fSSatish Balay 
38249b5e25fSSatish Balay   PetscFunctionBegin;
38349b5e25fSSatish Balay   /*
38449b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
38549b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
38649b5e25fSSatish Balay    rest should return immediately.
38749b5e25fSSatish Balay   */
38849b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
38949b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
39049b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
39149b5e25fSSatish Balay 
39249b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
39349b5e25fSSatish Balay 
394b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
395b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
39649b5e25fSSatish Balay 
39749b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
398b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
39949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
40049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
401c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
40249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
40349b5e25fSSatish Balay       aa = a->a + j*bs2;
40449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
40549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
40649b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
407b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
40849b5e25fSSatish Balay         }
40949b5e25fSSatish Balay       }
41049b5e25fSSatish Balay     }
41149b5e25fSSatish Balay   }
412b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
41349b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
41449b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
415c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
41649b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
41749b5e25fSSatish Balay       aa = a->a + j*bs2;
41849b5e25fSSatish Balay       for (k=0; k<bs; k++) {
41949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
42049b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
421b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
42249b5e25fSSatish Balay         }
42349b5e25fSSatish Balay       }
42449b5e25fSSatish Balay     }
42549b5e25fSSatish Balay   }
42649b5e25fSSatish Balay 
427b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
42849b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
42949b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
430c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
43149b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
43249b5e25fSSatish Balay       aa = a->a + j*bs2;
43349b5e25fSSatish Balay       for (k=0; k<bs; k++) {
43449b5e25fSSatish Balay         for (l=0; l<bs; l++) {
43549b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
436b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
43749b5e25fSSatish Balay         }
43849b5e25fSSatish Balay       }
43949b5e25fSSatish Balay     }
44049b5e25fSSatish Balay   }
44149b5e25fSSatish Balay   PetscFunctionReturn(0);
44249b5e25fSSatish Balay }
44349b5e25fSSatish Balay 
4444a2ae208SSatish Balay #undef __FUNCT__
4454a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
4466849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
44749b5e25fSSatish Balay {
448dfbe8321SBarry Smith   PetscErrorCode ierr;
44949b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
450b0a32e0cSBarry Smith   PetscDraw    draw;
45149b5e25fSSatish Balay   PetscTruth   isnull;
45249b5e25fSSatish Balay 
45349b5e25fSSatish Balay   PetscFunctionBegin;
45449b5e25fSSatish Balay 
455b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
456b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
45749b5e25fSSatish Balay 
45849b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
459c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
46049b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
461b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
462b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
46349b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
46449b5e25fSSatish Balay   PetscFunctionReturn(0);
46549b5e25fSSatish Balay }
46649b5e25fSSatish Balay 
4674a2ae208SSatish Balay #undef __FUNCT__
4684a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
469dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
47049b5e25fSSatish Balay {
471dfbe8321SBarry Smith   PetscErrorCode ierr;
47232077d6dSBarry Smith   PetscTruth iascii,isdraw;
47349b5e25fSSatish Balay 
47449b5e25fSSatish Balay   PetscFunctionBegin;
47532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
476fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
47732077d6dSBarry Smith   if (iascii){
47849b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
47949b5e25fSSatish Balay   } else if (isdraw) {
48049b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
48149b5e25fSSatish Balay   } else {
482a5e6ed63SBarry Smith     Mat B;
483a5e6ed63SBarry Smith     ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr);
484a5e6ed63SBarry Smith     ierr = MatView(B,viewer);CHKERRQ(ierr);
485a5e6ed63SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
48649b5e25fSSatish Balay   }
48749b5e25fSSatish Balay   PetscFunctionReturn(0);
48849b5e25fSSatish Balay }
48949b5e25fSSatish Balay 
49049b5e25fSSatish Balay 
4914a2ae208SSatish Balay #undef __FUNCT__
4924a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
493dfbe8321SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[])
49449b5e25fSSatish Balay {
495045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
49649b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
49749b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
49849b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
49949b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
50049b5e25fSSatish Balay 
50149b5e25fSSatish Balay   PetscFunctionBegin;
50249b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
50349b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
504590ac198SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
505590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
50649b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
50749b5e25fSSatish Balay     nrow = ailen[brow];
50849b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
509590ac198SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
510590ac198SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1);
51149b5e25fSSatish Balay       col  = in[l] ;
51249b5e25fSSatish Balay       bcol = col/bs;
51349b5e25fSSatish Balay       cidx = col%bs;
51449b5e25fSSatish Balay       ridx = row%bs;
51549b5e25fSSatish Balay       high = nrow;
51649b5e25fSSatish Balay       low  = 0; /* assume unsorted */
51749b5e25fSSatish Balay       while (high-low > 5) {
51849b5e25fSSatish Balay         t = (low+high)/2;
51949b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
52049b5e25fSSatish Balay         else             low  = t;
52149b5e25fSSatish Balay       }
52249b5e25fSSatish Balay       for (i=low; i<high; i++) {
52349b5e25fSSatish Balay         if (rp[i] > bcol) break;
52449b5e25fSSatish Balay         if (rp[i] == bcol) {
52549b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
52649b5e25fSSatish Balay           goto finished;
52749b5e25fSSatish Balay         }
52849b5e25fSSatish Balay       }
52949b5e25fSSatish Balay       *v++ = zero;
53049b5e25fSSatish Balay       finished:;
53149b5e25fSSatish Balay     }
53249b5e25fSSatish Balay   }
53349b5e25fSSatish Balay   PetscFunctionReturn(0);
53449b5e25fSSatish Balay }
53549b5e25fSSatish Balay 
53649b5e25fSSatish Balay 
5374a2ae208SSatish Balay #undef __FUNCT__
5384a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
539dfbe8321SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
54049b5e25fSSatish Balay {
5410880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
5426849ba73SBarry Smith   PetscErrorCode ierr;
5430880e062SHong Zhang   int             *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
5440880e062SHong Zhang   int             *imax=a->imax,*ai=a->i,*ailen=a->ilen;
5456849ba73SBarry Smith   int             *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
5460880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
547f15d580aSBarry Smith   const MatScalar *value = v;
548f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
5490880e062SHong Zhang 
55049b5e25fSSatish Balay   PetscFunctionBegin;
5510880e062SHong Zhang   if (roworiented) {
5520880e062SHong Zhang     stepval = (n-1)*bs;
5530880e062SHong Zhang   } else {
5540880e062SHong Zhang     stepval = (m-1)*bs;
5550880e062SHong Zhang   }
5560880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
5570880e062SHong Zhang     row  = im[k];
5580880e062SHong Zhang     if (row < 0) continue;
5590880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
560590ac198SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1);
5610880e062SHong Zhang #endif
5620880e062SHong Zhang     rp   = aj + ai[row];
5630880e062SHong Zhang     ap   = aa + bs2*ai[row];
5640880e062SHong Zhang     rmax = imax[row];
5650880e062SHong Zhang     nrow = ailen[row];
5660880e062SHong Zhang     low  = 0;
5670880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
5680880e062SHong Zhang       if (in[l] < 0) continue;
5690880e062SHong Zhang       col = in[l];
570b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g)
571b1823623SSatish Balay       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",col,a->nbs-1);
572b1823623SSatish Balay #endif
5730880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
5740880e062SHong Zhang       if (roworiented) {
5750880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
5760880e062SHong Zhang       } else {
5770880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
5780880e062SHong Zhang       }
5790880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
5800880e062SHong Zhang       while (high-low > 7) {
5810880e062SHong Zhang         t = (low+high)/2;
5820880e062SHong Zhang         if (rp[t] > col) high = t;
5830880e062SHong Zhang         else             low  = t;
5840880e062SHong Zhang       }
5850880e062SHong Zhang       for (i=low; i<high; i++) {
5860880e062SHong Zhang         if (rp[i] > col) break;
5870880e062SHong Zhang         if (rp[i] == col) {
5880880e062SHong Zhang           bap  = ap +  bs2*i;
5890880e062SHong Zhang           if (roworiented) {
5900880e062SHong Zhang             if (is == ADD_VALUES) {
5910880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
5920880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
5930880e062SHong Zhang                   bap[jj] += *value++;
5940880e062SHong Zhang                 }
5950880e062SHong Zhang               }
5960880e062SHong Zhang             } else {
5970880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
5980880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
5990880e062SHong Zhang                   bap[jj] = *value++;
6000880e062SHong Zhang                 }
6010880e062SHong Zhang               }
6020880e062SHong Zhang             }
6030880e062SHong Zhang           } else {
6040880e062SHong Zhang             if (is == ADD_VALUES) {
6050880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6060880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6070880e062SHong Zhang                   *bap++ += *value++;
6080880e062SHong Zhang                 }
6090880e062SHong Zhang               }
6100880e062SHong Zhang             } else {
6110880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6120880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6130880e062SHong Zhang                   *bap++  = *value++;
6140880e062SHong Zhang                 }
6150880e062SHong Zhang               }
6160880e062SHong Zhang             }
6170880e062SHong Zhang           }
6180880e062SHong Zhang           goto noinsert2;
6190880e062SHong Zhang         }
6200880e062SHong Zhang       }
6210880e062SHong Zhang       if (nonew == 1) goto noinsert2;
622a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6230880e062SHong Zhang       if (nrow >= rmax) {
6240880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
6250880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
6260880e062SHong Zhang         MatScalar *new_a;
6270880e062SHong Zhang 
628a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6290880e062SHong Zhang 
6300880e062SHong Zhang         /* malloc new storage space */
6310880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
6320880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
6330880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
6340880e062SHong Zhang         new_i   = new_j + new_nz;
6350880e062SHong Zhang 
6360880e062SHong Zhang         /* copy over old data into new slots */
6370880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
6380880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
6390880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
6400880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
6410880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
6420880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
6430880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
6440880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
6450880e062SHong Zhang         /* free up old matrix storage */
6460880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
6470880e062SHong Zhang         if (!a->singlemalloc) {
6480880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
6490880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
6500880e062SHong Zhang         }
6510880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
6520880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
6530880e062SHong Zhang 
6540880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
6550880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
6560880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
6576c6c5352SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
6580880e062SHong Zhang         a->reallocs++;
6596c6c5352SBarry Smith         a->nz++;
6600880e062SHong Zhang       }
6610880e062SHong Zhang       N = nrow++ - 1;
6620880e062SHong Zhang       /* shift up all the later entries in this row */
6630880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
6640880e062SHong Zhang         rp[ii+1] = rp[ii];
6650880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
6660880e062SHong Zhang       }
6670880e062SHong Zhang       if (N >= i) {
6680880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6690880e062SHong Zhang       }
6700880e062SHong Zhang       rp[i] = col;
6710880e062SHong Zhang       bap   = ap +  bs2*i;
6720880e062SHong Zhang       if (roworiented) {
6730880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6740880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
6750880e062SHong Zhang             bap[jj] = *value++;
6760880e062SHong Zhang           }
6770880e062SHong Zhang         }
6780880e062SHong Zhang       } else {
6790880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6800880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
6810880e062SHong Zhang             *bap++  = *value++;
6820880e062SHong Zhang           }
6830880e062SHong Zhang         }
6840880e062SHong Zhang       }
6850880e062SHong Zhang       noinsert2:;
6860880e062SHong Zhang       low = i;
6870880e062SHong Zhang     }
6880880e062SHong Zhang     ailen[row] = nrow;
6890880e062SHong Zhang   }
6900880e062SHong Zhang   PetscFunctionReturn(0);
69149b5e25fSSatish Balay }
69249b5e25fSSatish Balay 
6934a2ae208SSatish Balay #undef __FUNCT__
6944a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
695dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
69649b5e25fSSatish Balay {
69749b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
6986849ba73SBarry Smith   PetscErrorCode ierr;
69949b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
700c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
7016849ba73SBarry Smith   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0;
70249b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
70349b5e25fSSatish Balay 
70449b5e25fSSatish Balay   PetscFunctionBegin;
70549b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
70649b5e25fSSatish Balay 
70749b5e25fSSatish Balay   if (m) rmax = ailen[0];
70849b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
70949b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
71049b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
71149b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
71249b5e25fSSatish Balay     if (fshift) {
71349b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
71449b5e25fSSatish Balay       N = ailen[i];
71549b5e25fSSatish Balay       for (j=0; j<N; j++) {
71649b5e25fSSatish Balay         ip[j-fshift] = ip[j];
71749b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
71849b5e25fSSatish Balay       }
71949b5e25fSSatish Balay     }
72049b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
72149b5e25fSSatish Balay   }
72249b5e25fSSatish Balay   if (mbs) {
72349b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
72449b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
72549b5e25fSSatish Balay   }
72649b5e25fSSatish Balay   /* reset ilen and imax for each row */
72749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
72849b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
72949b5e25fSSatish Balay   }
7306c6c5352SBarry Smith   a->nz = ai[mbs];
73149b5e25fSSatish Balay 
732b424e231SHong Zhang   /* diagonals may have moved, reset it */
733b424e231SHong Zhang   if (a->diag) {
734b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
73549b5e25fSSatish Balay   }
736b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
7376c6c5352SBarry Smith            m,A->m,a->bs,fshift*bs2,a->nz*bs2);
738b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
73949b5e25fSSatish Balay            a->reallocs);
740b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
74149b5e25fSSatish Balay   a->reallocs          = 0;
74249b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
74349b5e25fSSatish Balay 
74449b5e25fSSatish Balay   PetscFunctionReturn(0);
74549b5e25fSSatish Balay }
74649b5e25fSSatish Balay 
74749b5e25fSSatish Balay /*
74849b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
74949b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
75049b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
75149b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
75249b5e25fSSatish Balay */
7534a2ae208SSatish Balay #undef __FUNCT__
7544a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
755dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
75649b5e25fSSatish Balay {
75749b5e25fSSatish Balay   int        i,j,k,row;
75849b5e25fSSatish Balay   PetscTruth flg;
75949b5e25fSSatish Balay 
76049b5e25fSSatish Balay   PetscFunctionBegin;
76149b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
76249b5e25fSSatish Balay     row = idx[i];
76349b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
76449b5e25fSSatish Balay       sizes[j] = 1;
76549b5e25fSSatish Balay       i++;
76649b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
76749b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
76849b5e25fSSatish Balay       i++;
76949b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
77049b5e25fSSatish Balay       flg = PETSC_TRUE;
77149b5e25fSSatish Balay       for (k=1; k<bs; k++) {
77249b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
77349b5e25fSSatish Balay           flg = PETSC_FALSE;
77449b5e25fSSatish Balay           break;
77549b5e25fSSatish Balay         }
77649b5e25fSSatish Balay       }
77749b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
77849b5e25fSSatish Balay         sizes[j] = bs;
77949b5e25fSSatish Balay         i+= bs;
78049b5e25fSSatish Balay       } else {
78149b5e25fSSatish Balay         sizes[j] = 1;
78249b5e25fSSatish Balay         i++;
78349b5e25fSSatish Balay       }
78449b5e25fSSatish Balay     }
78549b5e25fSSatish Balay   }
78649b5e25fSSatish Balay   *bs_max = j;
78749b5e25fSSatish Balay   PetscFunctionReturn(0);
78849b5e25fSSatish Balay }
78949b5e25fSSatish Balay 
7904a2ae208SSatish Balay #undef __FUNCT__
7914a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
792dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag)
79349b5e25fSSatish Balay {
79449b5e25fSSatish Balay   PetscFunctionBegin;
795c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
79649b5e25fSSatish Balay }
79749b5e25fSSatish Balay 
79849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
79949b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
80049b5e25fSSatish Balay */
80149b5e25fSSatish Balay 
8024a2ae208SSatish Balay #undef __FUNCT__
8034a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
804dfbe8321SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
80549b5e25fSSatish Balay {
80649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
8076849ba73SBarry Smith   PetscErrorCode ierr;
80849b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
80949b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
81049b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
8116849ba73SBarry Smith   int         ridx,cidx,bs2=a->bs2;
81249b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
81349b5e25fSSatish Balay 
81449b5e25fSSatish Balay   PetscFunctionBegin;
81549b5e25fSSatish Balay 
81649b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
81749b5e25fSSatish Balay     row  = im[k];       /* row number */
81849b5e25fSSatish Balay     brow = row/bs;      /* block row number */
81949b5e25fSSatish Balay     if (row < 0) continue;
82049b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
821590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
82249b5e25fSSatish Balay #endif
82349b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
82449b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
82549b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
82649b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
82749b5e25fSSatish Balay     low  = 0;
82849b5e25fSSatish Balay 
82949b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
83049b5e25fSSatish Balay       if (in[l] < 0) continue;
83149b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
832590ac198SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m-1);
83349b5e25fSSatish Balay #endif
83449b5e25fSSatish Balay       col = in[l];
83549b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
83649b5e25fSSatish Balay 
83749b5e25fSSatish Balay       if (brow <= bcol){
83849b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8398549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
84049b5e25fSSatish Balay           /* element value a(k,l) */
84149b5e25fSSatish Balay           if (roworiented) {
84249b5e25fSSatish Balay             value = v[l + k*n];
84349b5e25fSSatish Balay           } else {
84449b5e25fSSatish Balay             value = v[k + l*m];
84549b5e25fSSatish Balay           }
84649b5e25fSSatish Balay 
84749b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
84849b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
84949b5e25fSSatish Balay           while (high-low > 7) {
85049b5e25fSSatish Balay             t = (low+high)/2;
85149b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
85249b5e25fSSatish Balay             else              low  = t;
85349b5e25fSSatish Balay           }
85449b5e25fSSatish Balay           for (i=low; i<high; i++) {
85549b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
85649b5e25fSSatish Balay             if (rp[i] > bcol) break;
85749b5e25fSSatish Balay             if (rp[i] == bcol) {
85849b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
85949b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
86049b5e25fSSatish Balay               else                  *bap  = value;
8618549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8628549e402SHong Zhang               if (brow == bcol && ridx < cidx){
8638549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
8648549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
8658549e402SHong Zhang                 else                  *bap  = value;
8668549e402SHong Zhang               }
86749b5e25fSSatish Balay               goto noinsert1;
86849b5e25fSSatish Balay             }
86949b5e25fSSatish Balay           }
87049b5e25fSSatish Balay 
87149b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
872a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
87349b5e25fSSatish Balay       if (nrow >= rmax) {
87449b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
87549b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
87649b5e25fSSatish Balay         MatScalar *new_a;
87749b5e25fSSatish Balay 
878a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
87949b5e25fSSatish Balay 
88049b5e25fSSatish Balay         /* Malloc new storage space */
88149b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
88282502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
88349b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
88449b5e25fSSatish Balay         new_i = new_j + new_nz;
88549b5e25fSSatish Balay 
88649b5e25fSSatish Balay         /* copy over old data into new slots */
88749b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
88849b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
88949b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
89049b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
89149b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
89249b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
89349b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
89449b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
89549b5e25fSSatish Balay         /* free up old matrix storage */
89649b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
89749b5e25fSSatish Balay         if (!a->singlemalloc) {
89849b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
89949b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
90049b5e25fSSatish Balay         }
90149b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
90249b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
90349b5e25fSSatish Balay 
90449b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
90549b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
906b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
9076c6c5352SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
90849b5e25fSSatish Balay         a->reallocs++;
9096c6c5352SBarry Smith         a->nz++;
91049b5e25fSSatish Balay       }
91149b5e25fSSatish Balay 
91249b5e25fSSatish Balay       N = nrow++ - 1;
91349b5e25fSSatish Balay       /* shift up all the later entries in this row */
91449b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
91549b5e25fSSatish Balay         rp[ii+1] = rp[ii];
91649b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
91749b5e25fSSatish Balay       }
91849b5e25fSSatish Balay       if (N>=i) {
91949b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
92049b5e25fSSatish Balay       }
92149b5e25fSSatish Balay       rp[i]                      = bcol;
92249b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
92349b5e25fSSatish Balay       noinsert1:;
92449b5e25fSSatish Balay       low = i;
92549b5e25fSSatish Balay       /* } */
9268549e402SHong Zhang         }
92749b5e25fSSatish Balay       } /* end of if .. if..  */
92849b5e25fSSatish Balay     }                     /* end of loop over added columns */
92949b5e25fSSatish Balay     ailen[brow] = nrow;
93049b5e25fSSatish Balay   }                       /* end of loop over added rows */
93149b5e25fSSatish Balay 
93249b5e25fSSatish Balay   PetscFunctionReturn(0);
93349b5e25fSSatish Balay }
93449b5e25fSSatish Balay 
9356849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
9366849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
9376849ba73SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int);
9386849ba73SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
9396849ba73SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
9406849ba73SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
9416849ba73SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
9426849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat);
9436849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
9446849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
9456849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
9466849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
9476849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
9486849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
9496849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec);
95049b5e25fSSatish Balay 
9516849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
9526849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
9536849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
9546849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
9556849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
9566849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
9576849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
9586849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
9596849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
9606849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
9616849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
9626849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
9636849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
9646849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
9656849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
96649b5e25fSSatish Balay 
9676849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
968d59c15a7SBarry Smith 
9696849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
9706849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
9716849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
9726849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
9736849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
9746849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
9756849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
9766849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
97707f98182SSatish Balay 
9786849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
9796849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
9806849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
9816849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
9826849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
9836849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
9846849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
9856849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
9866849ba73SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
98749b5e25fSSatish Balay 
9886849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
9896849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
9906849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
9916849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
9926849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
9936849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
9946849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
9956849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
99649b5e25fSSatish Balay 
9976849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
9986849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
9996849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
10006849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
10016849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
10026849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
10036849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
10046849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
100549b5e25fSSatish Balay 
10064d101231SSatish Balay #ifdef HAVE_MatICCFactor
10076849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work!  */
10084a2ae208SSatish Balay #undef __FUNCT__
10094d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
1010dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
101149b5e25fSSatish Balay {
10124ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
101349b5e25fSSatish Balay   Mat         outA;
1014dfbe8321SBarry Smith   PetscErrorCode ierr;
101549b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
101649b5e25fSSatish Balay 
101749b5e25fSSatish Balay   PetscFunctionBegin;
10181a3463dfSHong Zhang   /*
101929bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
102049b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
102149b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
102249b5e25fSSatish Balay   if (!row_identity || !col_identity) {
102329bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
102449b5e25fSSatish Balay   }
10251a3463dfSHong Zhang   */
102649b5e25fSSatish Balay 
102749b5e25fSSatish Balay   outA          = inA;
10281260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
102949b5e25fSSatish Balay 
103049b5e25fSSatish Balay   if (!a->diag) {
10311a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
103249b5e25fSSatish Balay   }
103349b5e25fSSatish Balay   /*
103449b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
103549b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
103649b5e25fSSatish Balay   */
103749b5e25fSSatish Balay   switch (a->bs) {
103849b5e25fSSatish Balay   case 1:
10390fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
104007f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1041d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
10424d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
104349b5e25fSSatish Balay   case 2:
10441a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
10451a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10461a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10474d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
104849b5e25fSSatish Balay     break;
104949b5e25fSSatish Balay   case 3:
10501a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
10511a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10521a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10534d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
105449b5e25fSSatish Balay     break;
105549b5e25fSSatish Balay   case 4:
10561a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
10571a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10581a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10594d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
106049b5e25fSSatish Balay     break;
106149b5e25fSSatish Balay   case 5:
10621a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
10631a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10641a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10654d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
106649b5e25fSSatish Balay     break;
106749b5e25fSSatish Balay   case 6:
10681a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
10691a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10701a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10714d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
107249b5e25fSSatish Balay     break;
107349b5e25fSSatish Balay   case 7:
10741a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
10751a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10761a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10774d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
107849b5e25fSSatish Balay     break;
107949b5e25fSSatish Balay   default:
108049b5e25fSSatish Balay     a->row        = row;
10811a3463dfSHong Zhang     a->icol       = col;
108249b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
108349b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
108449b5e25fSSatish Balay 
108549b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
108649b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1087b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
108849b5e25fSSatish Balay 
108949b5e25fSSatish Balay     if (!a->solve_work) {
109087828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
109187828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
109249b5e25fSSatish Balay     }
109349b5e25fSSatish Balay   }
109449b5e25fSSatish Balay 
10951a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
109649b5e25fSSatish Balay 
109749b5e25fSSatish Balay   PetscFunctionReturn(0);
109849b5e25fSSatish Balay }
1099950f1e5bSHong Zhang #endif
1100950f1e5bSHong Zhang 
11014a2ae208SSatish Balay #undef __FUNCT__
11024a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
1103dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A)
110449b5e25fSSatish Balay {
110549b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
110649b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
1107dfbe8321SBarry Smith   PetscErrorCode ierr;
110849b5e25fSSatish Balay 
110949b5e25fSSatish Balay   PetscFunctionBegin;
111049b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
111149b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
111249b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
111349b5e25fSSatish Balay   PetscFunctionReturn(0);
111449b5e25fSSatish Balay }
111549b5e25fSSatish Balay 
111649b5e25fSSatish Balay EXTERN_C_BEGIN
11174a2ae208SSatish Balay #undef __FUNCT__
11184a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
1119dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
112049b5e25fSSatish Balay {
1121045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
112249b5e25fSSatish Balay   int         i,nz,n;
112349b5e25fSSatish Balay 
112449b5e25fSSatish Balay   PetscFunctionBegin;
11256c6c5352SBarry Smith   nz = baij->maxnz;
1126c464158bSHong Zhang   n  = mat->n;
112749b5e25fSSatish Balay   for (i=0; i<nz; i++) {
112849b5e25fSSatish Balay     baij->j[i] = indices[i];
112949b5e25fSSatish Balay   }
11306c6c5352SBarry Smith   baij->nz = nz;
113149b5e25fSSatish Balay   for (i=0; i<n; i++) {
113249b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
113349b5e25fSSatish Balay   }
113449b5e25fSSatish Balay 
113549b5e25fSSatish Balay   PetscFunctionReturn(0);
113649b5e25fSSatish Balay }
113749b5e25fSSatish Balay EXTERN_C_END
113849b5e25fSSatish Balay 
11394a2ae208SSatish Balay #undef __FUNCT__
11404a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
114149b5e25fSSatish Balay /*@
114219585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
114349b5e25fSSatish Balay        in the matrix.
114449b5e25fSSatish Balay 
114549b5e25fSSatish Balay   Input Parameters:
114619585528SSatish Balay +  mat     - the SeqSBAIJ matrix
114749b5e25fSSatish Balay -  indices - the column indices
114849b5e25fSSatish Balay 
114949b5e25fSSatish Balay   Level: advanced
115049b5e25fSSatish Balay 
115149b5e25fSSatish Balay   Notes:
115249b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
115349b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
115449b5e25fSSatish Balay   of the MatSetValues() operation.
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
115719585528SSatish Balay   MatCreateSeqSBAIJ().
115849b5e25fSSatish Balay 
1159ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
116049b5e25fSSatish Balay 
1161ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
116249b5e25fSSatish Balay @*/
1163dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
116449b5e25fSSatish Balay {
1165dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,int *);
116649b5e25fSSatish Balay 
116749b5e25fSSatish Balay   PetscFunctionBegin;
11684482741eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
11694482741eSBarry Smith   PetscValidPointer(indices,2);
1170c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
117149b5e25fSSatish Balay   if (f) {
117249b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
117349b5e25fSSatish Balay   } else {
117429bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
117549b5e25fSSatish Balay   }
117649b5e25fSSatish Balay   PetscFunctionReturn(0);
117749b5e25fSSatish Balay }
117849b5e25fSSatish Balay 
11794a2ae208SSatish Balay #undef __FUNCT__
11804a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1181dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A)
1182273d9f13SBarry Smith {
1183dfbe8321SBarry Smith   PetscErrorCode ierr;
1184273d9f13SBarry Smith 
1185273d9f13SBarry Smith   PetscFunctionBegin;
1186273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1187273d9f13SBarry Smith   PetscFunctionReturn(0);
1188273d9f13SBarry Smith }
1189273d9f13SBarry Smith 
1190a6ece127SHong Zhang #undef __FUNCT__
1191a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1192dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1193a6ece127SHong Zhang {
1194a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1195a6ece127SHong Zhang   PetscFunctionBegin;
1196a6ece127SHong Zhang   *array = a->a;
1197a6ece127SHong Zhang   PetscFunctionReturn(0);
1198a6ece127SHong Zhang }
1199a6ece127SHong Zhang 
1200a6ece127SHong Zhang #undef __FUNCT__
1201a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1202dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1203a6ece127SHong Zhang {
1204a6ece127SHong Zhang   PetscFunctionBegin;
1205a6ece127SHong Zhang   PetscFunctionReturn(0);
1206a6ece127SHong Zhang }
1207a6ece127SHong Zhang 
120842ee4b1aSHong Zhang #include "petscblaslapack.h"
120942ee4b1aSHong Zhang #undef __FUNCT__
121042ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1211dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
121242ee4b1aSHong Zhang {
121342ee4b1aSHong Zhang   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1214dfbe8321SBarry Smith   PetscErrorCode ierr;
1215*4ce68768SBarry Smith   int            i,bs=y->bs,bs2,j;
1216*4ce68768SBarry Smith   PetscBLASInt   bnz = (PetscBLASInt)x->nz,one = 1;
121742ee4b1aSHong Zhang 
121842ee4b1aSHong Zhang   PetscFunctionBegin;
121942ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1220*4ce68768SBarry Smith     BLaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one);
1221c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1222c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1223c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1224c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1225c537a176SHong Zhang     }
1226c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1227c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1228c4319e64SHong Zhang       y->XtoY = X;
1229c537a176SHong Zhang     }
1230c4319e64SHong Zhang     bs2 = bs*bs;
12316c6c5352SBarry Smith     for (i=0; i<x->nz; i++) {
1232c4319e64SHong Zhang       j = 0;
1233c4319e64SHong Zhang       while (j < bs2){
12346550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1235c4319e64SHong Zhang         j++;
1236c537a176SHong Zhang       }
1237c4319e64SHong Zhang     }
12386c6c5352SBarry Smith     PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %d/%d = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));
123942ee4b1aSHong Zhang   } else {
124042ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
124142ee4b1aSHong Zhang   }
124242ee4b1aSHong Zhang   PetscFunctionReturn(0);
124342ee4b1aSHong Zhang }
124442ee4b1aSHong Zhang 
1245efcf0fc3SBarry Smith #undef __FUNCT__
1246efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1247dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg)
1248efcf0fc3SBarry Smith {
1249efcf0fc3SBarry Smith   PetscFunctionBegin;
1250efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1251efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1252efcf0fc3SBarry Smith }
1253efcf0fc3SBarry Smith 
1254efcf0fc3SBarry Smith #undef __FUNCT__
1255efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1256dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1257efcf0fc3SBarry Smith {
1258efcf0fc3SBarry Smith   PetscFunctionBegin;
1259efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1260efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1261efcf0fc3SBarry Smith }
1262efcf0fc3SBarry Smith 
1263efcf0fc3SBarry Smith #undef __FUNCT__
1264efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1265dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg)
1266efcf0fc3SBarry Smith {
1267efcf0fc3SBarry Smith   PetscFunctionBegin;
1268efcf0fc3SBarry Smith   *flg = PETSC_FALSE;
1269efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1270efcf0fc3SBarry Smith }
1271efcf0fc3SBarry Smith 
127249b5e25fSSatish Balay /* -------------------------------------------------------------------*/
127349b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
127449b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
127549b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
127649b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
127797304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
127849b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
127949b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
128049b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
128149b5e25fSSatish Balay        0,
128249b5e25fSSatish Balay        0,
128397304618SKris Buschelman /*10*/ 0,
128449b5e25fSSatish Balay        0,
12855f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1286d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
128749b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
128897304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
128949b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
129049b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
129149b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
129249b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
129397304618SKris Buschelman /*20*/ 0,
129449b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
129549b5e25fSSatish Balay        0,
129649b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
129749b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
129897304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ,
129949b5e25fSSatish Balay        0,
130049b5e25fSSatish Balay        0,
13015f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13025f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
130397304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1304c464158bSHong Zhang        0,
13054d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1306a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1307a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
130897304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
130949b5e25fSSatish Balay        0,
131049b5e25fSSatish Balay        0,
131149b5e25fSSatish Balay        0,
1312950f1e5bSHong Zhang        0,
131397304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
131449b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
131549b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
131649b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
131749b5e25fSSatish Balay        0,
131897304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
131949b5e25fSSatish Balay        MatScale_SeqSBAIJ,
132049b5e25fSSatish Balay        0,
132149b5e25fSSatish Balay        0,
132249b5e25fSSatish Balay        0,
132397304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ,
132449b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
132549b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
132649b5e25fSSatish Balay        0,
132749b5e25fSSatish Balay        0,
132897304618SKris Buschelman /*55*/ 0,
132949b5e25fSSatish Balay        0,
133049b5e25fSSatish Balay        0,
133149b5e25fSSatish Balay        0,
133249b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
133397304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
133449b5e25fSSatish Balay        0,
133549b5e25fSSatish Balay        0,
13368a124369SBarry Smith        MatGetPetscMaps_Petsc,
1337d959ec07SHong Zhang        0,
133897304618SKris Buschelman /*65*/ 0,
1339d959ec07SHong Zhang        0,
1340d959ec07SHong Zhang        0,
1341d959ec07SHong Zhang        0,
1342d959ec07SHong Zhang        0,
134397304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
13443e0d88b5SBarry Smith        0,
13453e0d88b5SBarry Smith        0,
13463e0d88b5SBarry Smith        0,
13473e0d88b5SBarry Smith        0,
134897304618SKris Buschelman /*75*/ 0,
13493e0d88b5SBarry Smith        0,
13503e0d88b5SBarry Smith        0,
13513e0d88b5SBarry Smith        0,
13523e0d88b5SBarry Smith        0,
135397304618SKris Buschelman /*80*/ 0,
13543e0d88b5SBarry Smith        0,
13553e0d88b5SBarry Smith        0,
13563e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
135797304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
13583e0d88b5SBarry Smith #else
135997304618SKris Buschelman        0,
13603e0d88b5SBarry Smith #endif
136141acf15aSKris Buschelman /*84*/ MatLoad_SeqSBAIJ,
1362efcf0fc3SBarry Smith        MatIsSymmetric_SeqSBAIJ,
1363efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1364efcf0fc3SBarry Smith        MatIsHermitian_SeqSBAIJ
13653e0d88b5SBarry Smith };
136649b5e25fSSatish Balay 
136749b5e25fSSatish Balay EXTERN_C_BEGIN
13684a2ae208SSatish Balay #undef __FUNCT__
13694a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
1370dfbe8321SBarry Smith PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
137149b5e25fSSatish Balay {
13724afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1373c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1374dfbe8321SBarry Smith   PetscErrorCode ierr;
137549b5e25fSSatish Balay 
137649b5e25fSSatish Balay   PetscFunctionBegin;
137749b5e25fSSatish Balay   if (aij->nonew != 1) {
137829bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
137949b5e25fSSatish Balay   }
138049b5e25fSSatish Balay 
138149b5e25fSSatish Balay   /* allocate space for values if not already there */
138249b5e25fSSatish Balay   if (!aij->saved_values) {
138387828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
138449b5e25fSSatish Balay   }
138549b5e25fSSatish Balay 
138649b5e25fSSatish Balay   /* copy values over */
138787828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
138849b5e25fSSatish Balay   PetscFunctionReturn(0);
138949b5e25fSSatish Balay }
139049b5e25fSSatish Balay EXTERN_C_END
139149b5e25fSSatish Balay 
139249b5e25fSSatish Balay EXTERN_C_BEGIN
13934a2ae208SSatish Balay #undef __FUNCT__
13944a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
1395dfbe8321SBarry Smith PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
139649b5e25fSSatish Balay {
13974afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
13986849ba73SBarry Smith   PetscErrorCode ierr;
13996849ba73SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
140049b5e25fSSatish Balay 
140149b5e25fSSatish Balay   PetscFunctionBegin;
140249b5e25fSSatish Balay   if (aij->nonew != 1) {
140329bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
140449b5e25fSSatish Balay   }
140549b5e25fSSatish Balay   if (!aij->saved_values) {
140629bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
140749b5e25fSSatish Balay   }
140849b5e25fSSatish Balay 
140949b5e25fSSatish Balay   /* copy values over */
141087828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
141149b5e25fSSatish Balay   PetscFunctionReturn(0);
141249b5e25fSSatish Balay }
141349b5e25fSSatish Balay EXTERN_C_END
141449b5e25fSSatish Balay 
14158549e402SHong Zhang EXTERN_C_BEGIN
14164a2ae208SSatish Balay #undef __FUNCT__
1417a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1418dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz)
141949b5e25fSSatish Balay {
1420c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
14216849ba73SBarry Smith   PetscErrorCode ierr;
14226849ba73SBarry Smith   int          i,len,mbs,bs2;
142349b5e25fSSatish Balay   PetscTruth   flg;
142449b5e25fSSatish Balay 
142549b5e25fSSatish Balay   PetscFunctionBegin;
1426273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1427e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1428c464158bSHong Zhang   mbs  = B->m/bs;
142949b5e25fSSatish Balay   bs2  = bs*bs;
143049b5e25fSSatish Balay 
1431c464158bSHong Zhang   if (mbs*bs != B->m) {
143229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
143349b5e25fSSatish Balay   }
143449b5e25fSSatish Balay 
1435435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1436435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
143749b5e25fSSatish Balay   if (nnz) {
143849b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
143929bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
144080fe4e49SBarry Smith       if (nnz[i] > mbs) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %d value %d rowlength %d",i,nnz[i],mbs);
144149b5e25fSSatish Balay     }
144249b5e25fSSatish Balay   }
144349b5e25fSSatish Balay 
1444e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
144549b5e25fSSatish Balay   if (!flg) {
144649b5e25fSSatish Balay     switch (bs) {
144749b5e25fSSatish Balay     case 1:
14484ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
144949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1450d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
145149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
145249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
145349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
145449b5e25fSSatish Balay       break;
145549b5e25fSSatish Balay     case 2:
14564ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
145749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
145849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
145949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
146049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
146149b5e25fSSatish Balay       break;
146249b5e25fSSatish Balay     case 3:
14635f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
146449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
146549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
146649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
146749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
146849b5e25fSSatish Balay       break;
146949b5e25fSSatish Balay     case 4:
14705f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
147149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
147249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
147349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
147449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
147549b5e25fSSatish Balay       break;
147649b5e25fSSatish Balay     case 5:
14775f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
147849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
147949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
148049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
148149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
148249b5e25fSSatish Balay       break;
148349b5e25fSSatish Balay     case 6:
14845f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
148549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
148649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
148749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
148849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
148949b5e25fSSatish Balay       break;
149049b5e25fSSatish Balay     case 7:
1491de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
149249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
149349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1494de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
149549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
149649b5e25fSSatish Balay       break;
149749b5e25fSSatish Balay     }
149849b5e25fSSatish Balay   }
149949b5e25fSSatish Balay 
150049b5e25fSSatish Balay   b->mbs = mbs;
15014afc71dfSHong Zhang   b->nbs = mbs;
1502b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
150349b5e25fSSatish Balay   if (!nnz) {
1504435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
150549b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
150649b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
15078cef66ccSHong Zhang       b->imax[i] = nz;
150849b5e25fSSatish Balay     }
1509153ea458SHong Zhang     nz = nz*mbs; /* total nz */
151049b5e25fSSatish Balay   } else {
151149b5e25fSSatish Balay     nz = 0;
15128cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
151349b5e25fSSatish Balay   }
15146c6c5352SBarry Smith   /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
151549b5e25fSSatish Balay 
151649b5e25fSSatish Balay   /* allocate the matrix space */
15176c6c5352SBarry Smith   len  = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
151882502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
15196c6c5352SBarry Smith   ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
15206c6c5352SBarry Smith   b->j = (int*)(b->a + nz*bs2);
15216c6c5352SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
15226c6c5352SBarry Smith   b->i = b->j + nz;
152349b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
152449b5e25fSSatish Balay 
152549b5e25fSSatish Balay   /* pointer to beginning of each row */
152649b5e25fSSatish Balay   b->i[0] = 0;
152749b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
152849b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
152949b5e25fSSatish Balay   }
153049b5e25fSSatish Balay 
153149b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
15325033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1533b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
153449b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
153549b5e25fSSatish Balay 
153649b5e25fSSatish Balay   b->bs               = bs;
153749b5e25fSSatish Balay   b->bs2              = bs2;
15386c6c5352SBarry Smith   b->nz             = 0;
15396c6c5352SBarry Smith   b->maxnz          = nz*bs2;
1540153ea458SHong Zhang 
154116cdd363SHong Zhang   b->inew             = 0;
154216cdd363SHong Zhang   b->jnew             = 0;
154316cdd363SHong Zhang   b->anew             = 0;
154416cdd363SHong Zhang   b->a2anew           = 0;
15451a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1546c464158bSHong Zhang   PetscFunctionReturn(0);
1547c464158bSHong Zhang }
1548a23d5eceSKris Buschelman EXTERN_C_END
1549153ea458SHong Zhang 
1550d769727bSBarry Smith EXTERN_C_BEGIN
1551dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,Mat*);
1552dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,Mat*);
1553d769727bSBarry Smith EXTERN_C_END
1554d769727bSBarry Smith 
15550bad9183SKris Buschelman /*MC
1556fafad747SKris Buschelman    MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
15570bad9183SKris Buschelman    based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
15580bad9183SKris Buschelman 
15590bad9183SKris Buschelman    Options Database Keys:
15600bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
15610bad9183SKris Buschelman 
15620bad9183SKris Buschelman   Level: beginner
15630bad9183SKris Buschelman 
15640bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ
15650bad9183SKris Buschelman M*/
15660bad9183SKris Buschelman 
1567a23d5eceSKris Buschelman EXTERN_C_BEGIN
1568a23d5eceSKris Buschelman #undef __FUNCT__
1569a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1570dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1571a23d5eceSKris Buschelman {
1572a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
1573dfbe8321SBarry Smith   PetscErrorCode ierr;
1574dfbe8321SBarry Smith   int size;
1575a23d5eceSKris Buschelman 
1576a23d5eceSKris Buschelman   PetscFunctionBegin;
1577a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1578a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1579a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1580a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1581a23d5eceSKris Buschelman 
1582a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1583a23d5eceSKris Buschelman   B->data = (void*)b;
1584a23d5eceSKris Buschelman   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1585a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1586a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1587a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1588a23d5eceSKris Buschelman   B->factor           = 0;
1589a23d5eceSKris Buschelman   B->lupivotthreshold = 1.0;
1590a23d5eceSKris Buschelman   B->mapping          = 0;
1591a23d5eceSKris Buschelman   b->row              = 0;
1592a23d5eceSKris Buschelman   b->icol             = 0;
1593a23d5eceSKris Buschelman   b->reallocs         = 0;
1594a23d5eceSKris Buschelman   b->saved_values     = 0;
1595a23d5eceSKris Buschelman 
1596a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1597a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1598a23d5eceSKris Buschelman 
1599a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1600a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1601a23d5eceSKris Buschelman   b->nonew            = 0;
1602a23d5eceSKris Buschelman   b->diag             = 0;
1603a23d5eceSKris Buschelman   b->solve_work       = 0;
1604a23d5eceSKris Buschelman   b->mult_work        = 0;
1605a23d5eceSKris Buschelman   B->spptr            = 0;
1606a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1607a23d5eceSKris Buschelman   b->xtoy             = 0;
1608a23d5eceSKris Buschelman   b->XtoY             = 0;
1609a23d5eceSKris Buschelman 
1610a23d5eceSKris Buschelman   b->inew             = 0;
1611a23d5eceSKris Buschelman   b->jnew             = 0;
1612a23d5eceSKris Buschelman   b->anew             = 0;
1613a23d5eceSKris Buschelman   b->a2anew           = 0;
1614a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1615a23d5eceSKris Buschelman 
1616a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1617a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1618a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1619a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1620a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1621a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1622a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1623a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1624a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
16254e5e7fe4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",
16264e5e7fe4SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqAIJ",
16274e5e7fe4SHong Zhang                                       MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1628a0e1a404SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",
1629a0e1a404SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqBAIJ",
1630a0e1a404SHong Zhang                                       MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1631a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1632a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1633a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
163423ce1328SBarry Smith 
163523ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
163623ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
163723ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
163823ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
1639a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1640a23d5eceSKris Buschelman }
1641a23d5eceSKris Buschelman EXTERN_C_END
1642a23d5eceSKris Buschelman 
1643a23d5eceSKris Buschelman #undef __FUNCT__
1644a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1645a23d5eceSKris Buschelman /*@C
1646a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1647a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1648a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1649a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1650a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1651a23d5eceSKris Buschelman 
1652a23d5eceSKris Buschelman    Collective on Mat
1653a23d5eceSKris Buschelman 
1654a23d5eceSKris Buschelman    Input Parameters:
1655a23d5eceSKris Buschelman +  A - the symmetric matrix
1656a23d5eceSKris Buschelman .  bs - size of block
1657a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1658a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1659a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1660a23d5eceSKris Buschelman 
1661a23d5eceSKris Buschelman    Options Database Keys:
1662a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1663a23d5eceSKris Buschelman                      block calculations (much slower)
1664a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1665a23d5eceSKris Buschelman 
1666a23d5eceSKris Buschelman    Level: intermediate
1667a23d5eceSKris Buschelman 
1668a23d5eceSKris Buschelman    Notes:
1669a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1670a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1671a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1672a23d5eceSKris Buschelman    matrices.
1673a23d5eceSKris Buschelman 
1674a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1675a23d5eceSKris Buschelman @*/
1676dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) {
1677dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,int,int,const int[]);
1678a23d5eceSKris Buschelman 
1679a23d5eceSKris Buschelman   PetscFunctionBegin;
1680a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1681a23d5eceSKris Buschelman   if (f) {
1682a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1683a23d5eceSKris Buschelman   }
1684a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1685a23d5eceSKris Buschelman }
168649b5e25fSSatish Balay 
16874a2ae208SSatish Balay #undef __FUNCT__
16884a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1689c464158bSHong Zhang /*@C
1690c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1691c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1692c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1693c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1694c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
169549b5e25fSSatish Balay 
1696c464158bSHong Zhang    Collective on MPI_Comm
1697c464158bSHong Zhang 
1698c464158bSHong Zhang    Input Parameters:
1699c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1700c464158bSHong Zhang .  bs - size of block
1701c464158bSHong Zhang .  m - number of rows, or number of columns
1702c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1703744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1704744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1705c464158bSHong Zhang 
1706c464158bSHong Zhang    Output Parameter:
1707c464158bSHong Zhang .  A - the symmetric matrix
1708c464158bSHong Zhang 
1709c464158bSHong Zhang    Options Database Keys:
1710c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1711c464158bSHong Zhang                      block calculations (much slower)
1712c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1713c464158bSHong Zhang 
1714c464158bSHong Zhang    Level: intermediate
1715c464158bSHong Zhang 
1716c464158bSHong Zhang    Notes:
1717c464158bSHong Zhang 
1718c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1719c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1720c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1721c464158bSHong Zhang    matrices.
1722c464158bSHong Zhang 
1723c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1724c464158bSHong Zhang @*/
1725dfbe8321SBarry Smith PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A)
1726c464158bSHong Zhang {
1727dfbe8321SBarry Smith   PetscErrorCode ierr;
1728c464158bSHong Zhang 
1729c464158bSHong Zhang   PetscFunctionBegin;
1730c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1731c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1732273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
173349b5e25fSSatish Balay   PetscFunctionReturn(0);
173449b5e25fSSatish Balay }
173549b5e25fSSatish Balay 
17364a2ae208SSatish Balay #undef __FUNCT__
17374a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
1738dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
173949b5e25fSSatish Balay {
174049b5e25fSSatish Balay   Mat          C;
174149b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
17426849ba73SBarry Smith   PetscErrorCode ierr;
17436849ba73SBarry Smith   int          i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
174449b5e25fSSatish Balay 
174549b5e25fSSatish Balay   PetscFunctionBegin;
174629bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
174749b5e25fSSatish Balay 
174849b5e25fSSatish Balay   *B = 0;
1749692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1750be5d1d56SKris Buschelman   ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
17511d5dac46SHong Zhang   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1752692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1753692f9cbeSHong Zhang 
1754273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
175549b5e25fSSatish Balay   C->factor         = A->factor;
175649b5e25fSSatish Balay   c->row            = 0;
175749b5e25fSSatish Balay   c->icol           = 0;
175849b5e25fSSatish Balay   c->saved_values   = 0;
175949b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
176049b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
176149b5e25fSSatish Balay 
176282327fa8SHong Zhang   C->M    = A->M;
176382327fa8SHong Zhang   C->N    = A->N;
176449b5e25fSSatish Balay   c->bs   = a->bs;
176549b5e25fSSatish Balay   c->bs2  = a->bs2;
176649b5e25fSSatish Balay   c->mbs  = a->mbs;
176749b5e25fSSatish Balay   c->nbs  = a->nbs;
176849b5e25fSSatish Balay 
1769b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1770b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
177149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
177249b5e25fSSatish Balay     c->imax[i] = a->imax[i];
177349b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
177449b5e25fSSatish Balay   }
177549b5e25fSSatish Balay 
177649b5e25fSSatish Balay   /* allocate the matrix space */
177749b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
177849b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
177982502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
178049b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
178149b5e25fSSatish Balay   c->i = c->j + nz;
178249b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
178349b5e25fSSatish Balay   if (mbs > 0) {
178449b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
178549b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
178649b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178749b5e25fSSatish Balay     } else {
178849b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178949b5e25fSSatish Balay     }
179049b5e25fSSatish Balay   }
179149b5e25fSSatish Balay 
1792b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
179349b5e25fSSatish Balay   c->sorted      = a->sorted;
179449b5e25fSSatish Balay   c->roworiented = a->roworiented;
179549b5e25fSSatish Balay   c->nonew       = a->nonew;
179649b5e25fSSatish Balay 
179749b5e25fSSatish Balay   if (a->diag) {
1798b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1799b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
180049b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
180149b5e25fSSatish Balay       c->diag[i] = a->diag[i];
180249b5e25fSSatish Balay     }
180349b5e25fSSatish Balay   } else c->diag        = 0;
18046c6c5352SBarry Smith   c->nz               = a->nz;
18056c6c5352SBarry Smith   c->maxnz            = a->maxnz;
180649b5e25fSSatish Balay   c->solve_work         = 0;
180749b5e25fSSatish Balay   c->mult_work          = 0;
180849b5e25fSSatish Balay   *B = C;
1809b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
181049b5e25fSSatish Balay   PetscFunctionReturn(0);
181149b5e25fSSatish Balay }
181249b5e25fSSatish Balay 
18134a2ae208SSatish Balay #undef __FUNCT__
18144a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1815dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A)
181649b5e25fSSatish Balay {
181749b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
181849b5e25fSSatish Balay   Mat          B;
18196849ba73SBarry Smith   PetscErrorCode ierr;
18206849ba73SBarry Smith   int          i,nz,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18213f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
182249b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
182349b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
182487828ca2SBarry Smith   PetscScalar  *aa;
182549b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
182649b5e25fSSatish Balay 
182749b5e25fSSatish Balay   PetscFunctionBegin;
1828b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
182949b5e25fSSatish Balay   bs2  = bs*bs;
183049b5e25fSSatish Balay 
183149b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
183229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1833b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
183449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1835552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
183649b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
183749b5e25fSSatish Balay 
183849b5e25fSSatish Balay   if (header[3] < 0) {
183929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
184049b5e25fSSatish Balay   }
184149b5e25fSSatish Balay 
184229bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
184349b5e25fSSatish Balay 
184449b5e25fSSatish Balay   /*
184549b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
184649b5e25fSSatish Balay     divisible by the blocksize
184749b5e25fSSatish Balay   */
184849b5e25fSSatish Balay   mbs        = M/bs;
184949b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
185049b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
185149b5e25fSSatish Balay   else                  mbs++;
185249b5e25fSSatish Balay   if (extra_rows) {
1853b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
185449b5e25fSSatish Balay   }
185549b5e25fSSatish Balay 
185649b5e25fSSatish Balay   /* read in row lengths */
1857b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
185849b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
185949b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
186049b5e25fSSatish Balay 
186149b5e25fSSatish Balay   /* read in column indices */
1862b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
186349b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
186449b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
186549b5e25fSSatish Balay 
186649b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
186782502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1868a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
186982502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
187049b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
187149b5e25fSSatish Balay   masked   = mask + mbs;
187249b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
187349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
187449b5e25fSSatish Balay     nmask = 0;
187549b5e25fSSatish Balay     for (j=0; j<bs; j++) {
187649b5e25fSSatish Balay       kmax = rowlengths[rowcount];
187749b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18782d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
187903630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
188049b5e25fSSatish Balay       }
188149b5e25fSSatish Balay       rowcount++;
188249b5e25fSSatish Balay     }
1883574b2666SHong Zhang     s_browlengths[i] += nmask;
1884574b2666SHong Zhang 
188549b5e25fSSatish Balay     /* zero out the mask elements we set */
188649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
188749b5e25fSSatish Balay   }
188849b5e25fSSatish Balay 
188949b5e25fSSatish Balay   /* create our matrix */
18909abb65ffSKris Buschelman   ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr);
18919abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
189278473fd7SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr);
189349b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
189449b5e25fSSatish Balay 
189549b5e25fSSatish Balay   /* set matrix "i" values */
189649b5e25fSSatish Balay   a->i[0] = 0;
189749b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1898574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1899574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
190049b5e25fSSatish Balay   }
19016c6c5352SBarry Smith   a->nz = a->i[mbs];
190249b5e25fSSatish Balay 
190349b5e25fSSatish Balay   /* read in nonzero values */
190487828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
190549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
190649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
190749b5e25fSSatish Balay 
190849b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
190949b5e25fSSatish Balay   nzcount = 0; jcount = 0;
191049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
191149b5e25fSSatish Balay     nzcountb = nzcount;
191249b5e25fSSatish Balay     nmask    = 0;
191349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
191449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
191549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19162d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
191703630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19182d703238SHong Zhang       }
19192d703238SHong Zhang     }
19202d703238SHong Zhang     /* sort the masked values */
19212d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19222d703238SHong Zhang 
19232d703238SHong Zhang     /* set "j" values into matrix */
19242d703238SHong Zhang     maskcount = 1;
19252d703238SHong Zhang     for (j=0; j<nmask; j++) {
192649b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
192749b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
192849b5e25fSSatish Balay     }
1929574b2666SHong Zhang 
193049b5e25fSSatish Balay     /* set "a" values into matrix */
193149b5e25fSSatish Balay     ishift = bs2*a->i[i];
193249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
193349b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
193449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1935574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1936574b2666SHong Zhang         if (tmp >= i){
193749b5e25fSSatish Balay           block     = mask[tmp] - 1;
193849b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
193949b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1940574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1941574b2666SHong Zhang         }
1942574b2666SHong Zhang         nzcountb++;
194349b5e25fSSatish Balay       }
194449b5e25fSSatish Balay     }
194549b5e25fSSatish Balay     /* zero out the mask elements we set */
194649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
194749b5e25fSSatish Balay   }
19486c6c5352SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
194949b5e25fSSatish Balay 
195049b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1951574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
195249b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
195349b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
195449b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
195549b5e25fSSatish Balay 
19569abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19579abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
195849b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19599abb65ffSKris Buschelman   *A = B;
196049b5e25fSSatish Balay   PetscFunctionReturn(0);
196149b5e25fSSatish Balay }
1962574b2666SHong Zhang 
1963d06b337dSHong Zhang #undef __FUNCT__
1964d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1965dfbe8321SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1966d06b337dSHong Zhang {
1967d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1968d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1969d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
19706849ba73SBarry Smith   PetscErrorCode ierr;
19716849ba73SBarry Smith   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j;
1972d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1973d06b337dSHong Zhang 
1974d06b337dSHong Zhang   PetscFunctionBegin;
1975b965ef7fSBarry Smith   its = its*lits;
197691723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1977d06b337dSHong Zhang 
1978d06b337dSHong Zhang   if (bs > 1)
1979d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1980d06b337dSHong Zhang 
19811ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1982d06b337dSHong Zhang   if (xx != bb) {
19831ebc52fbSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1984d06b337dSHong Zhang   } else {
1985d06b337dSHong Zhang     b = x;
1986d06b337dSHong Zhang   }
1987d06b337dSHong Zhang 
1988d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1989d06b337dSHong Zhang 
1990d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19912798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1992d06b337dSHong Zhang       for (i=0; i<m; i++)
1993d06b337dSHong Zhang         t[i] = b[i];
1994d06b337dSHong Zhang 
1995d06b337dSHong Zhang       for (i=0; i<m; i++){
199644706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1997d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1998d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1999d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2000d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2001d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
200244706875SHong Zhang         PetscLogFlops(2*nz-1);
2003d06b337dSHong Zhang       }
2004d06b337dSHong Zhang     }
2005d06b337dSHong Zhang 
20062798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2007d06b337dSHong Zhang       for (i=0; i<m; i++)
2008d06b337dSHong Zhang         t[i] = b[i];
2009d06b337dSHong Zhang 
2010d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2011d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2012d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2013d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2014d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
201544706875SHong Zhang         PetscLogFlops(2*nz-1);
2016d06b337dSHong Zhang       }
2017d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2018d06b337dSHong Zhang         d  = *(aa + ai[i]);
2019d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2020d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2021d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2022d06b337dSHong Zhang         sum = t[i];
2023d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
202444706875SHong Zhang         PetscLogFlops(2*nz-1);
2025d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2026d06b337dSHong Zhang       }
2027d06b337dSHong Zhang     }
2028d06b337dSHong Zhang     its--;
2029d06b337dSHong Zhang   }
2030d06b337dSHong Zhang 
2031d06b337dSHong Zhang   while (its--) {
203244706875SHong Zhang     /*
203344706875SHong Zhang        forward sweep:
203444706875SHong Zhang        for i=0,...,m-1:
203544706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
203644706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
203744706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2038d06b337dSHong Zhang 
203944706875SHong Zhang     */
20402798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2041d06b337dSHong Zhang       for (i=0; i<m; i++)
2042d06b337dSHong Zhang         t[i] = b[i];
2043d06b337dSHong Zhang 
2044d06b337dSHong Zhang       for (i=0; i<m; i++){
204544706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2046d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2047d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2048d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2049d06b337dSHong Zhang         sum = t[i];
2050d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2051d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2052d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
205344706875SHong Zhang         PetscLogFlops(4*nz-2);
2054d06b337dSHong Zhang       }
2055d06b337dSHong Zhang     }
2056d06b337dSHong Zhang 
20572798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
205844706875SHong Zhang       /*
205944706875SHong Zhang        backward sweep:
206044706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
206144706875SHong Zhang        for i=m-1,...,0:
206244706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
206344706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
206444706875SHong Zhang       */
2065d06b337dSHong Zhang       for (i=0; i<m; i++)
2066d06b337dSHong Zhang         t[i] = b[i];
2067d06b337dSHong Zhang 
2068d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2069d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2070d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2071d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2072d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
207344706875SHong Zhang         PetscLogFlops(2*nz-1);
2074d06b337dSHong Zhang       }
2075d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2076d06b337dSHong Zhang         d  = *(aa + ai[i]);
2077d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2078d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2079d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2080d06b337dSHong Zhang         sum = t[i];
2081d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
208244706875SHong Zhang         PetscLogFlops(2*nz-1);
2083d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2084d06b337dSHong Zhang       }
2085d06b337dSHong Zhang     }
2086d06b337dSHong Zhang   }
2087d06b337dSHong Zhang 
2088d06b337dSHong Zhang   ierr = PetscFree(t);CHKERRQ(ierr);
20891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2090d06b337dSHong Zhang   if (bb != xx) {
20911ebc52fbSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2092d06b337dSHong Zhang   }
20932798e883SHong Zhang 
2094d06b337dSHong Zhang   PetscFunctionReturn(0);
2095d06b337dSHong Zhang }
2096d06b337dSHong Zhang 
2097d06b337dSHong Zhang 
2098d06b337dSHong Zhang 
2099d06b337dSHong Zhang 
210049b5e25fSSatish Balay 
210149b5e25fSSatish Balay 
2102