xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 901853e05a144ee1b7dab90d28db377fef86bcf8)
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);
135*901853e0SKris Buschelman 
136*901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr);
137*901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr);
138*901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr);
139*901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
140*901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr);
141*901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
14249b5e25fSSatish Balay   PetscFunctionReturn(0);
14349b5e25fSSatish Balay }
14449b5e25fSSatish Balay 
1454a2ae208SSatish Balay #undef __FUNCT__
1464a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
147dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14849b5e25fSSatish Balay {
149045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
15049b5e25fSSatish Balay 
15149b5e25fSSatish Balay   PetscFunctionBegin;
1524d9d31abSKris Buschelman   switch (op) {
1534d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1544d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1554d9d31abSKris Buschelman     break;
1564d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1574d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1584d9d31abSKris Buschelman     break;
1594d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1604d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1614d9d31abSKris Buschelman     break;
1624d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1634d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1644d9d31abSKris Buschelman     break;
1654d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1664d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1674d9d31abSKris Buschelman     break;
1684d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1694d9d31abSKris Buschelman     a->nonew          = 1;
1704d9d31abSKris Buschelman     break;
1714d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1724d9d31abSKris Buschelman     a->nonew          = -1;
1734d9d31abSKris Buschelman     break;
1744d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1754d9d31abSKris Buschelman     a->nonew          = -2;
1764d9d31abSKris Buschelman     break;
1774d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1784d9d31abSKris Buschelman     a->nonew          = 0;
1794d9d31abSKris Buschelman     break;
1804d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1814d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1824d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1834d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1844d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
185b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1864d9d31abSKris Buschelman     break;
1874d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
18829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1899a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
1909a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1919a4540c5SBarry Smith   case MAT_HERMITIAN:
1929a4540c5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
19377e54ba9SKris Buschelman   case MAT_SYMMETRIC:
19477e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
1959a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
1969a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1979a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
19877e54ba9SKris Buschelman     break;
1994d9d31abSKris Buschelman   default:
20029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
20149b5e25fSSatish Balay   }
20249b5e25fSSatish Balay   PetscFunctionReturn(0);
20349b5e25fSSatish Balay }
20449b5e25fSSatish Balay 
2054a2ae208SSatish Balay #undef __FUNCT__
2064a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
207dfbe8321SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
20849b5e25fSSatish Balay {
20949b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
2106849ba73SBarry Smith   PetscErrorCode ierr;
2116849ba73SBarry Smith   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
21249b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
21387828ca2SBarry Smith   PetscScalar  *v_i;
21449b5e25fSSatish Balay 
21549b5e25fSSatish Balay   PetscFunctionBegin;
21649b5e25fSSatish Balay   bs  = a->bs;
21749b5e25fSSatish Balay   ai  = a->i;
21849b5e25fSSatish Balay   aj  = a->j;
21949b5e25fSSatish Balay   aa  = a->a;
22049b5e25fSSatish Balay   bs2 = a->bs2;
22149b5e25fSSatish Balay 
222a45adfd6SMatthew Knepley   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %d out of range", row);
22349b5e25fSSatish Balay 
22449b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
22549b5e25fSSatish Balay   bp  = row % bs; /* Block position */
22649b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
22749b5e25fSSatish Balay   *ncols = bs*M;
22849b5e25fSSatish Balay 
22949b5e25fSSatish Balay   if (v) {
23049b5e25fSSatish Balay     *v = 0;
23149b5e25fSSatish Balay     if (*ncols) {
23287828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
23349b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23449b5e25fSSatish Balay         v_i  = *v + i*bs;
23549b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
23649b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
23749b5e25fSSatish Balay       }
23849b5e25fSSatish Balay     }
23949b5e25fSSatish Balay   }
24049b5e25fSSatish Balay 
24149b5e25fSSatish Balay   if (cols) {
24249b5e25fSSatish Balay     *cols = 0;
24349b5e25fSSatish Balay     if (*ncols) {
24482502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
24549b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
24649b5e25fSSatish Balay         cols_i = *cols + i*bs;
24749b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24849b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24949b5e25fSSatish Balay       }
25049b5e25fSSatish Balay     }
25149b5e25fSSatish Balay   }
25249b5e25fSSatish Balay 
25349b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2545ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2555ddb2528SHong Zhang #ifdef column_search
25649b5e25fSSatish Balay   v_i    = *v    + M*bs;
25749b5e25fSSatish Balay   cols_i = *cols + M*bs;
25849b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25949b5e25fSSatish Balay     M = ai[i+1] - ai[i];
26049b5e25fSSatish Balay     for (j=0; j<M; j++){
26149b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
26249b5e25fSSatish Balay       if (itmp == bn){
26349b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
26449b5e25fSSatish Balay         for (k=0; k<bs; k++) {
26549b5e25fSSatish Balay           *cols_i++ = i*bs+k;
26649b5e25fSSatish Balay           *v_i++    = aa_i[k];
26749b5e25fSSatish Balay         }
26849b5e25fSSatish Balay         *ncols += bs;
26949b5e25fSSatish Balay         break;
27049b5e25fSSatish Balay       }
27149b5e25fSSatish Balay     }
27249b5e25fSSatish Balay   }
2735ddb2528SHong Zhang #endif
27449b5e25fSSatish Balay 
27549b5e25fSSatish Balay   PetscFunctionReturn(0);
27649b5e25fSSatish Balay }
27749b5e25fSSatish Balay 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
280dfbe8321SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
28149b5e25fSSatish Balay {
282dfbe8321SBarry Smith   PetscErrorCode ierr;
28349b5e25fSSatish Balay 
28449b5e25fSSatish Balay   PetscFunctionBegin;
28549b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
28649b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
28749b5e25fSSatish Balay   PetscFunctionReturn(0);
28849b5e25fSSatish Balay }
28949b5e25fSSatish Balay 
2904a2ae208SSatish Balay #undef __FUNCT__
2914a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
292dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B)
29349b5e25fSSatish Balay {
294dfbe8321SBarry Smith   PetscErrorCode ierr;
29549b5e25fSSatish Balay   PetscFunctionBegin;
296999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2978115998fSBarry Smith   PetscFunctionReturn(0);
29849b5e25fSSatish Balay }
29949b5e25fSSatish Balay 
3004a2ae208SSatish Balay #undef __FUNCT__
3014a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
3026849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
30349b5e25fSSatish Balay {
30449b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
305dfbe8321SBarry Smith   PetscErrorCode ierr;
306dfbe8321SBarry Smith   int i,j,bs = a->bs,k,l,bs2=a->bs2;
307fb9695e5SSatish Balay   char              *name;
308f3ef73ceSBarry Smith   PetscViewerFormat format;
30949b5e25fSSatish Balay 
31049b5e25fSSatish Balay   PetscFunctionBegin;
31180fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
312b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
313456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
314b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
315fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
31629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
317fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
318b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
31949b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
32049b5e25fSSatish Balay       for (j=0; j<bs; j++) {
321b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
32249b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
32349b5e25fSSatish Balay           for (l=0; l<bs; l++) {
32449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
32549b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
32662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
32749b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
32849b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
32962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
33049b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
33149b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
33262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
33349b5e25fSSatish Balay             }
33449b5e25fSSatish Balay #else
33549b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
33662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
33749b5e25fSSatish Balay             }
33849b5e25fSSatish Balay #endif
33949b5e25fSSatish Balay           }
34049b5e25fSSatish Balay         }
341b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
34249b5e25fSSatish Balay       }
34349b5e25fSSatish Balay     }
344b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
34549b5e25fSSatish Balay   } else {
346b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
34749b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
34849b5e25fSSatish Balay       for (j=0; j<bs; j++) {
349b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
35049b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
35149b5e25fSSatish Balay           for (l=0; l<bs; l++) {
35249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
35349b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
35462b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
35549b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35649b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
35762b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
35849b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35949b5e25fSSatish Balay             } else {
36062b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
36149b5e25fSSatish Balay             }
36249b5e25fSSatish Balay #else
363b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
36449b5e25fSSatish Balay #endif
36549b5e25fSSatish Balay           }
36649b5e25fSSatish Balay         }
367b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
36849b5e25fSSatish Balay       }
36949b5e25fSSatish Balay     }
370b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
37149b5e25fSSatish Balay   }
372b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
37349b5e25fSSatish Balay   PetscFunctionReturn(0);
37449b5e25fSSatish Balay }
37549b5e25fSSatish Balay 
3764a2ae208SSatish Balay #undef __FUNCT__
3774a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
3786849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
37949b5e25fSSatish Balay {
38049b5e25fSSatish Balay   Mat          A = (Mat) Aa;
38149b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
3826849ba73SBarry Smith   PetscErrorCode ierr;
3836849ba73SBarry Smith   int          row,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
38449b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
38549b5e25fSSatish Balay   MatScalar    *aa;
38649b5e25fSSatish Balay   MPI_Comm     comm;
387b0a32e0cSBarry Smith   PetscViewer  viewer;
38849b5e25fSSatish Balay 
38949b5e25fSSatish Balay   PetscFunctionBegin;
39049b5e25fSSatish Balay   /*
39149b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
39249b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
39349b5e25fSSatish Balay    rest should return immediately.
39449b5e25fSSatish Balay   */
39549b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
39649b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
39749b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
39849b5e25fSSatish Balay 
39949b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
40049b5e25fSSatish Balay 
401b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
402b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
40349b5e25fSSatish Balay 
40449b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
405b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
40649b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
40749b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
408c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
40949b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
41049b5e25fSSatish Balay       aa = a->a + j*bs2;
41149b5e25fSSatish Balay       for (k=0; k<bs; k++) {
41249b5e25fSSatish Balay         for (l=0; l<bs; l++) {
41349b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
414b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
41549b5e25fSSatish Balay         }
41649b5e25fSSatish Balay       }
41749b5e25fSSatish Balay     }
41849b5e25fSSatish Balay   }
419b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
42049b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
42149b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
422c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
42349b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
42449b5e25fSSatish Balay       aa = a->a + j*bs2;
42549b5e25fSSatish Balay       for (k=0; k<bs; k++) {
42649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
42749b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
428b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
42949b5e25fSSatish Balay         }
43049b5e25fSSatish Balay       }
43149b5e25fSSatish Balay     }
43249b5e25fSSatish Balay   }
43349b5e25fSSatish Balay 
434b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
43549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
43649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
437c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
43849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
43949b5e25fSSatish Balay       aa = a->a + j*bs2;
44049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
44149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
44249b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
443b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
44449b5e25fSSatish Balay         }
44549b5e25fSSatish Balay       }
44649b5e25fSSatish Balay     }
44749b5e25fSSatish Balay   }
44849b5e25fSSatish Balay   PetscFunctionReturn(0);
44949b5e25fSSatish Balay }
45049b5e25fSSatish Balay 
4514a2ae208SSatish Balay #undef __FUNCT__
4524a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
4536849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
45449b5e25fSSatish Balay {
455dfbe8321SBarry Smith   PetscErrorCode ierr;
45649b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
457b0a32e0cSBarry Smith   PetscDraw    draw;
45849b5e25fSSatish Balay   PetscTruth   isnull;
45949b5e25fSSatish Balay 
46049b5e25fSSatish Balay   PetscFunctionBegin;
46149b5e25fSSatish Balay 
462b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
463b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
46449b5e25fSSatish Balay 
46549b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
466c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
46749b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
468b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
469b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
47049b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
47149b5e25fSSatish Balay   PetscFunctionReturn(0);
47249b5e25fSSatish Balay }
47349b5e25fSSatish Balay 
4744a2ae208SSatish Balay #undef __FUNCT__
4754a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
476dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
47749b5e25fSSatish Balay {
478dfbe8321SBarry Smith   PetscErrorCode ierr;
47932077d6dSBarry Smith   PetscTruth iascii,isdraw;
48049b5e25fSSatish Balay 
48149b5e25fSSatish Balay   PetscFunctionBegin;
48232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
483fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
48432077d6dSBarry Smith   if (iascii){
48549b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
48649b5e25fSSatish Balay   } else if (isdraw) {
48749b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
48849b5e25fSSatish Balay   } else {
489a5e6ed63SBarry Smith     Mat B;
490a5e6ed63SBarry Smith     ierr = MatConvert(A,MATSEQAIJ,&B);CHKERRQ(ierr);
491a5e6ed63SBarry Smith     ierr = MatView(B,viewer);CHKERRQ(ierr);
492a5e6ed63SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
49349b5e25fSSatish Balay   }
49449b5e25fSSatish Balay   PetscFunctionReturn(0);
49549b5e25fSSatish Balay }
49649b5e25fSSatish Balay 
49749b5e25fSSatish Balay 
4984a2ae208SSatish Balay #undef __FUNCT__
4994a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
500dfbe8321SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[])
50149b5e25fSSatish Balay {
502045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
50349b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
50449b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
50549b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
50649b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
50749b5e25fSSatish Balay 
50849b5e25fSSatish Balay   PetscFunctionBegin;
50949b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
51049b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
511590ac198SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
512590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
51349b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
51449b5e25fSSatish Balay     nrow = ailen[brow];
51549b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
516590ac198SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
517590ac198SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1);
51849b5e25fSSatish Balay       col  = in[l] ;
51949b5e25fSSatish Balay       bcol = col/bs;
52049b5e25fSSatish Balay       cidx = col%bs;
52149b5e25fSSatish Balay       ridx = row%bs;
52249b5e25fSSatish Balay       high = nrow;
52349b5e25fSSatish Balay       low  = 0; /* assume unsorted */
52449b5e25fSSatish Balay       while (high-low > 5) {
52549b5e25fSSatish Balay         t = (low+high)/2;
52649b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
52749b5e25fSSatish Balay         else             low  = t;
52849b5e25fSSatish Balay       }
52949b5e25fSSatish Balay       for (i=low; i<high; i++) {
53049b5e25fSSatish Balay         if (rp[i] > bcol) break;
53149b5e25fSSatish Balay         if (rp[i] == bcol) {
53249b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
53349b5e25fSSatish Balay           goto finished;
53449b5e25fSSatish Balay         }
53549b5e25fSSatish Balay       }
53649b5e25fSSatish Balay       *v++ = zero;
53749b5e25fSSatish Balay       finished:;
53849b5e25fSSatish Balay     }
53949b5e25fSSatish Balay   }
54049b5e25fSSatish Balay   PetscFunctionReturn(0);
54149b5e25fSSatish Balay }
54249b5e25fSSatish Balay 
54349b5e25fSSatish Balay 
5444a2ae208SSatish Balay #undef __FUNCT__
5454a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
546dfbe8321SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
54749b5e25fSSatish Balay {
5480880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
5496849ba73SBarry Smith   PetscErrorCode ierr;
5500880e062SHong Zhang   int             *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
5510880e062SHong Zhang   int             *imax=a->imax,*ai=a->i,*ailen=a->ilen;
5526849ba73SBarry Smith   int             *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval;
5530880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
554f15d580aSBarry Smith   const MatScalar *value = v;
555f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
5560880e062SHong Zhang 
55749b5e25fSSatish Balay   PetscFunctionBegin;
5580880e062SHong Zhang   if (roworiented) {
5590880e062SHong Zhang     stepval = (n-1)*bs;
5600880e062SHong Zhang   } else {
5610880e062SHong Zhang     stepval = (m-1)*bs;
5620880e062SHong Zhang   }
5630880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
5640880e062SHong Zhang     row  = im[k];
5650880e062SHong Zhang     if (row < 0) continue;
5660880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
567590ac198SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1);
5680880e062SHong Zhang #endif
5690880e062SHong Zhang     rp   = aj + ai[row];
5700880e062SHong Zhang     ap   = aa + bs2*ai[row];
5710880e062SHong Zhang     rmax = imax[row];
5720880e062SHong Zhang     nrow = ailen[row];
5730880e062SHong Zhang     low  = 0;
5740880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
5750880e062SHong Zhang       if (in[l] < 0) continue;
5760880e062SHong Zhang       col = in[l];
577b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g)
578b1823623SSatish Balay       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",col,a->nbs-1);
579b1823623SSatish Balay #endif
5800880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
5810880e062SHong Zhang       if (roworiented) {
5820880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
5830880e062SHong Zhang       } else {
5840880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
5850880e062SHong Zhang       }
5860880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
5870880e062SHong Zhang       while (high-low > 7) {
5880880e062SHong Zhang         t = (low+high)/2;
5890880e062SHong Zhang         if (rp[t] > col) high = t;
5900880e062SHong Zhang         else             low  = t;
5910880e062SHong Zhang       }
5920880e062SHong Zhang       for (i=low; i<high; i++) {
5930880e062SHong Zhang         if (rp[i] > col) break;
5940880e062SHong Zhang         if (rp[i] == col) {
5950880e062SHong Zhang           bap  = ap +  bs2*i;
5960880e062SHong Zhang           if (roworiented) {
5970880e062SHong Zhang             if (is == ADD_VALUES) {
5980880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
5990880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6000880e062SHong Zhang                   bap[jj] += *value++;
6010880e062SHong Zhang                 }
6020880e062SHong Zhang               }
6030880e062SHong Zhang             } else {
6040880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6050880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6060880e062SHong Zhang                   bap[jj] = *value++;
6070880e062SHong Zhang                 }
6080880e062SHong Zhang               }
6090880e062SHong Zhang             }
6100880e062SHong Zhang           } else {
6110880e062SHong Zhang             if (is == ADD_VALUES) {
6120880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6130880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6140880e062SHong Zhang                   *bap++ += *value++;
6150880e062SHong Zhang                 }
6160880e062SHong Zhang               }
6170880e062SHong Zhang             } else {
6180880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6190880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6200880e062SHong Zhang                   *bap++  = *value++;
6210880e062SHong Zhang                 }
6220880e062SHong Zhang               }
6230880e062SHong Zhang             }
6240880e062SHong Zhang           }
6250880e062SHong Zhang           goto noinsert2;
6260880e062SHong Zhang         }
6270880e062SHong Zhang       }
6280880e062SHong Zhang       if (nonew == 1) goto noinsert2;
629a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6300880e062SHong Zhang       if (nrow >= rmax) {
6310880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
6320880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
6330880e062SHong Zhang         MatScalar *new_a;
6340880e062SHong Zhang 
635a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6360880e062SHong Zhang 
6370880e062SHong Zhang         /* malloc new storage space */
6380880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
6390880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
6400880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
6410880e062SHong Zhang         new_i   = new_j + new_nz;
6420880e062SHong Zhang 
6430880e062SHong Zhang         /* copy over old data into new slots */
6440880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
6450880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
6460880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
6470880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
6480880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
6490880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
6500880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
6510880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
6520880e062SHong Zhang         /* free up old matrix storage */
6530880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
6540880e062SHong Zhang         if (!a->singlemalloc) {
6550880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
6560880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
6570880e062SHong Zhang         }
6580880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
6590880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
6600880e062SHong Zhang 
6610880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
6620880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
6630880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
6646c6c5352SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
6650880e062SHong Zhang         a->reallocs++;
6666c6c5352SBarry Smith         a->nz++;
6670880e062SHong Zhang       }
6680880e062SHong Zhang       N = nrow++ - 1;
6690880e062SHong Zhang       /* shift up all the later entries in this row */
6700880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
6710880e062SHong Zhang         rp[ii+1] = rp[ii];
6720880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
6730880e062SHong Zhang       }
6740880e062SHong Zhang       if (N >= i) {
6750880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6760880e062SHong Zhang       }
6770880e062SHong Zhang       rp[i] = col;
6780880e062SHong Zhang       bap   = ap +  bs2*i;
6790880e062SHong Zhang       if (roworiented) {
6800880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6810880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
6820880e062SHong Zhang             bap[jj] = *value++;
6830880e062SHong Zhang           }
6840880e062SHong Zhang         }
6850880e062SHong Zhang       } else {
6860880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6870880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
6880880e062SHong Zhang             *bap++  = *value++;
6890880e062SHong Zhang           }
6900880e062SHong Zhang         }
6910880e062SHong Zhang       }
6920880e062SHong Zhang       noinsert2:;
6930880e062SHong Zhang       low = i;
6940880e062SHong Zhang     }
6950880e062SHong Zhang     ailen[row] = nrow;
6960880e062SHong Zhang   }
6970880e062SHong Zhang   PetscFunctionReturn(0);
69849b5e25fSSatish Balay }
69949b5e25fSSatish Balay 
7004a2ae208SSatish Balay #undef __FUNCT__
7014a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
702dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
70349b5e25fSSatish Balay {
70449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
7056849ba73SBarry Smith   PetscErrorCode ierr;
70649b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
707c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
7086849ba73SBarry Smith   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0;
70949b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
71049b5e25fSSatish Balay 
71149b5e25fSSatish Balay   PetscFunctionBegin;
71249b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
71349b5e25fSSatish Balay 
71449b5e25fSSatish Balay   if (m) rmax = ailen[0];
71549b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
71649b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
71749b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
71849b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
71949b5e25fSSatish Balay     if (fshift) {
72049b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
72149b5e25fSSatish Balay       N = ailen[i];
72249b5e25fSSatish Balay       for (j=0; j<N; j++) {
72349b5e25fSSatish Balay         ip[j-fshift] = ip[j];
72449b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
72549b5e25fSSatish Balay       }
72649b5e25fSSatish Balay     }
72749b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
72849b5e25fSSatish Balay   }
72949b5e25fSSatish Balay   if (mbs) {
73049b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
73149b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
73249b5e25fSSatish Balay   }
73349b5e25fSSatish Balay   /* reset ilen and imax for each row */
73449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
73549b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
73649b5e25fSSatish Balay   }
7376c6c5352SBarry Smith   a->nz = ai[mbs];
73849b5e25fSSatish Balay 
739b424e231SHong Zhang   /* diagonals may have moved, reset it */
740b424e231SHong Zhang   if (a->diag) {
741b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
74249b5e25fSSatish Balay   }
743b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
7446c6c5352SBarry Smith            m,A->m,a->bs,fshift*bs2,a->nz*bs2);
745b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
74649b5e25fSSatish Balay            a->reallocs);
747b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
74849b5e25fSSatish Balay   a->reallocs          = 0;
74949b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
75049b5e25fSSatish Balay 
75149b5e25fSSatish Balay   PetscFunctionReturn(0);
75249b5e25fSSatish Balay }
75349b5e25fSSatish Balay 
75449b5e25fSSatish Balay /*
75549b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
75649b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
75749b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
75849b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
75949b5e25fSSatish Balay */
7604a2ae208SSatish Balay #undef __FUNCT__
7614a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
762dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
76349b5e25fSSatish Balay {
76449b5e25fSSatish Balay   int        i,j,k,row;
76549b5e25fSSatish Balay   PetscTruth flg;
76649b5e25fSSatish Balay 
76749b5e25fSSatish Balay   PetscFunctionBegin;
76849b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
76949b5e25fSSatish Balay     row = idx[i];
77049b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
77149b5e25fSSatish Balay       sizes[j] = 1;
77249b5e25fSSatish Balay       i++;
77349b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
77449b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
77549b5e25fSSatish Balay       i++;
77649b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
77749b5e25fSSatish Balay       flg = PETSC_TRUE;
77849b5e25fSSatish Balay       for (k=1; k<bs; k++) {
77949b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
78049b5e25fSSatish Balay           flg = PETSC_FALSE;
78149b5e25fSSatish Balay           break;
78249b5e25fSSatish Balay         }
78349b5e25fSSatish Balay       }
78449b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
78549b5e25fSSatish Balay         sizes[j] = bs;
78649b5e25fSSatish Balay         i+= bs;
78749b5e25fSSatish Balay       } else {
78849b5e25fSSatish Balay         sizes[j] = 1;
78949b5e25fSSatish Balay         i++;
79049b5e25fSSatish Balay       }
79149b5e25fSSatish Balay     }
79249b5e25fSSatish Balay   }
79349b5e25fSSatish Balay   *bs_max = j;
79449b5e25fSSatish Balay   PetscFunctionReturn(0);
79549b5e25fSSatish Balay }
79649b5e25fSSatish Balay 
7974a2ae208SSatish Balay #undef __FUNCT__
7984a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
799dfbe8321SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag)
80049b5e25fSSatish Balay {
80149b5e25fSSatish Balay   PetscFunctionBegin;
802c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
80349b5e25fSSatish Balay }
80449b5e25fSSatish Balay 
80549b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
80649b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
80749b5e25fSSatish Balay */
80849b5e25fSSatish Balay 
8094a2ae208SSatish Balay #undef __FUNCT__
8104a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
811dfbe8321SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
81249b5e25fSSatish Balay {
81349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
8146849ba73SBarry Smith   PetscErrorCode ierr;
81549b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
81649b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
81749b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
8186849ba73SBarry Smith   int         ridx,cidx,bs2=a->bs2;
81949b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
82049b5e25fSSatish Balay 
82149b5e25fSSatish Balay   PetscFunctionBegin;
82249b5e25fSSatish Balay 
82349b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
82449b5e25fSSatish Balay     row  = im[k];       /* row number */
82549b5e25fSSatish Balay     brow = row/bs;      /* block row number */
82649b5e25fSSatish Balay     if (row < 0) continue;
82749b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
828590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
82949b5e25fSSatish Balay #endif
83049b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
83149b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
83249b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
83349b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
83449b5e25fSSatish Balay     low  = 0;
83549b5e25fSSatish Balay 
83649b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
83749b5e25fSSatish Balay       if (in[l] < 0) continue;
83849b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
839590ac198SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m-1);
84049b5e25fSSatish Balay #endif
84149b5e25fSSatish Balay       col = in[l];
84249b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
84349b5e25fSSatish Balay 
84449b5e25fSSatish Balay       if (brow <= bcol){
84549b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8468549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
84749b5e25fSSatish Balay           /* element value a(k,l) */
84849b5e25fSSatish Balay           if (roworiented) {
84949b5e25fSSatish Balay             value = v[l + k*n];
85049b5e25fSSatish Balay           } else {
85149b5e25fSSatish Balay             value = v[k + l*m];
85249b5e25fSSatish Balay           }
85349b5e25fSSatish Balay 
85449b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
85549b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
85649b5e25fSSatish Balay           while (high-low > 7) {
85749b5e25fSSatish Balay             t = (low+high)/2;
85849b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
85949b5e25fSSatish Balay             else              low  = t;
86049b5e25fSSatish Balay           }
86149b5e25fSSatish Balay           for (i=low; i<high; i++) {
86249b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
86349b5e25fSSatish Balay             if (rp[i] > bcol) break;
86449b5e25fSSatish Balay             if (rp[i] == bcol) {
86549b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
86649b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
86749b5e25fSSatish Balay               else                  *bap  = value;
8688549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8698549e402SHong Zhang               if (brow == bcol && ridx < cidx){
8708549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
8718549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
8728549e402SHong Zhang                 else                  *bap  = value;
8738549e402SHong Zhang               }
87449b5e25fSSatish Balay               goto noinsert1;
87549b5e25fSSatish Balay             }
87649b5e25fSSatish Balay           }
87749b5e25fSSatish Balay 
87849b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
879a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
88049b5e25fSSatish Balay       if (nrow >= rmax) {
88149b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
88249b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
88349b5e25fSSatish Balay         MatScalar *new_a;
88449b5e25fSSatish Balay 
885a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
88649b5e25fSSatish Balay 
88749b5e25fSSatish Balay         /* Malloc new storage space */
88849b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
88982502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
89049b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
89149b5e25fSSatish Balay         new_i = new_j + new_nz;
89249b5e25fSSatish Balay 
89349b5e25fSSatish Balay         /* copy over old data into new slots */
89449b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
89549b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
89649b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
89749b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
89849b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
89949b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
90049b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
90149b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
90249b5e25fSSatish Balay         /* free up old matrix storage */
90349b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
90449b5e25fSSatish Balay         if (!a->singlemalloc) {
90549b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
90649b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
90749b5e25fSSatish Balay         }
90849b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
90949b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
91049b5e25fSSatish Balay 
91149b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
91249b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
913b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
9146c6c5352SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
91549b5e25fSSatish Balay         a->reallocs++;
9166c6c5352SBarry Smith         a->nz++;
91749b5e25fSSatish Balay       }
91849b5e25fSSatish Balay 
91949b5e25fSSatish Balay       N = nrow++ - 1;
92049b5e25fSSatish Balay       /* shift up all the later entries in this row */
92149b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
92249b5e25fSSatish Balay         rp[ii+1] = rp[ii];
92349b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
92449b5e25fSSatish Balay       }
92549b5e25fSSatish Balay       if (N>=i) {
92649b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
92749b5e25fSSatish Balay       }
92849b5e25fSSatish Balay       rp[i]                      = bcol;
92949b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
93049b5e25fSSatish Balay       noinsert1:;
93149b5e25fSSatish Balay       low = i;
93249b5e25fSSatish Balay       /* } */
9338549e402SHong Zhang         }
93449b5e25fSSatish Balay       } /* end of if .. if..  */
93549b5e25fSSatish Balay     }                     /* end of loop over added columns */
93649b5e25fSSatish Balay     ailen[brow] = nrow;
93749b5e25fSSatish Balay   }                       /* end of loop over added rows */
93849b5e25fSSatish Balay 
93949b5e25fSSatish Balay   PetscFunctionReturn(0);
94049b5e25fSSatish Balay }
94149b5e25fSSatish Balay 
9426849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
9436849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
9446849ba73SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int);
9456849ba73SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
9466849ba73SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
9476849ba73SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
9486849ba73SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
9496849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat);
9506849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
9516849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
9526849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
9536849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
9546849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
9556849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
9566849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec);
95749b5e25fSSatish Balay 
9586849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
9596849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
9606849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
9616849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
9626849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
9636849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
9646849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
9656849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
9666849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
9676849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
9686849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
9696849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
9706849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
9716849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
9726849ba73SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
97349b5e25fSSatish Balay 
9746849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
975d59c15a7SBarry Smith 
9766849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
9776849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
9786849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
9796849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
9806849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
9816849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
9826849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
9836849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
98407f98182SSatish Balay 
9856849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
9866849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
9876849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
9886849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
9896849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
9906849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
9916849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
9926849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
9936849ba73SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
99449b5e25fSSatish Balay 
9956849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
9966849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
9976849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
9986849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
9996849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
10006849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
10016849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
10026849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
100349b5e25fSSatish Balay 
10046849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
10056849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
10066849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
10076849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
10086849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
10096849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
10106849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
10116849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
101249b5e25fSSatish Balay 
10134d101231SSatish Balay #ifdef HAVE_MatICCFactor
10146849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work!  */
10154a2ae208SSatish Balay #undef __FUNCT__
10164d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
1017dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
101849b5e25fSSatish Balay {
10194ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
102049b5e25fSSatish Balay   Mat         outA;
1021dfbe8321SBarry Smith   PetscErrorCode ierr;
102249b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
102349b5e25fSSatish Balay 
102449b5e25fSSatish Balay   PetscFunctionBegin;
10251a3463dfSHong Zhang   /*
102629bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
102749b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
102849b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
102949b5e25fSSatish Balay   if (!row_identity || !col_identity) {
103029bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
103149b5e25fSSatish Balay   }
10321a3463dfSHong Zhang   */
103349b5e25fSSatish Balay 
103449b5e25fSSatish Balay   outA          = inA;
10351260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
103649b5e25fSSatish Balay 
103749b5e25fSSatish Balay   if (!a->diag) {
10381a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
103949b5e25fSSatish Balay   }
104049b5e25fSSatish Balay   /*
104149b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
104249b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
104349b5e25fSSatish Balay   */
104449b5e25fSSatish Balay   switch (a->bs) {
104549b5e25fSSatish Balay   case 1:
10460fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
104707f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1048d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
10494d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
105049b5e25fSSatish Balay   case 2:
10511a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
10521a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10531a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10544d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
105549b5e25fSSatish Balay     break;
105649b5e25fSSatish Balay   case 3:
10571a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
10581a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10591a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10604d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
106149b5e25fSSatish Balay     break;
106249b5e25fSSatish Balay   case 4:
10631a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
10641a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10651a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10664d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
106749b5e25fSSatish Balay     break;
106849b5e25fSSatish Balay   case 5:
10691a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
10701a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10711a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10724d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
107349b5e25fSSatish Balay     break;
107449b5e25fSSatish Balay   case 6:
10751a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
10761a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10771a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10784d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
107949b5e25fSSatish Balay     break;
108049b5e25fSSatish Balay   case 7:
10811a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
10821a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10831a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10844d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
108549b5e25fSSatish Balay     break;
108649b5e25fSSatish Balay   default:
108749b5e25fSSatish Balay     a->row        = row;
10881a3463dfSHong Zhang     a->icol       = col;
108949b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
109049b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
109149b5e25fSSatish Balay 
109249b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
109349b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1094b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
109549b5e25fSSatish Balay 
109649b5e25fSSatish Balay     if (!a->solve_work) {
109787828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
109887828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
109949b5e25fSSatish Balay     }
110049b5e25fSSatish Balay   }
110149b5e25fSSatish Balay 
11021a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
110349b5e25fSSatish Balay 
110449b5e25fSSatish Balay   PetscFunctionReturn(0);
110549b5e25fSSatish Balay }
1106950f1e5bSHong Zhang #endif
1107950f1e5bSHong Zhang 
11084a2ae208SSatish Balay #undef __FUNCT__
11094a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
1110dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A)
111149b5e25fSSatish Balay {
111249b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
111349b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
1114dfbe8321SBarry Smith   PetscErrorCode ierr;
111549b5e25fSSatish Balay 
111649b5e25fSSatish Balay   PetscFunctionBegin;
111749b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
111849b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
111949b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
112049b5e25fSSatish Balay   PetscFunctionReturn(0);
112149b5e25fSSatish Balay }
112249b5e25fSSatish Balay 
112349b5e25fSSatish Balay EXTERN_C_BEGIN
11244a2ae208SSatish Balay #undef __FUNCT__
11254a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
1126dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
112749b5e25fSSatish Balay {
1128045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
112949b5e25fSSatish Balay   int         i,nz,n;
113049b5e25fSSatish Balay 
113149b5e25fSSatish Balay   PetscFunctionBegin;
11326c6c5352SBarry Smith   nz = baij->maxnz;
1133c464158bSHong Zhang   n  = mat->n;
113449b5e25fSSatish Balay   for (i=0; i<nz; i++) {
113549b5e25fSSatish Balay     baij->j[i] = indices[i];
113649b5e25fSSatish Balay   }
11376c6c5352SBarry Smith   baij->nz = nz;
113849b5e25fSSatish Balay   for (i=0; i<n; i++) {
113949b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
114049b5e25fSSatish Balay   }
114149b5e25fSSatish Balay 
114249b5e25fSSatish Balay   PetscFunctionReturn(0);
114349b5e25fSSatish Balay }
114449b5e25fSSatish Balay EXTERN_C_END
114549b5e25fSSatish Balay 
11464a2ae208SSatish Balay #undef __FUNCT__
11474a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
114849b5e25fSSatish Balay /*@
114919585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
115049b5e25fSSatish Balay        in the matrix.
115149b5e25fSSatish Balay 
115249b5e25fSSatish Balay   Input Parameters:
115319585528SSatish Balay +  mat     - the SeqSBAIJ matrix
115449b5e25fSSatish Balay -  indices - the column indices
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay   Level: advanced
115749b5e25fSSatish Balay 
115849b5e25fSSatish Balay   Notes:
115949b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
116049b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
116149b5e25fSSatish Balay   of the MatSetValues() operation.
116249b5e25fSSatish Balay 
116349b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
116419585528SSatish Balay   MatCreateSeqSBAIJ().
116549b5e25fSSatish Balay 
1166ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
116749b5e25fSSatish Balay 
1168ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
116949b5e25fSSatish Balay @*/
1170dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
117149b5e25fSSatish Balay {
1172dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,int *);
117349b5e25fSSatish Balay 
117449b5e25fSSatish Balay   PetscFunctionBegin;
11754482741eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
11764482741eSBarry Smith   PetscValidPointer(indices,2);
1177c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
117849b5e25fSSatish Balay   if (f) {
117949b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
118049b5e25fSSatish Balay   } else {
118129bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
118249b5e25fSSatish Balay   }
118349b5e25fSSatish Balay   PetscFunctionReturn(0);
118449b5e25fSSatish Balay }
118549b5e25fSSatish Balay 
11864a2ae208SSatish Balay #undef __FUNCT__
11874a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1188dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A)
1189273d9f13SBarry Smith {
1190dfbe8321SBarry Smith   PetscErrorCode ierr;
1191273d9f13SBarry Smith 
1192273d9f13SBarry Smith   PetscFunctionBegin;
1193273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1194273d9f13SBarry Smith   PetscFunctionReturn(0);
1195273d9f13SBarry Smith }
1196273d9f13SBarry Smith 
1197a6ece127SHong Zhang #undef __FUNCT__
1198a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1199dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1200a6ece127SHong Zhang {
1201a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1202a6ece127SHong Zhang   PetscFunctionBegin;
1203a6ece127SHong Zhang   *array = a->a;
1204a6ece127SHong Zhang   PetscFunctionReturn(0);
1205a6ece127SHong Zhang }
1206a6ece127SHong Zhang 
1207a6ece127SHong Zhang #undef __FUNCT__
1208a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1209dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1210a6ece127SHong Zhang {
1211a6ece127SHong Zhang   PetscFunctionBegin;
1212a6ece127SHong Zhang   PetscFunctionReturn(0);
1213a6ece127SHong Zhang }
1214a6ece127SHong Zhang 
121542ee4b1aSHong Zhang #include "petscblaslapack.h"
121642ee4b1aSHong Zhang #undef __FUNCT__
121742ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1218dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
121942ee4b1aSHong Zhang {
122042ee4b1aSHong Zhang   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1221dfbe8321SBarry Smith   PetscErrorCode ierr;
12224ce68768SBarry Smith   int            i,bs=y->bs,bs2,j;
12234ce68768SBarry Smith   PetscBLASInt   bnz = (PetscBLASInt)x->nz,one = 1;
122442ee4b1aSHong Zhang 
122542ee4b1aSHong Zhang   PetscFunctionBegin;
122642ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
12274ce68768SBarry Smith     BLaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one);
1228c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1229c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1230c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1231c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1232c537a176SHong Zhang     }
1233c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1234c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1235c4319e64SHong Zhang       y->XtoY = X;
1236c537a176SHong Zhang     }
1237c4319e64SHong Zhang     bs2 = bs*bs;
12386c6c5352SBarry Smith     for (i=0; i<x->nz; i++) {
1239c4319e64SHong Zhang       j = 0;
1240c4319e64SHong Zhang       while (j < bs2){
12416550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1242c4319e64SHong Zhang         j++;
1243c537a176SHong Zhang       }
1244c4319e64SHong Zhang     }
12456c6c5352SBarry 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));
124642ee4b1aSHong Zhang   } else {
124742ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
124842ee4b1aSHong Zhang   }
124942ee4b1aSHong Zhang   PetscFunctionReturn(0);
125042ee4b1aSHong Zhang }
125142ee4b1aSHong Zhang 
1252efcf0fc3SBarry Smith #undef __FUNCT__
1253efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1254dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg)
1255efcf0fc3SBarry Smith {
1256efcf0fc3SBarry Smith   PetscFunctionBegin;
1257efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1258efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1259efcf0fc3SBarry Smith }
1260efcf0fc3SBarry Smith 
1261efcf0fc3SBarry Smith #undef __FUNCT__
1262efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1263dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1264efcf0fc3SBarry Smith {
1265efcf0fc3SBarry Smith   PetscFunctionBegin;
1266efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1267efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1268efcf0fc3SBarry Smith }
1269efcf0fc3SBarry Smith 
1270efcf0fc3SBarry Smith #undef __FUNCT__
1271efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1272dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg)
1273efcf0fc3SBarry Smith {
1274efcf0fc3SBarry Smith   PetscFunctionBegin;
1275efcf0fc3SBarry Smith   *flg = PETSC_FALSE;
1276efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1277efcf0fc3SBarry Smith }
1278efcf0fc3SBarry Smith 
127949b5e25fSSatish Balay /* -------------------------------------------------------------------*/
128049b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
128149b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
128249b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
128349b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
128497304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
128549b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
128649b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
128749b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
128849b5e25fSSatish Balay        0,
128949b5e25fSSatish Balay        0,
129097304618SKris Buschelman /*10*/ 0,
129149b5e25fSSatish Balay        0,
12925f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1293d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
129449b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
129597304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
129649b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
129749b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
129849b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
129949b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
130097304618SKris Buschelman /*20*/ 0,
130149b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
130249b5e25fSSatish Balay        0,
130349b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
130449b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
130597304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ,
130649b5e25fSSatish Balay        0,
130749b5e25fSSatish Balay        0,
13085f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13095f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
131097304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1311c464158bSHong Zhang        0,
13124d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1313a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1314a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
131597304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
131649b5e25fSSatish Balay        0,
131749b5e25fSSatish Balay        0,
131849b5e25fSSatish Balay        0,
1319950f1e5bSHong Zhang        0,
132097304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
132149b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
132249b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
132349b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
132449b5e25fSSatish Balay        0,
132597304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
132649b5e25fSSatish Balay        MatScale_SeqSBAIJ,
132749b5e25fSSatish Balay        0,
132849b5e25fSSatish Balay        0,
132949b5e25fSSatish Balay        0,
133097304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ,
133149b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
133249b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
133349b5e25fSSatish Balay        0,
133449b5e25fSSatish Balay        0,
133597304618SKris Buschelman /*55*/ 0,
133649b5e25fSSatish Balay        0,
133749b5e25fSSatish Balay        0,
133849b5e25fSSatish Balay        0,
133949b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
134097304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
134149b5e25fSSatish Balay        0,
134249b5e25fSSatish Balay        0,
13438a124369SBarry Smith        MatGetPetscMaps_Petsc,
1344d959ec07SHong Zhang        0,
134597304618SKris Buschelman /*65*/ 0,
1346d959ec07SHong Zhang        0,
1347d959ec07SHong Zhang        0,
1348d959ec07SHong Zhang        0,
1349d959ec07SHong Zhang        0,
135097304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
13513e0d88b5SBarry Smith        0,
13523e0d88b5SBarry Smith        0,
13533e0d88b5SBarry Smith        0,
13543e0d88b5SBarry Smith        0,
135597304618SKris Buschelman /*75*/ 0,
13563e0d88b5SBarry Smith        0,
13573e0d88b5SBarry Smith        0,
13583e0d88b5SBarry Smith        0,
13593e0d88b5SBarry Smith        0,
136097304618SKris Buschelman /*80*/ 0,
13613e0d88b5SBarry Smith        0,
13623e0d88b5SBarry Smith        0,
13633e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
136497304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
13653e0d88b5SBarry Smith #else
136697304618SKris Buschelman        0,
13673e0d88b5SBarry Smith #endif
1368865e5f61SKris Buschelman        MatLoad_SeqSBAIJ,
1369865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ,
1370865e5f61SKris Buschelman        MatIsHermitian_SeqSBAIJ,
1371efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1372865e5f61SKris Buschelman        0,
1373865e5f61SKris Buschelman        0,
1374865e5f61SKris Buschelman /*90*/ 0,
1375865e5f61SKris Buschelman        0,
1376865e5f61SKris Buschelman        0,
1377865e5f61SKris Buschelman        0,
1378865e5f61SKris Buschelman        0,
1379865e5f61SKris Buschelman /*95*/ 0,
1380865e5f61SKris Buschelman        0,
1381865e5f61SKris Buschelman        0,
1382865e5f61SKris Buschelman        0};
138349b5e25fSSatish Balay 
138449b5e25fSSatish Balay EXTERN_C_BEGIN
13854a2ae208SSatish Balay #undef __FUNCT__
13864a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
1387dfbe8321SBarry Smith PetscErrorCode MatStoreValues_SeqSBAIJ(Mat mat)
138849b5e25fSSatish Balay {
13894afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1390c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
1391dfbe8321SBarry Smith   PetscErrorCode ierr;
139249b5e25fSSatish Balay 
139349b5e25fSSatish Balay   PetscFunctionBegin;
139449b5e25fSSatish Balay   if (aij->nonew != 1) {
139529bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
139649b5e25fSSatish Balay   }
139749b5e25fSSatish Balay 
139849b5e25fSSatish Balay   /* allocate space for values if not already there */
139949b5e25fSSatish Balay   if (!aij->saved_values) {
140087828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
140149b5e25fSSatish Balay   }
140249b5e25fSSatish Balay 
140349b5e25fSSatish Balay   /* copy values over */
140487828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
140549b5e25fSSatish Balay   PetscFunctionReturn(0);
140649b5e25fSSatish Balay }
140749b5e25fSSatish Balay EXTERN_C_END
140849b5e25fSSatish Balay 
140949b5e25fSSatish Balay EXTERN_C_BEGIN
14104a2ae208SSatish Balay #undef __FUNCT__
14114a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
1412dfbe8321SBarry Smith PetscErrorCode MatRetrieveValues_SeqSBAIJ(Mat mat)
141349b5e25fSSatish Balay {
14144afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
14156849ba73SBarry Smith   PetscErrorCode ierr;
14166849ba73SBarry Smith   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
141749b5e25fSSatish Balay 
141849b5e25fSSatish Balay   PetscFunctionBegin;
141949b5e25fSSatish Balay   if (aij->nonew != 1) {
142029bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
142149b5e25fSSatish Balay   }
142249b5e25fSSatish Balay   if (!aij->saved_values) {
142329bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
142449b5e25fSSatish Balay   }
142549b5e25fSSatish Balay 
142649b5e25fSSatish Balay   /* copy values over */
142787828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
142849b5e25fSSatish Balay   PetscFunctionReturn(0);
142949b5e25fSSatish Balay }
143049b5e25fSSatish Balay EXTERN_C_END
143149b5e25fSSatish Balay 
14328549e402SHong Zhang EXTERN_C_BEGIN
14334a2ae208SSatish Balay #undef __FUNCT__
1434a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1435dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz)
143649b5e25fSSatish Balay {
1437c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
14386849ba73SBarry Smith   PetscErrorCode ierr;
14396849ba73SBarry Smith   int          i,len,mbs,bs2;
144049b5e25fSSatish Balay   PetscTruth   flg;
144149b5e25fSSatish Balay 
144249b5e25fSSatish Balay   PetscFunctionBegin;
1443273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1444e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1445c464158bSHong Zhang   mbs  = B->m/bs;
144649b5e25fSSatish Balay   bs2  = bs*bs;
144749b5e25fSSatish Balay 
1448c464158bSHong Zhang   if (mbs*bs != B->m) {
144929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
145049b5e25fSSatish Balay   }
145149b5e25fSSatish Balay 
1452435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1453435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
145449b5e25fSSatish Balay   if (nnz) {
145549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
145629bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
145780fe4e49SBarry 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);
145849b5e25fSSatish Balay     }
145949b5e25fSSatish Balay   }
146049b5e25fSSatish Balay 
1461e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
146249b5e25fSSatish Balay   if (!flg) {
146349b5e25fSSatish Balay     switch (bs) {
146449b5e25fSSatish Balay     case 1:
14654ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
146649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1467d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
146849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
146949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
147049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
147149b5e25fSSatish Balay       break;
147249b5e25fSSatish Balay     case 2:
14734ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
147449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
147549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
147649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
147749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
147849b5e25fSSatish Balay       break;
147949b5e25fSSatish Balay     case 3:
14805f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
148149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
148249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
148349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
148449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
148549b5e25fSSatish Balay       break;
148649b5e25fSSatish Balay     case 4:
14875f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
148849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
148949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
149049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
149149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
149249b5e25fSSatish Balay       break;
149349b5e25fSSatish Balay     case 5:
14945f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
149549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
149649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
149749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
149849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
149949b5e25fSSatish Balay       break;
150049b5e25fSSatish Balay     case 6:
15015f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
150249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
150349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
150449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
150549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
150649b5e25fSSatish Balay       break;
150749b5e25fSSatish Balay     case 7:
1508de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
150949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
151049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1511de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
151249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
151349b5e25fSSatish Balay       break;
151449b5e25fSSatish Balay     }
151549b5e25fSSatish Balay   }
151649b5e25fSSatish Balay 
151749b5e25fSSatish Balay   b->mbs = mbs;
15184afc71dfSHong Zhang   b->nbs = mbs;
1519b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
152049b5e25fSSatish Balay   if (!nnz) {
1521435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
152249b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
152349b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
15248cef66ccSHong Zhang       b->imax[i] = nz;
152549b5e25fSSatish Balay     }
1526153ea458SHong Zhang     nz = nz*mbs; /* total nz */
152749b5e25fSSatish Balay   } else {
152849b5e25fSSatish Balay     nz = 0;
15298cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
153049b5e25fSSatish Balay   }
15316c6c5352SBarry Smith   /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
153249b5e25fSSatish Balay 
153349b5e25fSSatish Balay   /* allocate the matrix space */
15346c6c5352SBarry Smith   len  = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
153582502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
15366c6c5352SBarry Smith   ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
15376c6c5352SBarry Smith   b->j = (int*)(b->a + nz*bs2);
15386c6c5352SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
15396c6c5352SBarry Smith   b->i = b->j + nz;
154049b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
154149b5e25fSSatish Balay 
154249b5e25fSSatish Balay   /* pointer to beginning of each row */
154349b5e25fSSatish Balay   b->i[0] = 0;
154449b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
154549b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
154649b5e25fSSatish Balay   }
154749b5e25fSSatish Balay 
154849b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
15495033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1550b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
155149b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
155249b5e25fSSatish Balay 
155349b5e25fSSatish Balay   b->bs               = bs;
155449b5e25fSSatish Balay   b->bs2              = bs2;
15556c6c5352SBarry Smith   b->nz             = 0;
15566c6c5352SBarry Smith   b->maxnz          = nz*bs2;
1557153ea458SHong Zhang 
155816cdd363SHong Zhang   b->inew             = 0;
155916cdd363SHong Zhang   b->jnew             = 0;
156016cdd363SHong Zhang   b->anew             = 0;
156116cdd363SHong Zhang   b->a2anew           = 0;
15621a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1563c464158bSHong Zhang   PetscFunctionReturn(0);
1564c464158bSHong Zhang }
1565a23d5eceSKris Buschelman EXTERN_C_END
1566153ea458SHong Zhang 
1567d769727bSBarry Smith EXTERN_C_BEGIN
1568dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,Mat*);
1569dfbe8321SBarry Smith EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,Mat*);
1570d769727bSBarry Smith EXTERN_C_END
1571d769727bSBarry Smith 
15720bad9183SKris Buschelman /*MC
1573fafad747SKris Buschelman    MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
15740bad9183SKris Buschelman    based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
15750bad9183SKris Buschelman 
15760bad9183SKris Buschelman    Options Database Keys:
15770bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
15780bad9183SKris Buschelman 
15790bad9183SKris Buschelman   Level: beginner
15800bad9183SKris Buschelman 
15810bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ
15820bad9183SKris Buschelman M*/
15830bad9183SKris Buschelman 
1584a23d5eceSKris Buschelman EXTERN_C_BEGIN
1585a23d5eceSKris Buschelman #undef __FUNCT__
1586a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1587dfbe8321SBarry Smith PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1588a23d5eceSKris Buschelman {
1589a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
1590dfbe8321SBarry Smith   PetscErrorCode ierr;
1591dfbe8321SBarry Smith   int size;
1592a23d5eceSKris Buschelman 
1593a23d5eceSKris Buschelman   PetscFunctionBegin;
1594a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1595a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1596a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1597a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1598a23d5eceSKris Buschelman 
1599a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1600a23d5eceSKris Buschelman   B->data = (void*)b;
1601a23d5eceSKris Buschelman   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1602a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1603a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1604a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1605a23d5eceSKris Buschelman   B->factor           = 0;
1606a23d5eceSKris Buschelman   B->lupivotthreshold = 1.0;
1607a23d5eceSKris Buschelman   B->mapping          = 0;
1608a23d5eceSKris Buschelman   b->row              = 0;
1609a23d5eceSKris Buschelman   b->icol             = 0;
1610a23d5eceSKris Buschelman   b->reallocs         = 0;
1611a23d5eceSKris Buschelman   b->saved_values     = 0;
1612a23d5eceSKris Buschelman 
1613a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1614a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1615a23d5eceSKris Buschelman 
1616a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1617a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1618a23d5eceSKris Buschelman   b->nonew            = 0;
1619a23d5eceSKris Buschelman   b->diag             = 0;
1620a23d5eceSKris Buschelman   b->solve_work       = 0;
1621a23d5eceSKris Buschelman   b->mult_work        = 0;
1622a23d5eceSKris Buschelman   B->spptr            = 0;
1623a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1624a23d5eceSKris Buschelman   b->xtoy             = 0;
1625a23d5eceSKris Buschelman   b->XtoY             = 0;
1626a23d5eceSKris Buschelman 
1627a23d5eceSKris Buschelman   b->inew             = 0;
1628a23d5eceSKris Buschelman   b->jnew             = 0;
1629a23d5eceSKris Buschelman   b->anew             = 0;
1630a23d5eceSKris Buschelman   b->a2anew           = 0;
1631a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1632a23d5eceSKris Buschelman 
1633a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1634a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1635a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1636a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1637a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1638a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1639a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1640a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1641a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
16424e5e7fe4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",
16434e5e7fe4SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqAIJ",
16444e5e7fe4SHong Zhang                                       MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1645a0e1a404SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",
1646a0e1a404SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqBAIJ",
1647a0e1a404SHong Zhang                                       MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1648a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1649a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1650a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
165123ce1328SBarry Smith 
165223ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
165323ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
165423ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
165523ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
1656a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1657a23d5eceSKris Buschelman }
1658a23d5eceSKris Buschelman EXTERN_C_END
1659a23d5eceSKris Buschelman 
1660a23d5eceSKris Buschelman #undef __FUNCT__
1661a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1662a23d5eceSKris Buschelman /*@C
1663a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1664a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1665a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1666a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1667a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1668a23d5eceSKris Buschelman 
1669a23d5eceSKris Buschelman    Collective on Mat
1670a23d5eceSKris Buschelman 
1671a23d5eceSKris Buschelman    Input Parameters:
1672a23d5eceSKris Buschelman +  A - the symmetric matrix
1673a23d5eceSKris Buschelman .  bs - size of block
1674a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1675a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1676a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1677a23d5eceSKris Buschelman 
1678a23d5eceSKris Buschelman    Options Database Keys:
1679a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1680a23d5eceSKris Buschelman                      block calculations (much slower)
1681a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1682a23d5eceSKris Buschelman 
1683a23d5eceSKris Buschelman    Level: intermediate
1684a23d5eceSKris Buschelman 
1685a23d5eceSKris Buschelman    Notes:
1686a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1687a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1688a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1689a23d5eceSKris Buschelman    matrices.
1690a23d5eceSKris Buschelman 
1691a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1692a23d5eceSKris Buschelman @*/
1693dfbe8321SBarry Smith PetscErrorCode MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) {
1694dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,int,int,const int[]);
1695a23d5eceSKris Buschelman 
1696a23d5eceSKris Buschelman   PetscFunctionBegin;
1697a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1698a23d5eceSKris Buschelman   if (f) {
1699a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1700a23d5eceSKris Buschelman   }
1701a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1702a23d5eceSKris Buschelman }
170349b5e25fSSatish Balay 
17044a2ae208SSatish Balay #undef __FUNCT__
17054a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1706c464158bSHong Zhang /*@C
1707c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1708c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1709c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1710c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1711c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
171249b5e25fSSatish Balay 
1713c464158bSHong Zhang    Collective on MPI_Comm
1714c464158bSHong Zhang 
1715c464158bSHong Zhang    Input Parameters:
1716c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1717c464158bSHong Zhang .  bs - size of block
1718c464158bSHong Zhang .  m - number of rows, or number of columns
1719c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1720744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1721744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1722c464158bSHong Zhang 
1723c464158bSHong Zhang    Output Parameter:
1724c464158bSHong Zhang .  A - the symmetric matrix
1725c464158bSHong Zhang 
1726c464158bSHong Zhang    Options Database Keys:
1727c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1728c464158bSHong Zhang                      block calculations (much slower)
1729c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1730c464158bSHong Zhang 
1731c464158bSHong Zhang    Level: intermediate
1732c464158bSHong Zhang 
1733c464158bSHong Zhang    Notes:
1734c464158bSHong Zhang 
1735c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1736c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1737c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1738c464158bSHong Zhang    matrices.
1739c464158bSHong Zhang 
1740c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1741c464158bSHong Zhang @*/
1742dfbe8321SBarry Smith PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A)
1743c464158bSHong Zhang {
1744dfbe8321SBarry Smith   PetscErrorCode ierr;
1745c464158bSHong Zhang 
1746c464158bSHong Zhang   PetscFunctionBegin;
1747c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1748c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1749273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
175049b5e25fSSatish Balay   PetscFunctionReturn(0);
175149b5e25fSSatish Balay }
175249b5e25fSSatish Balay 
17534a2ae208SSatish Balay #undef __FUNCT__
17544a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
1755dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
175649b5e25fSSatish Balay {
175749b5e25fSSatish Balay   Mat          C;
175849b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
17596849ba73SBarry Smith   PetscErrorCode ierr;
17606849ba73SBarry Smith   int          i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
176149b5e25fSSatish Balay 
176249b5e25fSSatish Balay   PetscFunctionBegin;
176329bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
176449b5e25fSSatish Balay 
176549b5e25fSSatish Balay   *B = 0;
1766692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1767be5d1d56SKris Buschelman   ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
17681d5dac46SHong Zhang   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1769692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1770692f9cbeSHong Zhang 
1771273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
177249b5e25fSSatish Balay   C->factor         = A->factor;
177349b5e25fSSatish Balay   c->row            = 0;
177449b5e25fSSatish Balay   c->icol           = 0;
177549b5e25fSSatish Balay   c->saved_values   = 0;
177649b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
177749b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
177849b5e25fSSatish Balay 
177982327fa8SHong Zhang   C->M    = A->M;
178082327fa8SHong Zhang   C->N    = A->N;
178149b5e25fSSatish Balay   c->bs   = a->bs;
178249b5e25fSSatish Balay   c->bs2  = a->bs2;
178349b5e25fSSatish Balay   c->mbs  = a->mbs;
178449b5e25fSSatish Balay   c->nbs  = a->nbs;
178549b5e25fSSatish Balay 
1786b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1787b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
178849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
178949b5e25fSSatish Balay     c->imax[i] = a->imax[i];
179049b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
179149b5e25fSSatish Balay   }
179249b5e25fSSatish Balay 
179349b5e25fSSatish Balay   /* allocate the matrix space */
179449b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
179549b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
179682502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
179749b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
179849b5e25fSSatish Balay   c->i = c->j + nz;
179949b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
180049b5e25fSSatish Balay   if (mbs > 0) {
180149b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
180249b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
180349b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
180449b5e25fSSatish Balay     } else {
180549b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
180649b5e25fSSatish Balay     }
180749b5e25fSSatish Balay   }
180849b5e25fSSatish Balay 
1809b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
181049b5e25fSSatish Balay   c->sorted      = a->sorted;
181149b5e25fSSatish Balay   c->roworiented = a->roworiented;
181249b5e25fSSatish Balay   c->nonew       = a->nonew;
181349b5e25fSSatish Balay 
181449b5e25fSSatish Balay   if (a->diag) {
1815b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1816b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
181749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
181849b5e25fSSatish Balay       c->diag[i] = a->diag[i];
181949b5e25fSSatish Balay     }
182049b5e25fSSatish Balay   } else c->diag        = 0;
18216c6c5352SBarry Smith   c->nz               = a->nz;
18226c6c5352SBarry Smith   c->maxnz            = a->maxnz;
182349b5e25fSSatish Balay   c->solve_work         = 0;
182449b5e25fSSatish Balay   c->mult_work          = 0;
182549b5e25fSSatish Balay   *B = C;
1826b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
182749b5e25fSSatish Balay   PetscFunctionReturn(0);
182849b5e25fSSatish Balay }
182949b5e25fSSatish Balay 
18304a2ae208SSatish Balay #undef __FUNCT__
18314a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1832dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A)
183349b5e25fSSatish Balay {
183449b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
183549b5e25fSSatish Balay   Mat          B;
18366849ba73SBarry Smith   PetscErrorCode ierr;
18376849ba73SBarry Smith   int          i,nz,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18383f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
183949b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
184049b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
184187828ca2SBarry Smith   PetscScalar  *aa;
184249b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
184349b5e25fSSatish Balay 
184449b5e25fSSatish Balay   PetscFunctionBegin;
1845b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
184649b5e25fSSatish Balay   bs2  = bs*bs;
184749b5e25fSSatish Balay 
184849b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
184929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1850b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
185149b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1852552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
185349b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
185449b5e25fSSatish Balay 
185549b5e25fSSatish Balay   if (header[3] < 0) {
185629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
185749b5e25fSSatish Balay   }
185849b5e25fSSatish Balay 
185929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
186049b5e25fSSatish Balay 
186149b5e25fSSatish Balay   /*
186249b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
186349b5e25fSSatish Balay     divisible by the blocksize
186449b5e25fSSatish Balay   */
186549b5e25fSSatish Balay   mbs        = M/bs;
186649b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
186749b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
186849b5e25fSSatish Balay   else                  mbs++;
186949b5e25fSSatish Balay   if (extra_rows) {
1870b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
187149b5e25fSSatish Balay   }
187249b5e25fSSatish Balay 
187349b5e25fSSatish Balay   /* read in row lengths */
1874b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
187549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
187649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
187749b5e25fSSatish Balay 
187849b5e25fSSatish Balay   /* read in column indices */
1879b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
188049b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
188149b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
188249b5e25fSSatish Balay 
188349b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
188482502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1885a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
188682502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
188749b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
188849b5e25fSSatish Balay   masked   = mask + mbs;
188949b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
189049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
189149b5e25fSSatish Balay     nmask = 0;
189249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
189349b5e25fSSatish Balay       kmax = rowlengths[rowcount];
189449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18952d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
189603630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
189749b5e25fSSatish Balay       }
189849b5e25fSSatish Balay       rowcount++;
189949b5e25fSSatish Balay     }
1900574b2666SHong Zhang     s_browlengths[i] += nmask;
1901574b2666SHong Zhang 
190249b5e25fSSatish Balay     /* zero out the mask elements we set */
190349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
190449b5e25fSSatish Balay   }
190549b5e25fSSatish Balay 
190649b5e25fSSatish Balay   /* create our matrix */
19079abb65ffSKris Buschelman   ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr);
19089abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
190978473fd7SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr);
191049b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
191149b5e25fSSatish Balay 
191249b5e25fSSatish Balay   /* set matrix "i" values */
191349b5e25fSSatish Balay   a->i[0] = 0;
191449b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1915574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1916574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
191749b5e25fSSatish Balay   }
19186c6c5352SBarry Smith   a->nz = a->i[mbs];
191949b5e25fSSatish Balay 
192049b5e25fSSatish Balay   /* read in nonzero values */
192187828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
192249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
192349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
192449b5e25fSSatish Balay 
192549b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
192649b5e25fSSatish Balay   nzcount = 0; jcount = 0;
192749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
192849b5e25fSSatish Balay     nzcountb = nzcount;
192949b5e25fSSatish Balay     nmask    = 0;
193049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
193149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
193249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19332d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
193403630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19352d703238SHong Zhang       }
19362d703238SHong Zhang     }
19372d703238SHong Zhang     /* sort the masked values */
19382d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19392d703238SHong Zhang 
19402d703238SHong Zhang     /* set "j" values into matrix */
19412d703238SHong Zhang     maskcount = 1;
19422d703238SHong Zhang     for (j=0; j<nmask; j++) {
194349b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
194449b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
194549b5e25fSSatish Balay     }
1946574b2666SHong Zhang 
194749b5e25fSSatish Balay     /* set "a" values into matrix */
194849b5e25fSSatish Balay     ishift = bs2*a->i[i];
194949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
195049b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
195149b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1952574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1953574b2666SHong Zhang         if (tmp >= i){
195449b5e25fSSatish Balay           block     = mask[tmp] - 1;
195549b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
195649b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1957574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1958574b2666SHong Zhang         }
1959574b2666SHong Zhang         nzcountb++;
196049b5e25fSSatish Balay       }
196149b5e25fSSatish Balay     }
196249b5e25fSSatish Balay     /* zero out the mask elements we set */
196349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
196449b5e25fSSatish Balay   }
19656c6c5352SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
196649b5e25fSSatish Balay 
196749b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1968574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
196949b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
197049b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
197149b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
197249b5e25fSSatish Balay 
19739abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19749abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
197549b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19769abb65ffSKris Buschelman   *A = B;
197749b5e25fSSatish Balay   PetscFunctionReturn(0);
197849b5e25fSSatish Balay }
1979574b2666SHong Zhang 
1980d06b337dSHong Zhang #undef __FUNCT__
1981d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1982dfbe8321SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1983d06b337dSHong Zhang {
1984d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1985d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1986d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
19876849ba73SBarry Smith   PetscErrorCode ierr;
19886849ba73SBarry Smith   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j;
1989d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1990d06b337dSHong Zhang 
1991d06b337dSHong Zhang   PetscFunctionBegin;
1992b965ef7fSBarry Smith   its = its*lits;
199391723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1994d06b337dSHong Zhang 
1995d06b337dSHong Zhang   if (bs > 1)
1996d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1997d06b337dSHong Zhang 
19981ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1999d06b337dSHong Zhang   if (xx != bb) {
20001ebc52fbSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
2001d06b337dSHong Zhang   } else {
2002d06b337dSHong Zhang     b = x;
2003d06b337dSHong Zhang   }
2004d06b337dSHong Zhang 
2005d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
2006d06b337dSHong Zhang 
2007d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
20082798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2009d06b337dSHong Zhang       for (i=0; i<m; i++)
2010d06b337dSHong Zhang         t[i] = b[i];
2011d06b337dSHong Zhang 
2012d06b337dSHong Zhang       for (i=0; i<m; i++){
201344706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2014d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2015d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2016d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2017d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2018d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
201944706875SHong Zhang         PetscLogFlops(2*nz-1);
2020d06b337dSHong Zhang       }
2021d06b337dSHong Zhang     }
2022d06b337dSHong Zhang 
20232798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2024d06b337dSHong Zhang       for (i=0; i<m; i++)
2025d06b337dSHong Zhang         t[i] = b[i];
2026d06b337dSHong Zhang 
2027d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2028d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2029d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2030d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2031d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
203244706875SHong Zhang         PetscLogFlops(2*nz-1);
2033d06b337dSHong Zhang       }
2034d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2035d06b337dSHong Zhang         d  = *(aa + ai[i]);
2036d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2037d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2038d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2039d06b337dSHong Zhang         sum = t[i];
2040d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
204144706875SHong Zhang         PetscLogFlops(2*nz-1);
2042d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2043d06b337dSHong Zhang       }
2044d06b337dSHong Zhang     }
2045d06b337dSHong Zhang     its--;
2046d06b337dSHong Zhang   }
2047d06b337dSHong Zhang 
2048d06b337dSHong Zhang   while (its--) {
204944706875SHong Zhang     /*
205044706875SHong Zhang        forward sweep:
205144706875SHong Zhang        for i=0,...,m-1:
205244706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
205344706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
205444706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2055d06b337dSHong Zhang 
205644706875SHong Zhang     */
20572798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2058d06b337dSHong Zhang       for (i=0; i<m; i++)
2059d06b337dSHong Zhang         t[i] = b[i];
2060d06b337dSHong Zhang 
2061d06b337dSHong Zhang       for (i=0; i<m; i++){
206244706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2063d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2064d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2065d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2066d06b337dSHong Zhang         sum = t[i];
2067d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2068d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2069d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
207044706875SHong Zhang         PetscLogFlops(4*nz-2);
2071d06b337dSHong Zhang       }
2072d06b337dSHong Zhang     }
2073d06b337dSHong Zhang 
20742798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
207544706875SHong Zhang       /*
207644706875SHong Zhang        backward sweep:
207744706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
207844706875SHong Zhang        for i=m-1,...,0:
207944706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
208044706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
208144706875SHong Zhang       */
2082d06b337dSHong Zhang       for (i=0; i<m; i++)
2083d06b337dSHong Zhang         t[i] = b[i];
2084d06b337dSHong Zhang 
2085d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2086d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2087d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2088d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2089d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
209044706875SHong Zhang         PetscLogFlops(2*nz-1);
2091d06b337dSHong Zhang       }
2092d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2093d06b337dSHong Zhang         d  = *(aa + ai[i]);
2094d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2095d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2096d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2097d06b337dSHong Zhang         sum = t[i];
2098d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
209944706875SHong Zhang         PetscLogFlops(2*nz-1);
2100d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2101d06b337dSHong Zhang       }
2102d06b337dSHong Zhang     }
2103d06b337dSHong Zhang   }
2104d06b337dSHong Zhang 
2105d06b337dSHong Zhang   ierr = PetscFree(t);CHKERRQ(ierr);
21061ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2107d06b337dSHong Zhang   if (bb != xx) {
21081ebc52fbSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2109d06b337dSHong Zhang   }
21102798e883SHong Zhang 
2111d06b337dSHong Zhang   PetscFunctionReturn(0);
2112d06b337dSHong Zhang }
2113d06b337dSHong Zhang 
2114d06b337dSHong Zhang 
2115d06b337dSHong Zhang 
2116d06b337dSHong Zhang 
211749b5e25fSSatish Balay 
211849b5e25fSSatish Balay 
2119