xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision efcf0fc353b705cc10d937cc67ff6a314bf622fa)
173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
4a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h"         /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/vec/vecimpl.h"
949b5e25fSSatish Balay #include "src/inline/spops.h"
10aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1149b5e25fSSatish Balay 
1249b5e25fSSatish Balay #define CHUNKSIZE  10
1349b5e25fSSatish Balay 
1449b5e25fSSatish Balay /*
1549b5e25fSSatish Balay      Checks for missing diagonals
1649b5e25fSSatish Balay */
174a2ae208SSatish Balay #undef __FUNCT__
184a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
1949b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A)
2049b5e25fSSatish Balay {
21045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2249b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2349b5e25fSSatish Balay 
2449b5e25fSSatish Balay   PetscFunctionBegin;
25045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2649b5e25fSSatish Balay   diag = a->diag;
2749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
28e11b0e14SHong Zhang     if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i);
2949b5e25fSSatish Balay   }
3049b5e25fSSatish Balay   PetscFunctionReturn(0);
3149b5e25fSSatish Balay }
3249b5e25fSSatish Balay 
334a2ae208SSatish Balay #undef __FUNCT__
344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
3549b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A)
3649b5e25fSSatish Balay {
37045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
38b424e231SHong Zhang   int          i,mbs = a->mbs,ierr;
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"
514e7234bfSBarry Smith static int 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 */
62a6ece127SHong Zhang     int s_nz = a->i[n];
63a6ece127SHong Zhang     for (i=0; i<s_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"
744e7234bfSBarry Smith static int 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) {
83a6ece127SHong Zhang     int s_nz = a->i[n]-1;
84a6ece127SHong Zhang     for (i=0; i<s_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"
9249b5e25fSSatish Balay int 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"
10349b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
10449b5e25fSSatish Balay {
10549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10649b5e25fSSatish Balay   int         ierr;
10749b5e25fSSatish Balay 
10849b5e25fSSatish Balay   PetscFunctionBegin;
10949b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
110b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz);
11149b5e25fSSatish Balay #endif
11249b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
11349b5e25fSSatish Balay   if (!a->singlemalloc) {
11449b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
11563c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
11649b5e25fSSatish Balay   }
11749b5e25fSSatish Balay   if (a->row) {
11849b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
11949b5e25fSSatish Balay   }
12049b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
12149b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
12249b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
12349b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
124d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
125d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
126d59c15a7SBarry Smith   if (a->mult_work)   {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
12749b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
128c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
1291a3463dfSHong Zhang 
1301a3463dfSHong Zhang   if (a->inew){
1311a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1321a3463dfSHong Zhang     a->inew = 0;
1331a3463dfSHong Zhang   }
13449b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13549b5e25fSSatish Balay   PetscFunctionReturn(0);
13649b5e25fSSatish Balay }
13749b5e25fSSatish Balay 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
14049b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14149b5e25fSSatish Balay {
142045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14349b5e25fSSatish Balay 
14449b5e25fSSatish Balay   PetscFunctionBegin;
1454d9d31abSKris Buschelman   switch (op) {
1464d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1474d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1484d9d31abSKris Buschelman     break;
1494d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1504d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1514d9d31abSKris Buschelman     break;
1524d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1534d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1564d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1594d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1624d9d31abSKris Buschelman     a->nonew          = 1;
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1654d9d31abSKris Buschelman     a->nonew          = -1;
1664d9d31abSKris Buschelman     break;
1674d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1684d9d31abSKris Buschelman     a->nonew          = -2;
1694d9d31abSKris Buschelman     break;
1704d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1714d9d31abSKris Buschelman     a->nonew          = 0;
1724d9d31abSKris Buschelman     break;
1734d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1744d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1754d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1764d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1774d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
178b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1794d9d31abSKris Buschelman     break;
1804d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
18129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
18277e54ba9SKris Buschelman   case MAT_SYMMETRIC:
18377e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
18477e54ba9SKris Buschelman     break;
1854d9d31abSKris Buschelman   default:
18629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
18749b5e25fSSatish Balay   }
18849b5e25fSSatish Balay   PetscFunctionReturn(0);
18949b5e25fSSatish Balay }
19049b5e25fSSatish Balay 
1914a2ae208SSatish Balay #undef __FUNCT__
1924a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19387828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
19449b5e25fSSatish Balay {
19549b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
19682502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
19749b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
19887828ca2SBarry Smith   PetscScalar  *v_i;
19949b5e25fSSatish Balay 
20049b5e25fSSatish Balay   PetscFunctionBegin;
20149b5e25fSSatish Balay   bs  = a->bs;
20249b5e25fSSatish Balay   ai  = a->i;
20349b5e25fSSatish Balay   aj  = a->j;
20449b5e25fSSatish Balay   aa  = a->a;
20549b5e25fSSatish Balay   bs2 = a->bs2;
20649b5e25fSSatish Balay 
207a45adfd6SMatthew Knepley   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %d out of range", row);
20849b5e25fSSatish Balay 
20949b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21049b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21149b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
21249b5e25fSSatish Balay   *ncols = bs*M;
21349b5e25fSSatish Balay 
21449b5e25fSSatish Balay   if (v) {
21549b5e25fSSatish Balay     *v = 0;
21649b5e25fSSatish Balay     if (*ncols) {
21787828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
21849b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
21949b5e25fSSatish Balay         v_i  = *v + i*bs;
22049b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22149b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
22249b5e25fSSatish Balay       }
22349b5e25fSSatish Balay     }
22449b5e25fSSatish Balay   }
22549b5e25fSSatish Balay 
22649b5e25fSSatish Balay   if (cols) {
22749b5e25fSSatish Balay     *cols = 0;
22849b5e25fSSatish Balay     if (*ncols) {
22982502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
23049b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23149b5e25fSSatish Balay         cols_i = *cols + i*bs;
23249b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
23349b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
23449b5e25fSSatish Balay       }
23549b5e25fSSatish Balay     }
23649b5e25fSSatish Balay   }
23749b5e25fSSatish Balay 
23849b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2395ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2405ddb2528SHong Zhang #ifdef column_search
24149b5e25fSSatish Balay   v_i    = *v    + M*bs;
24249b5e25fSSatish Balay   cols_i = *cols + M*bs;
24349b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
24449b5e25fSSatish Balay     M = ai[i+1] - ai[i];
24549b5e25fSSatish Balay     for (j=0; j<M; j++){
24649b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
24749b5e25fSSatish Balay       if (itmp == bn){
24849b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
24949b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25049b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25149b5e25fSSatish Balay           *v_i++    = aa_i[k];
25249b5e25fSSatish Balay         }
25349b5e25fSSatish Balay         *ncols += bs;
25449b5e25fSSatish Balay         break;
25549b5e25fSSatish Balay       }
25649b5e25fSSatish Balay     }
25749b5e25fSSatish Balay   }
2585ddb2528SHong Zhang #endif
25949b5e25fSSatish Balay 
26049b5e25fSSatish Balay   PetscFunctionReturn(0);
26149b5e25fSSatish Balay }
26249b5e25fSSatish Balay 
2634a2ae208SSatish Balay #undef __FUNCT__
2644a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
26587828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
26649b5e25fSSatish Balay {
26749b5e25fSSatish Balay   int ierr;
26849b5e25fSSatish Balay 
26949b5e25fSSatish Balay   PetscFunctionBegin;
27049b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27149b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
27249b5e25fSSatish Balay   PetscFunctionReturn(0);
27349b5e25fSSatish Balay }
27449b5e25fSSatish Balay 
2754a2ae208SSatish Balay #undef __FUNCT__
2764a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
27749b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
27849b5e25fSSatish Balay {
2798115998fSBarry Smith   int ierr;
28049b5e25fSSatish Balay   PetscFunctionBegin;
281999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2828115998fSBarry Smith   PetscFunctionReturn(0);
28349b5e25fSSatish Balay }
28449b5e25fSSatish Balay 
2854a2ae208SSatish Balay #undef __FUNCT__
2864a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary"
287b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
28849b5e25fSSatish Balay {
28949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
29049b5e25fSSatish Balay   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
29187828ca2SBarry Smith   PetscScalar  *aa;
29249b5e25fSSatish Balay   FILE         *file;
29349b5e25fSSatish Balay 
29449b5e25fSSatish Balay   PetscFunctionBegin;
295b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
29682502324SSatish Balay   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
297552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
29849b5e25fSSatish Balay 
299c464158bSHong Zhang   col_lens[1] = A->m;
300c464158bSHong Zhang   col_lens[2] = A->m;
30149b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
30249b5e25fSSatish Balay 
30349b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
30449b5e25fSSatish Balay   count = 0;
30549b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
30649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
30749b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
30849b5e25fSSatish Balay     }
30949b5e25fSSatish Balay   }
310c464158bSHong Zhang   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
31149b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
31249b5e25fSSatish Balay 
31349b5e25fSSatish Balay   /* store column indices (zero start index) */
31482502324SSatish Balay   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
31549b5e25fSSatish Balay   count = 0;
31649b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31749b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31849b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
31949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
32049b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
32149b5e25fSSatish Balay         }
32249b5e25fSSatish Balay       }
32349b5e25fSSatish Balay     }
32449b5e25fSSatish Balay   }
32549b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
32649b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
32749b5e25fSSatish Balay 
32849b5e25fSSatish Balay   /* store nonzero values */
32987828ca2SBarry Smith   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
33049b5e25fSSatish Balay   count = 0;
33149b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
33249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
33349b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
33449b5e25fSSatish Balay         for (l=0; l<bs; l++) {
33549b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
33649b5e25fSSatish Balay         }
33749b5e25fSSatish Balay       }
33849b5e25fSSatish Balay     }
33949b5e25fSSatish Balay   }
34049b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
34149b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
34249b5e25fSSatish Balay 
343b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
34449b5e25fSSatish Balay   if (file) {
34549b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
34649b5e25fSSatish Balay   }
34749b5e25fSSatish Balay   PetscFunctionReturn(0);
34849b5e25fSSatish Balay }
34949b5e25fSSatish Balay 
3504a2ae208SSatish Balay #undef __FUNCT__
3514a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
352b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
35349b5e25fSSatish Balay {
35449b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
355fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
356fb9695e5SSatish Balay   char              *name;
357f3ef73ceSBarry Smith   PetscViewerFormat format;
35849b5e25fSSatish Balay 
35949b5e25fSSatish Balay   PetscFunctionBegin;
36080fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
361b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
362456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
363b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
364fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
36529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
366fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
367b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
36849b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
36949b5e25fSSatish Balay       for (j=0; j<bs; j++) {
370b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
37149b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
37249b5e25fSSatish Balay           for (l=0; l<bs; l++) {
37349b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
37449b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
37649b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37749b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37862b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
37949b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38049b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38162b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38249b5e25fSSatish Balay             }
38349b5e25fSSatish Balay #else
38449b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
38562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
38649b5e25fSSatish Balay             }
38749b5e25fSSatish Balay #endif
38849b5e25fSSatish Balay           }
38949b5e25fSSatish Balay         }
390b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39149b5e25fSSatish Balay       }
39249b5e25fSSatish Balay     }
393b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
39449b5e25fSSatish Balay   } else {
395b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
39649b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
39749b5e25fSSatish Balay       for (j=0; j<bs; j++) {
398b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
39949b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
40049b5e25fSSatish Balay           for (l=0; l<bs; l++) {
40149b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
40249b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
40362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
40449b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40549b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
40662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
40749b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40849b5e25fSSatish Balay             } else {
40962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41049b5e25fSSatish Balay             }
41149b5e25fSSatish Balay #else
412b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
41349b5e25fSSatish Balay #endif
41449b5e25fSSatish Balay           }
41549b5e25fSSatish Balay         }
416b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
41749b5e25fSSatish Balay       }
41849b5e25fSSatish Balay     }
419b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
42049b5e25fSSatish Balay   }
421b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
42249b5e25fSSatish Balay   PetscFunctionReturn(0);
42349b5e25fSSatish Balay }
42449b5e25fSSatish Balay 
4254a2ae208SSatish Balay #undef __FUNCT__
4264a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
427b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
42849b5e25fSSatish Balay {
42949b5e25fSSatish Balay   Mat          A = (Mat) Aa;
43049b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
43149b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
43249b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
43349b5e25fSSatish Balay   MatScalar    *aa;
43449b5e25fSSatish Balay   MPI_Comm     comm;
435b0a32e0cSBarry Smith   PetscViewer  viewer;
43649b5e25fSSatish Balay 
43749b5e25fSSatish Balay   PetscFunctionBegin;
43849b5e25fSSatish Balay   /*
43949b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
44049b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
44149b5e25fSSatish Balay    rest should return immediately.
44249b5e25fSSatish Balay   */
44349b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
44449b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
44549b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
44649b5e25fSSatish Balay 
44749b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
44849b5e25fSSatish Balay 
449b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
450b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
45149b5e25fSSatish Balay 
45249b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
453b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
45449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
456c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
45749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
45849b5e25fSSatish Balay       aa = a->a + j*bs2;
45949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
46049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
46149b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
462b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
46349b5e25fSSatish Balay         }
46449b5e25fSSatish Balay       }
46549b5e25fSSatish Balay     }
46649b5e25fSSatish Balay   }
467b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
46849b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46949b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
470c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
47149b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
47249b5e25fSSatish Balay       aa = a->a + j*bs2;
47349b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47449b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47549b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
476b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47749b5e25fSSatish Balay         }
47849b5e25fSSatish Balay       }
47949b5e25fSSatish Balay     }
48049b5e25fSSatish Balay   }
48149b5e25fSSatish Balay 
482b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
48349b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
48449b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
485c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
48649b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
48749b5e25fSSatish Balay       aa = a->a + j*bs2;
48849b5e25fSSatish Balay       for (k=0; k<bs; k++) {
48949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
49049b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
491b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
49249b5e25fSSatish Balay         }
49349b5e25fSSatish Balay       }
49449b5e25fSSatish Balay     }
49549b5e25fSSatish Balay   }
49649b5e25fSSatish Balay   PetscFunctionReturn(0);
49749b5e25fSSatish Balay }
49849b5e25fSSatish Balay 
4994a2ae208SSatish Balay #undef __FUNCT__
5004a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
501b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
50249b5e25fSSatish Balay {
50349b5e25fSSatish Balay   int          ierr;
50449b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
505b0a32e0cSBarry Smith   PetscDraw    draw;
50649b5e25fSSatish Balay   PetscTruth   isnull;
50749b5e25fSSatish Balay 
50849b5e25fSSatish Balay   PetscFunctionBegin;
50949b5e25fSSatish Balay 
510b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
511b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
51249b5e25fSSatish Balay 
51349b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
514c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
51549b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
516b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
517b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
51849b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
51949b5e25fSSatish Balay   PetscFunctionReturn(0);
52049b5e25fSSatish Balay }
52149b5e25fSSatish Balay 
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
524b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
52549b5e25fSSatish Balay {
52649b5e25fSSatish Balay   int        ierr;
52749b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
52849b5e25fSSatish Balay 
52949b5e25fSSatish Balay   PetscFunctionBegin;
530b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
531b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
532fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
533fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
53449b5e25fSSatish Balay   if (issocket) {
53529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
53649b5e25fSSatish Balay   } else if (isascii){
53749b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
53849b5e25fSSatish Balay   } else if (isbinary) {
53949b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
54049b5e25fSSatish Balay   } else if (isdraw) {
54149b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
54249b5e25fSSatish Balay   } else {
543b7aaefc3SHong Zhang     SETERRQ1(1,"Viewer type %s not supported by SeqSBAIJ matrices",((PetscObject)viewer)->type_name);
54449b5e25fSSatish Balay   }
54549b5e25fSSatish Balay   PetscFunctionReturn(0);
54649b5e25fSSatish Balay }
54749b5e25fSSatish Balay 
54849b5e25fSSatish Balay 
5494a2ae208SSatish Balay #undef __FUNCT__
5504a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
551f15d580aSBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[])
55249b5e25fSSatish Balay {
553045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
55449b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
55549b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
55649b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
55749b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
55849b5e25fSSatish Balay 
55949b5e25fSSatish Balay   PetscFunctionBegin;
56049b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
56149b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
562590ac198SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
563590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
56449b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
56549b5e25fSSatish Balay     nrow = ailen[brow];
56649b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
567590ac198SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
568590ac198SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1);
56949b5e25fSSatish Balay       col  = in[l] ;
57049b5e25fSSatish Balay       bcol = col/bs;
57149b5e25fSSatish Balay       cidx = col%bs;
57249b5e25fSSatish Balay       ridx = row%bs;
57349b5e25fSSatish Balay       high = nrow;
57449b5e25fSSatish Balay       low  = 0; /* assume unsorted */
57549b5e25fSSatish Balay       while (high-low > 5) {
57649b5e25fSSatish Balay         t = (low+high)/2;
57749b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
57849b5e25fSSatish Balay         else             low  = t;
57949b5e25fSSatish Balay       }
58049b5e25fSSatish Balay       for (i=low; i<high; i++) {
58149b5e25fSSatish Balay         if (rp[i] > bcol) break;
58249b5e25fSSatish Balay         if (rp[i] == bcol) {
58349b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
58449b5e25fSSatish Balay           goto finished;
58549b5e25fSSatish Balay         }
58649b5e25fSSatish Balay       }
58749b5e25fSSatish Balay       *v++ = zero;
58849b5e25fSSatish Balay       finished:;
58949b5e25fSSatish Balay     }
59049b5e25fSSatish Balay   }
59149b5e25fSSatish Balay   PetscFunctionReturn(0);
59249b5e25fSSatish Balay }
59349b5e25fSSatish Balay 
59449b5e25fSSatish Balay 
5954a2ae208SSatish Balay #undef __FUNCT__
5964a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
597f15d580aSBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
59849b5e25fSSatish Balay {
5990880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
6000880e062SHong Zhang   int             *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
6010880e062SHong Zhang   int             *imax=a->imax,*ai=a->i,*ailen=a->ilen;
6020880e062SHong Zhang   int             *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
6030880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
604f15d580aSBarry Smith   const MatScalar *value = v;
605f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
6060880e062SHong Zhang 
60749b5e25fSSatish Balay   PetscFunctionBegin;
6080880e062SHong Zhang   if (roworiented) {
6090880e062SHong Zhang     stepval = (n-1)*bs;
6100880e062SHong Zhang   } else {
6110880e062SHong Zhang     stepval = (m-1)*bs;
6120880e062SHong Zhang   }
6130880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
6140880e062SHong Zhang     row  = im[k];
6150880e062SHong Zhang     if (row < 0) continue;
6160880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
617590ac198SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1);
6180880e062SHong Zhang #endif
6190880e062SHong Zhang     rp   = aj + ai[row];
6200880e062SHong Zhang     ap   = aa + bs2*ai[row];
6210880e062SHong Zhang     rmax = imax[row];
6220880e062SHong Zhang     nrow = ailen[row];
6230880e062SHong Zhang     low  = 0;
6240880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
6250880e062SHong Zhang       if (in[l] < 0) continue;
6260880e062SHong Zhang       col = in[l];
627b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g)
628b1823623SSatish Balay       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",col,a->nbs-1);
629b1823623SSatish Balay #endif
6300880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
6310880e062SHong Zhang       if (roworiented) {
6320880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
6330880e062SHong Zhang       } else {
6340880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
6350880e062SHong Zhang       }
6360880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
6370880e062SHong Zhang       while (high-low > 7) {
6380880e062SHong Zhang         t = (low+high)/2;
6390880e062SHong Zhang         if (rp[t] > col) high = t;
6400880e062SHong Zhang         else             low  = t;
6410880e062SHong Zhang       }
6420880e062SHong Zhang       for (i=low; i<high; i++) {
6430880e062SHong Zhang         if (rp[i] > col) break;
6440880e062SHong Zhang         if (rp[i] == col) {
6450880e062SHong Zhang           bap  = ap +  bs2*i;
6460880e062SHong Zhang           if (roworiented) {
6470880e062SHong Zhang             if (is == ADD_VALUES) {
6480880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6490880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6500880e062SHong Zhang                   bap[jj] += *value++;
6510880e062SHong Zhang                 }
6520880e062SHong Zhang               }
6530880e062SHong Zhang             } else {
6540880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6550880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6560880e062SHong Zhang                   bap[jj] = *value++;
6570880e062SHong Zhang                 }
6580880e062SHong Zhang               }
6590880e062SHong Zhang             }
6600880e062SHong Zhang           } else {
6610880e062SHong Zhang             if (is == ADD_VALUES) {
6620880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6630880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6640880e062SHong Zhang                   *bap++ += *value++;
6650880e062SHong Zhang                 }
6660880e062SHong Zhang               }
6670880e062SHong Zhang             } else {
6680880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6690880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6700880e062SHong Zhang                   *bap++  = *value++;
6710880e062SHong Zhang                 }
6720880e062SHong Zhang               }
6730880e062SHong Zhang             }
6740880e062SHong Zhang           }
6750880e062SHong Zhang           goto noinsert2;
6760880e062SHong Zhang         }
6770880e062SHong Zhang       }
6780880e062SHong Zhang       if (nonew == 1) goto noinsert2;
679a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6800880e062SHong Zhang       if (nrow >= rmax) {
6810880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
6820880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
6830880e062SHong Zhang         MatScalar *new_a;
6840880e062SHong Zhang 
685a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6860880e062SHong Zhang 
6870880e062SHong Zhang         /* malloc new storage space */
6880880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
6890880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
6900880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
6910880e062SHong Zhang         new_i   = new_j + new_nz;
6920880e062SHong Zhang 
6930880e062SHong Zhang         /* copy over old data into new slots */
6940880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
6950880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
6960880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
6970880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
6980880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
6990880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
7000880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
7010880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
7020880e062SHong Zhang         /* free up old matrix storage */
7030880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
7040880e062SHong Zhang         if (!a->singlemalloc) {
7050880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
7060880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
7070880e062SHong Zhang         }
7080880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
7090880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
7100880e062SHong Zhang 
7110880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
7120880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
7130880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
7140880e062SHong Zhang         a->s_maxnz += bs2*CHUNKSIZE;
7150880e062SHong Zhang         a->reallocs++;
7160880e062SHong Zhang         a->s_nz++;
7170880e062SHong Zhang       }
7180880e062SHong Zhang       N = nrow++ - 1;
7190880e062SHong Zhang       /* shift up all the later entries in this row */
7200880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
7210880e062SHong Zhang         rp[ii+1] = rp[ii];
7220880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
7230880e062SHong Zhang       }
7240880e062SHong Zhang       if (N >= i) {
7250880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
7260880e062SHong Zhang       }
7270880e062SHong Zhang       rp[i] = col;
7280880e062SHong Zhang       bap   = ap +  bs2*i;
7290880e062SHong Zhang       if (roworiented) {
7300880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7310880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
7320880e062SHong Zhang             bap[jj] = *value++;
7330880e062SHong Zhang           }
7340880e062SHong Zhang         }
7350880e062SHong Zhang       } else {
7360880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7370880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
7380880e062SHong Zhang             *bap++  = *value++;
7390880e062SHong Zhang           }
7400880e062SHong Zhang         }
7410880e062SHong Zhang       }
7420880e062SHong Zhang       noinsert2:;
7430880e062SHong Zhang       low = i;
7440880e062SHong Zhang     }
7450880e062SHong Zhang     ailen[row] = nrow;
7460880e062SHong Zhang   }
7470880e062SHong Zhang   PetscFunctionReturn(0);
74849b5e25fSSatish Balay }
74949b5e25fSSatish Balay 
7504a2ae208SSatish Balay #undef __FUNCT__
7514a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
75249b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
75349b5e25fSSatish Balay {
75449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
75549b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
756c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
75749b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
75849b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
75949b5e25fSSatish Balay 
76049b5e25fSSatish Balay   PetscFunctionBegin;
76149b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
76249b5e25fSSatish Balay 
76349b5e25fSSatish Balay   if (m) rmax = ailen[0];
76449b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
76549b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
76649b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
76749b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
76849b5e25fSSatish Balay     if (fshift) {
76949b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
77049b5e25fSSatish Balay       N = ailen[i];
77149b5e25fSSatish Balay       for (j=0; j<N; j++) {
77249b5e25fSSatish Balay         ip[j-fshift] = ip[j];
77349b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
77449b5e25fSSatish Balay       }
77549b5e25fSSatish Balay     }
77649b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
77749b5e25fSSatish Balay   }
77849b5e25fSSatish Balay   if (mbs) {
77949b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
78049b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
78149b5e25fSSatish Balay   }
78249b5e25fSSatish Balay   /* reset ilen and imax for each row */
78349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
78449b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
78549b5e25fSSatish Balay   }
78649b5e25fSSatish Balay   a->s_nz = ai[mbs];
78749b5e25fSSatish Balay 
788b424e231SHong Zhang   /* diagonals may have moved, reset it */
789b424e231SHong Zhang   if (a->diag) {
790b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
79149b5e25fSSatish Balay   }
792b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
793c464158bSHong Zhang            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
794b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
79549b5e25fSSatish Balay            a->reallocs);
796b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
79749b5e25fSSatish Balay   a->reallocs          = 0;
79849b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
79949b5e25fSSatish Balay 
80049b5e25fSSatish Balay   PetscFunctionReturn(0);
80149b5e25fSSatish Balay }
80249b5e25fSSatish Balay 
80349b5e25fSSatish Balay /*
80449b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
80549b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
80649b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
80749b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
80849b5e25fSSatish Balay */
8094a2ae208SSatish Balay #undef __FUNCT__
8104a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
811db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
81249b5e25fSSatish Balay {
81349b5e25fSSatish Balay   int        i,j,k,row;
81449b5e25fSSatish Balay   PetscTruth flg;
81549b5e25fSSatish Balay 
81649b5e25fSSatish Balay   PetscFunctionBegin;
81749b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
81849b5e25fSSatish Balay     row = idx[i];
81949b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
82049b5e25fSSatish Balay       sizes[j] = 1;
82149b5e25fSSatish Balay       i++;
82249b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
82349b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
82449b5e25fSSatish Balay       i++;
82549b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
82649b5e25fSSatish Balay       flg = PETSC_TRUE;
82749b5e25fSSatish Balay       for (k=1; k<bs; k++) {
82849b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
82949b5e25fSSatish Balay           flg = PETSC_FALSE;
83049b5e25fSSatish Balay           break;
83149b5e25fSSatish Balay         }
83249b5e25fSSatish Balay       }
83349b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
83449b5e25fSSatish Balay         sizes[j] = bs;
83549b5e25fSSatish Balay         i+= bs;
83649b5e25fSSatish Balay       } else {
83749b5e25fSSatish Balay         sizes[j] = 1;
83849b5e25fSSatish Balay         i++;
83949b5e25fSSatish Balay       }
84049b5e25fSSatish Balay     }
84149b5e25fSSatish Balay   }
84249b5e25fSSatish Balay   *bs_max = j;
84349b5e25fSSatish Balay   PetscFunctionReturn(0);
84449b5e25fSSatish Balay }
84549b5e25fSSatish Balay 
8464a2ae208SSatish Balay #undef __FUNCT__
8474a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
848268466fbSBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag)
84949b5e25fSSatish Balay {
85049b5e25fSSatish Balay   PetscFunctionBegin;
851c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
85249b5e25fSSatish Balay }
85349b5e25fSSatish Balay 
85449b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
85549b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
85649b5e25fSSatish Balay */
85749b5e25fSSatish Balay 
8584a2ae208SSatish Balay #undef __FUNCT__
8594a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
860f15d580aSBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
86149b5e25fSSatish Balay {
86249b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
86349b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
86449b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
86549b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
86649b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
86749b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
86849b5e25fSSatish Balay 
86949b5e25fSSatish Balay   PetscFunctionBegin;
87049b5e25fSSatish Balay 
87149b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
87249b5e25fSSatish Balay     row  = im[k];       /* row number */
87349b5e25fSSatish Balay     brow = row/bs;      /* block row number */
87449b5e25fSSatish Balay     if (row < 0) continue;
87549b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
876590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
87749b5e25fSSatish Balay #endif
87849b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
87949b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
88049b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
88149b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
88249b5e25fSSatish Balay     low  = 0;
88349b5e25fSSatish Balay 
88449b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
88549b5e25fSSatish Balay       if (in[l] < 0) continue;
88649b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
887590ac198SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m-1);
88849b5e25fSSatish Balay #endif
88949b5e25fSSatish Balay       col = in[l];
89049b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
89149b5e25fSSatish Balay 
89249b5e25fSSatish Balay       if (brow <= bcol){
89349b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8948549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
89549b5e25fSSatish Balay           /* element value a(k,l) */
89649b5e25fSSatish Balay           if (roworiented) {
89749b5e25fSSatish Balay             value = v[l + k*n];
89849b5e25fSSatish Balay           } else {
89949b5e25fSSatish Balay             value = v[k + l*m];
90049b5e25fSSatish Balay           }
90149b5e25fSSatish Balay 
90249b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
90349b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
90449b5e25fSSatish Balay           while (high-low > 7) {
90549b5e25fSSatish Balay             t = (low+high)/2;
90649b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
90749b5e25fSSatish Balay             else              low  = t;
90849b5e25fSSatish Balay           }
90949b5e25fSSatish Balay           for (i=low; i<high; i++) {
91049b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
91149b5e25fSSatish Balay             if (rp[i] > bcol) break;
91249b5e25fSSatish Balay             if (rp[i] == bcol) {
91349b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
91449b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
91549b5e25fSSatish Balay               else                  *bap  = value;
9168549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9178549e402SHong Zhang               if (brow == bcol && ridx < cidx){
9188549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
9198549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
9208549e402SHong Zhang                 else                  *bap  = value;
9218549e402SHong Zhang               }
92249b5e25fSSatish Balay               goto noinsert1;
92349b5e25fSSatish Balay             }
92449b5e25fSSatish Balay           }
92549b5e25fSSatish Balay 
92649b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
927a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
92849b5e25fSSatish Balay       if (nrow >= rmax) {
92949b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
93049b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
93149b5e25fSSatish Balay         MatScalar *new_a;
93249b5e25fSSatish Balay 
933a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
93449b5e25fSSatish Balay 
93549b5e25fSSatish Balay         /* Malloc new storage space */
93649b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
93782502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
93849b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
93949b5e25fSSatish Balay         new_i = new_j + new_nz;
94049b5e25fSSatish Balay 
94149b5e25fSSatish Balay         /* copy over old data into new slots */
94249b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
94349b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
94449b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
94549b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
94649b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
94749b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
94849b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
94949b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
95049b5e25fSSatish Balay         /* free up old matrix storage */
95149b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
95249b5e25fSSatish Balay         if (!a->singlemalloc) {
95349b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
95449b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
95549b5e25fSSatish Balay         }
95649b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
95749b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
95849b5e25fSSatish Balay 
95949b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
96049b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
961b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
96249b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
96349b5e25fSSatish Balay         a->reallocs++;
96449b5e25fSSatish Balay         a->s_nz++;
96549b5e25fSSatish Balay       }
96649b5e25fSSatish Balay 
96749b5e25fSSatish Balay       N = nrow++ - 1;
96849b5e25fSSatish Balay       /* shift up all the later entries in this row */
96949b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
97049b5e25fSSatish Balay         rp[ii+1] = rp[ii];
97149b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
97249b5e25fSSatish Balay       }
97349b5e25fSSatish Balay       if (N>=i) {
97449b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
97549b5e25fSSatish Balay       }
97649b5e25fSSatish Balay       rp[i]                      = bcol;
97749b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
97849b5e25fSSatish Balay       noinsert1:;
97949b5e25fSSatish Balay       low = i;
98049b5e25fSSatish Balay       /* } */
9818549e402SHong Zhang         }
98249b5e25fSSatish Balay       } /* end of if .. if..  */
98349b5e25fSSatish Balay     }                     /* end of loop over added columns */
98449b5e25fSSatish Balay     ailen[brow] = nrow;
98549b5e25fSSatish Balay   }                       /* end of loop over added rows */
98649b5e25fSSatish Balay 
98749b5e25fSSatish Balay   PetscFunctionReturn(0);
98849b5e25fSSatish Balay }
98949b5e25fSSatish Balay 
99015e8a5b3SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
99115e8a5b3SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
99252ebccd6SSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int);
99349b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
994268466fbSBarry Smith extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
99549b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
99649b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
997268466fbSBarry Smith extern int MatScale_SeqSBAIJ(const PetscScalar*,Mat);
99849b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
99949b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
100049b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
100149b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
100249b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
100349b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
100413a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
100549b5e25fSSatish Balay 
100649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
100749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
100849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
100949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
101049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
101149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
101249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
101349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
101449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
101549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
101649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
101749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
101849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
101949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
102049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
102149b5e25fSSatish Balay 
1022d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1023d59c15a7SBarry Smith 
102407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
102507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
102607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
102707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
102807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
102907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
103007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
103107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
103207f98182SSatish Balay 
10335f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
10345f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
10355f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
10365f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
10375f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
10385f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
10395f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
104057d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
10413e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
104249b5e25fSSatish Balay 
104349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
104449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
104549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
104649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
104749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
104849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
104949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
105049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
105149b5e25fSSatish Balay 
105249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
105349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
105449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
105549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
105649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
105749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
105849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
105949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
106049b5e25fSSatish Balay 
10614d101231SSatish Balay #ifdef HAVE_MatICCFactor
10621a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
10634a2ae208SSatish Balay #undef __FUNCT__
10644d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
106515e8a5b3SHong Zhang int MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
106649b5e25fSSatish Balay {
10674ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
106849b5e25fSSatish Balay   Mat         outA;
106949b5e25fSSatish Balay   int         ierr;
107049b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
107149b5e25fSSatish Balay 
107249b5e25fSSatish Balay   PetscFunctionBegin;
10731a3463dfSHong Zhang   /*
107429bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
107549b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
107649b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
107749b5e25fSSatish Balay   if (!row_identity || !col_identity) {
107829bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
107949b5e25fSSatish Balay   }
10801a3463dfSHong Zhang   */
108149b5e25fSSatish Balay 
108249b5e25fSSatish Balay   outA          = inA;
10831260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
108449b5e25fSSatish Balay 
108549b5e25fSSatish Balay   if (!a->diag) {
10861a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
108749b5e25fSSatish Balay   }
108849b5e25fSSatish Balay   /*
108949b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
109049b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
109149b5e25fSSatish Balay   */
109249b5e25fSSatish Balay   switch (a->bs) {
109349b5e25fSSatish Balay   case 1:
10940fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
109507f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1096d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
10974d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
109849b5e25fSSatish Balay   case 2:
10991a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
11001a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11011a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11024d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
110349b5e25fSSatish Balay     break;
110449b5e25fSSatish Balay   case 3:
11051a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
11061a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11071a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11084d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
110949b5e25fSSatish Balay     break;
111049b5e25fSSatish Balay   case 4:
11111a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
11121a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11131a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11144d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
111549b5e25fSSatish Balay     break;
111649b5e25fSSatish Balay   case 5:
11171a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
11181a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11191a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11204d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
112149b5e25fSSatish Balay     break;
112249b5e25fSSatish Balay   case 6:
11231a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
11241a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11251a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11264d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
112749b5e25fSSatish Balay     break;
112849b5e25fSSatish Balay   case 7:
11291a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
11301a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11311a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11324d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
113349b5e25fSSatish Balay     break;
113449b5e25fSSatish Balay   default:
113549b5e25fSSatish Balay     a->row        = row;
11361a3463dfSHong Zhang     a->icol       = col;
113749b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
113849b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
113949b5e25fSSatish Balay 
114049b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
114149b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1142b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
114349b5e25fSSatish Balay 
114449b5e25fSSatish Balay     if (!a->solve_work) {
114587828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
114687828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
114749b5e25fSSatish Balay     }
114849b5e25fSSatish Balay   }
114949b5e25fSSatish Balay 
11501a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
115149b5e25fSSatish Balay 
115249b5e25fSSatish Balay   PetscFunctionReturn(0);
115349b5e25fSSatish Balay }
1154950f1e5bSHong Zhang #endif
1155950f1e5bSHong Zhang 
11564a2ae208SSatish Balay #undef __FUNCT__
11574a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
115849b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
115949b5e25fSSatish Balay {
116049b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
116149b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
116249b5e25fSSatish Balay   int               ierr;
116349b5e25fSSatish Balay 
116449b5e25fSSatish Balay   PetscFunctionBegin;
116549b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
116649b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
116749b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
116849b5e25fSSatish Balay   PetscFunctionReturn(0);
116949b5e25fSSatish Balay }
117049b5e25fSSatish Balay 
117149b5e25fSSatish Balay EXTERN_C_BEGIN
11724a2ae208SSatish Balay #undef __FUNCT__
11734a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
117449b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
117549b5e25fSSatish Balay {
1176045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
117749b5e25fSSatish Balay   int         i,nz,n;
117849b5e25fSSatish Balay 
117949b5e25fSSatish Balay   PetscFunctionBegin;
1180045c9aa0SHong Zhang   nz = baij->s_maxnz;
1181c464158bSHong Zhang   n  = mat->n;
118249b5e25fSSatish Balay   for (i=0; i<nz; i++) {
118349b5e25fSSatish Balay     baij->j[i] = indices[i];
118449b5e25fSSatish Balay   }
1185045c9aa0SHong Zhang   baij->s_nz = nz;
118649b5e25fSSatish Balay   for (i=0; i<n; i++) {
118749b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
118849b5e25fSSatish Balay   }
118949b5e25fSSatish Balay 
119049b5e25fSSatish Balay   PetscFunctionReturn(0);
119149b5e25fSSatish Balay }
119249b5e25fSSatish Balay EXTERN_C_END
119349b5e25fSSatish Balay 
11944a2ae208SSatish Balay #undef __FUNCT__
11954a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
119649b5e25fSSatish Balay /*@
119719585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
119849b5e25fSSatish Balay        in the matrix.
119949b5e25fSSatish Balay 
120049b5e25fSSatish Balay   Input Parameters:
120119585528SSatish Balay +  mat     - the SeqSBAIJ matrix
120249b5e25fSSatish Balay -  indices - the column indices
120349b5e25fSSatish Balay 
120449b5e25fSSatish Balay   Level: advanced
120549b5e25fSSatish Balay 
120649b5e25fSSatish Balay   Notes:
120749b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
120849b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
120949b5e25fSSatish Balay   of the MatSetValues() operation.
121049b5e25fSSatish Balay 
121149b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
121219585528SSatish Balay   MatCreateSeqSBAIJ().
121349b5e25fSSatish Balay 
1214ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
121549b5e25fSSatish Balay 
1216ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
121749b5e25fSSatish Balay @*/
121849b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
121949b5e25fSSatish Balay {
122049b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
122149b5e25fSSatish Balay 
122249b5e25fSSatish Balay   PetscFunctionBegin;
122349b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1224c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
122549b5e25fSSatish Balay   if (f) {
122649b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
122749b5e25fSSatish Balay   } else {
122829bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
122949b5e25fSSatish Balay   }
123049b5e25fSSatish Balay   PetscFunctionReturn(0);
123149b5e25fSSatish Balay }
123249b5e25fSSatish Balay 
12334a2ae208SSatish Balay #undef __FUNCT__
12344a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1235273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1236273d9f13SBarry Smith {
1237273d9f13SBarry Smith   int        ierr;
1238273d9f13SBarry Smith 
1239273d9f13SBarry Smith   PetscFunctionBegin;
1240273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1241273d9f13SBarry Smith   PetscFunctionReturn(0);
1242273d9f13SBarry Smith }
1243273d9f13SBarry Smith 
1244a6ece127SHong Zhang #undef __FUNCT__
1245a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
12464e7234bfSBarry Smith int MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1247a6ece127SHong Zhang {
1248a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1249a6ece127SHong Zhang   PetscFunctionBegin;
1250a6ece127SHong Zhang   *array = a->a;
1251a6ece127SHong Zhang   PetscFunctionReturn(0);
1252a6ece127SHong Zhang }
1253a6ece127SHong Zhang 
1254a6ece127SHong Zhang #undef __FUNCT__
1255a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
12564e7234bfSBarry Smith int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1257a6ece127SHong Zhang {
1258a6ece127SHong Zhang   PetscFunctionBegin;
1259a6ece127SHong Zhang   PetscFunctionReturn(0);
1260a6ece127SHong Zhang }
1261a6ece127SHong Zhang 
126242ee4b1aSHong Zhang #include "petscblaslapack.h"
126342ee4b1aSHong Zhang #undef __FUNCT__
126442ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1265268466fbSBarry Smith int MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
126642ee4b1aSHong Zhang {
126742ee4b1aSHong Zhang   Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1268c4319e64SHong Zhang   int          ierr,one=1,i,bs=y->bs,bs2,j;
126942ee4b1aSHong Zhang 
127042ee4b1aSHong Zhang   PetscFunctionBegin;
127142ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1272268466fbSBarry Smith     BLaxpy_(&x->s_nz,(PetscScalar*)a,x->a,&one,y->a,&one);
1273c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1274c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1275c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1276c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1277c537a176SHong Zhang     }
1278c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1279c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1280c4319e64SHong Zhang       y->XtoY = X;
1281c537a176SHong Zhang     }
1282c4319e64SHong Zhang     bs2 = bs*bs;
1283c537a176SHong Zhang     for (i=0; i<x->s_nz; i++) {
1284c4319e64SHong Zhang       j = 0;
1285c4319e64SHong Zhang       while (j < bs2){
12866550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1287c4319e64SHong Zhang         j++;
1288c537a176SHong Zhang       }
1289c4319e64SHong Zhang     }
1290c4319e64SHong Zhang     PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %d/%d = %g\n",bs2*x->s_nz,bs2*y->s_nz,(PetscReal)(bs2*x->s_nz)/(bs2*y->s_nz));
129142ee4b1aSHong Zhang   } else {
129242ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
129342ee4b1aSHong Zhang   }
129442ee4b1aSHong Zhang   PetscFunctionReturn(0);
129542ee4b1aSHong Zhang }
129642ee4b1aSHong Zhang 
1297*efcf0fc3SBarry Smith #undef __FUNCT__
1298*efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1299*efcf0fc3SBarry Smith int MatIsSymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1300*efcf0fc3SBarry Smith {
1301*efcf0fc3SBarry Smith   PetscFunctionBegin;
1302*efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1303*efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1304*efcf0fc3SBarry Smith }
1305*efcf0fc3SBarry Smith 
1306*efcf0fc3SBarry Smith #undef __FUNCT__
1307*efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1308*efcf0fc3SBarry Smith int MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1309*efcf0fc3SBarry Smith {
1310*efcf0fc3SBarry Smith   PetscFunctionBegin;
1311*efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1312*efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1313*efcf0fc3SBarry Smith }
1314*efcf0fc3SBarry Smith 
1315*efcf0fc3SBarry Smith #undef __FUNCT__
1316*efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1317*efcf0fc3SBarry Smith int MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg)
1318*efcf0fc3SBarry Smith {
1319*efcf0fc3SBarry Smith   PetscFunctionBegin;
1320*efcf0fc3SBarry Smith   *flg = PETSC_FALSE;
1321*efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1322*efcf0fc3SBarry Smith }
1323*efcf0fc3SBarry Smith 
132449b5e25fSSatish Balay /* -------------------------------------------------------------------*/
132549b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
132649b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
132749b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
132849b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
132997304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
133049b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
133149b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
133249b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
133349b5e25fSSatish Balay        0,
133449b5e25fSSatish Balay        0,
133597304618SKris Buschelman /*10*/ 0,
133649b5e25fSSatish Balay        0,
13375f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1338d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
133949b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
134097304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
134149b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
134249b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
134349b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
134449b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
134597304618SKris Buschelman /*20*/ 0,
134649b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
134749b5e25fSSatish Balay        0,
134849b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
134949b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
135097304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ,
135149b5e25fSSatish Balay        0,
135249b5e25fSSatish Balay        0,
13535f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13545f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
135597304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1356c464158bSHong Zhang        0,
13574d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1358a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1359a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
136097304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
136149b5e25fSSatish Balay        0,
136249b5e25fSSatish Balay        0,
136349b5e25fSSatish Balay        0,
1364950f1e5bSHong Zhang        0,
136597304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
136649b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
136749b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
136849b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
136949b5e25fSSatish Balay        0,
137097304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
137149b5e25fSSatish Balay        MatScale_SeqSBAIJ,
137249b5e25fSSatish Balay        0,
137349b5e25fSSatish Balay        0,
137449b5e25fSSatish Balay        0,
137597304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ,
137649b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
137749b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
137849b5e25fSSatish Balay        0,
137949b5e25fSSatish Balay        0,
138097304618SKris Buschelman /*55*/ 0,
138149b5e25fSSatish Balay        0,
138249b5e25fSSatish Balay        0,
138349b5e25fSSatish Balay        0,
138449b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
138597304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
138649b5e25fSSatish Balay        0,
138749b5e25fSSatish Balay        0,
13888a124369SBarry Smith        MatGetPetscMaps_Petsc,
1389d959ec07SHong Zhang        0,
139097304618SKris Buschelman /*65*/ 0,
1391d959ec07SHong Zhang        0,
1392d959ec07SHong Zhang        0,
1393d959ec07SHong Zhang        0,
1394d959ec07SHong Zhang        0,
139597304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
13963e0d88b5SBarry Smith        0,
13973e0d88b5SBarry Smith        0,
13983e0d88b5SBarry Smith        0,
13993e0d88b5SBarry Smith        0,
140097304618SKris Buschelman /*75*/ 0,
14013e0d88b5SBarry Smith        0,
14023e0d88b5SBarry Smith        0,
14033e0d88b5SBarry Smith        0,
14043e0d88b5SBarry Smith        0,
140597304618SKris Buschelman /*80*/ 0,
14063e0d88b5SBarry Smith        0,
14073e0d88b5SBarry Smith        0,
14083e0d88b5SBarry Smith        0,
14093e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
141097304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
14113e0d88b5SBarry Smith #else
141297304618SKris Buschelman        0,
14133e0d88b5SBarry Smith #endif
1414*efcf0fc3SBarry Smith /*85*/ MatLoad_SeqSBAIJ,
1415*efcf0fc3SBarry Smith        MatIsSymmetric_SeqSBAIJ,
1416*efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1417*efcf0fc3SBarry Smith        MatIsHermitian_SeqSBAIJ
14183e0d88b5SBarry Smith };
141949b5e25fSSatish Balay 
142049b5e25fSSatish Balay EXTERN_C_BEGIN
14214a2ae208SSatish Balay #undef __FUNCT__
14224a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
142349b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
142449b5e25fSSatish Balay {
14254afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1426c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
142749b5e25fSSatish Balay   int         ierr;
142849b5e25fSSatish Balay 
142949b5e25fSSatish Balay   PetscFunctionBegin;
143049b5e25fSSatish Balay   if (aij->nonew != 1) {
143129bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
143249b5e25fSSatish Balay   }
143349b5e25fSSatish Balay 
143449b5e25fSSatish Balay   /* allocate space for values if not already there */
143549b5e25fSSatish Balay   if (!aij->saved_values) {
143687828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
143749b5e25fSSatish Balay   }
143849b5e25fSSatish Balay 
143949b5e25fSSatish Balay   /* copy values over */
144087828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
144149b5e25fSSatish Balay   PetscFunctionReturn(0);
144249b5e25fSSatish Balay }
144349b5e25fSSatish Balay EXTERN_C_END
144449b5e25fSSatish Balay 
144549b5e25fSSatish Balay EXTERN_C_BEGIN
14464a2ae208SSatish Balay #undef __FUNCT__
14474a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
144849b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
144949b5e25fSSatish Balay {
14504afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1451c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
145249b5e25fSSatish Balay 
145349b5e25fSSatish Balay   PetscFunctionBegin;
145449b5e25fSSatish Balay   if (aij->nonew != 1) {
145529bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
145649b5e25fSSatish Balay   }
145749b5e25fSSatish Balay   if (!aij->saved_values) {
145829bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
145949b5e25fSSatish Balay   }
146049b5e25fSSatish Balay 
146149b5e25fSSatish Balay   /* copy values over */
146287828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
146349b5e25fSSatish Balay   PetscFunctionReturn(0);
146449b5e25fSSatish Balay }
146549b5e25fSSatish Balay EXTERN_C_END
146649b5e25fSSatish Balay 
14678549e402SHong Zhang EXTERN_C_BEGIN
14684a2ae208SSatish Balay #undef __FUNCT__
1469a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1470a23d5eceSKris Buschelman int MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz)
147149b5e25fSSatish Balay {
1472c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
14730c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
147449b5e25fSSatish Balay   PetscTruth   flg;
14754afc71dfSHong Zhang   int          s_nz;
147649b5e25fSSatish Balay 
147749b5e25fSSatish Balay   PetscFunctionBegin;
1478273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1479e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1480c464158bSHong Zhang   mbs  = B->m/bs;
148149b5e25fSSatish Balay   bs2  = bs*bs;
148249b5e25fSSatish Balay 
1483c464158bSHong Zhang   if (mbs*bs != B->m) {
148429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
148549b5e25fSSatish Balay   }
148649b5e25fSSatish Balay 
1487435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1488435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
148949b5e25fSSatish Balay   if (nnz) {
149049b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
149129bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
149280fe4e49SBarry 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);
149349b5e25fSSatish Balay     }
149449b5e25fSSatish Balay   }
149549b5e25fSSatish Balay 
1496e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
149749b5e25fSSatish Balay   if (!flg) {
149849b5e25fSSatish Balay     switch (bs) {
149949b5e25fSSatish Balay     case 1:
15004ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
150149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1502d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
150349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
150449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
150549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
150649b5e25fSSatish Balay       break;
150749b5e25fSSatish Balay     case 2:
15084ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
150949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
151049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
151149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
151249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
151349b5e25fSSatish Balay       break;
151449b5e25fSSatish Balay     case 3:
15155f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
151649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
151749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
151849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
151949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
152049b5e25fSSatish Balay       break;
152149b5e25fSSatish Balay     case 4:
15225f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
152349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
152449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
152549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
152649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
152749b5e25fSSatish Balay       break;
152849b5e25fSSatish Balay     case 5:
15295f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
153049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
153149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
153249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
153349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
153449b5e25fSSatish Balay       break;
153549b5e25fSSatish Balay     case 6:
15365f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
153749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
153849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
153949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
154049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
154149b5e25fSSatish Balay       break;
154249b5e25fSSatish Balay     case 7:
1543de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
154449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
154549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1546de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
154749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
154849b5e25fSSatish Balay       break;
154949b5e25fSSatish Balay     }
155049b5e25fSSatish Balay   }
155149b5e25fSSatish Balay 
155249b5e25fSSatish Balay   b->mbs = mbs;
15534afc71dfSHong Zhang   b->nbs = mbs;
1554b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
155549b5e25fSSatish Balay   if (!nnz) {
1556435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
155749b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
155849b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
15598cef66ccSHong Zhang       b->imax[i] = nz;
156049b5e25fSSatish Balay     }
1561153ea458SHong Zhang     nz = nz*mbs; /* total nz */
156249b5e25fSSatish Balay   } else {
156349b5e25fSSatish Balay     nz = 0;
15648cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
156549b5e25fSSatish Balay   }
15668cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
15678cef66ccSHong Zhang   s_nz = nz;
156849b5e25fSSatish Balay 
156949b5e25fSSatish Balay   /* allocate the matrix space */
1570c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
157182502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
157249b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
157349b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
157449b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
157549b5e25fSSatish Balay   b->i = b->j + s_nz;
157649b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
157749b5e25fSSatish Balay 
157849b5e25fSSatish Balay   /* pointer to beginning of each row */
157949b5e25fSSatish Balay   b->i[0] = 0;
158049b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
158149b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
158249b5e25fSSatish Balay   }
158349b5e25fSSatish Balay 
158449b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
15855033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1586b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
158749b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
158849b5e25fSSatish Balay 
158949b5e25fSSatish Balay   b->bs               = bs;
159049b5e25fSSatish Balay   b->bs2              = bs2;
159149b5e25fSSatish Balay   b->s_nz             = 0;
159249b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1593153ea458SHong Zhang 
159416cdd363SHong Zhang   b->inew             = 0;
159516cdd363SHong Zhang   b->jnew             = 0;
159616cdd363SHong Zhang   b->anew             = 0;
159716cdd363SHong Zhang   b->a2anew           = 0;
15981a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1599c464158bSHong Zhang   PetscFunctionReturn(0);
1600c464158bSHong Zhang }
1601a23d5eceSKris Buschelman EXTERN_C_END
1602153ea458SHong Zhang 
16030bad9183SKris Buschelman /*MC
1604fafad747SKris Buschelman    MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
16050bad9183SKris Buschelman    based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
16060bad9183SKris Buschelman 
16070bad9183SKris Buschelman    Options Database Keys:
16080bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
16090bad9183SKris Buschelman 
16100bad9183SKris Buschelman   Level: beginner
16110bad9183SKris Buschelman 
16120bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ
16130bad9183SKris Buschelman M*/
16140bad9183SKris Buschelman 
1615a23d5eceSKris Buschelman EXTERN_C_BEGIN
1616a23d5eceSKris Buschelman #undef __FUNCT__
1617a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1618a23d5eceSKris Buschelman int MatCreate_SeqSBAIJ(Mat B)
1619a23d5eceSKris Buschelman {
1620a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
1621a23d5eceSKris Buschelman   int          ierr,size;
1622a23d5eceSKris Buschelman 
1623a23d5eceSKris Buschelman   PetscFunctionBegin;
1624a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1625a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1626a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1627a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1628a23d5eceSKris Buschelman 
1629a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1630a23d5eceSKris Buschelman   B->data = (void*)b;
1631a23d5eceSKris Buschelman   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1632a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1633a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1634a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1635a23d5eceSKris Buschelman   B->factor           = 0;
1636a23d5eceSKris Buschelman   B->lupivotthreshold = 1.0;
1637a23d5eceSKris Buschelman   B->mapping          = 0;
1638a23d5eceSKris Buschelman   b->row              = 0;
1639a23d5eceSKris Buschelman   b->icol             = 0;
1640a23d5eceSKris Buschelman   b->reallocs         = 0;
1641a23d5eceSKris Buschelman   b->saved_values     = 0;
1642a23d5eceSKris Buschelman 
1643a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1644a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1645a23d5eceSKris Buschelman 
1646a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1647a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1648a23d5eceSKris Buschelman   b->nonew            = 0;
1649a23d5eceSKris Buschelman   b->diag             = 0;
1650a23d5eceSKris Buschelman   b->solve_work       = 0;
1651a23d5eceSKris Buschelman   b->mult_work        = 0;
1652a23d5eceSKris Buschelman   B->spptr            = 0;
1653a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1654a23d5eceSKris Buschelman   b->xtoy             = 0;
1655a23d5eceSKris Buschelman   b->XtoY             = 0;
1656a23d5eceSKris Buschelman 
1657a23d5eceSKris Buschelman   b->inew             = 0;
1658a23d5eceSKris Buschelman   b->jnew             = 0;
1659a23d5eceSKris Buschelman   b->anew             = 0;
1660a23d5eceSKris Buschelman   b->a2anew           = 0;
1661a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1662a23d5eceSKris Buschelman 
1663a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1664a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1665a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1666a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1667a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1668a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1669a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1670a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1671a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1672a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1673a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1674a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
1675a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1676a23d5eceSKris Buschelman }
1677a23d5eceSKris Buschelman EXTERN_C_END
1678a23d5eceSKris Buschelman 
1679a23d5eceSKris Buschelman #undef __FUNCT__
1680a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1681a23d5eceSKris Buschelman /*@C
1682a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1683a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1684a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1685a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1686a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1687a23d5eceSKris Buschelman 
1688a23d5eceSKris Buschelman    Collective on Mat
1689a23d5eceSKris Buschelman 
1690a23d5eceSKris Buschelman    Input Parameters:
1691a23d5eceSKris Buschelman +  A - the symmetric matrix
1692a23d5eceSKris Buschelman .  bs - size of block
1693a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1694a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1695a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1696a23d5eceSKris Buschelman 
1697a23d5eceSKris Buschelman    Options Database Keys:
1698a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1699a23d5eceSKris Buschelman                      block calculations (much slower)
1700a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1701a23d5eceSKris Buschelman 
1702a23d5eceSKris Buschelman    Level: intermediate
1703a23d5eceSKris Buschelman 
1704a23d5eceSKris Buschelman    Notes:
1705a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1706a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1707a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1708a23d5eceSKris Buschelman    matrices.
1709a23d5eceSKris Buschelman 
1710a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1711a23d5eceSKris Buschelman @*/
1712ca01db9bSBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) {
1713ca01db9bSBarry Smith   int ierr,(*f)(Mat,int,int,const int[]);
1714a23d5eceSKris Buschelman 
1715a23d5eceSKris Buschelman   PetscFunctionBegin;
1716a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1717a23d5eceSKris Buschelman   if (f) {
1718a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1719a23d5eceSKris Buschelman   }
1720a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1721a23d5eceSKris Buschelman }
172249b5e25fSSatish Balay 
17234a2ae208SSatish Balay #undef __FUNCT__
17244a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1725c464158bSHong Zhang /*@C
1726c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1727c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1728c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1729c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1730c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
173149b5e25fSSatish Balay 
1732c464158bSHong Zhang    Collective on MPI_Comm
1733c464158bSHong Zhang 
1734c464158bSHong Zhang    Input Parameters:
1735c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1736c464158bSHong Zhang .  bs - size of block
1737c464158bSHong Zhang .  m - number of rows, or number of columns
1738c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1739744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1740744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1741c464158bSHong Zhang 
1742c464158bSHong Zhang    Output Parameter:
1743c464158bSHong Zhang .  A - the symmetric matrix
1744c464158bSHong Zhang 
1745c464158bSHong Zhang    Options Database Keys:
1746c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1747c464158bSHong Zhang                      block calculations (much slower)
1748c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1749c464158bSHong Zhang 
1750c464158bSHong Zhang    Level: intermediate
1751c464158bSHong Zhang 
1752c464158bSHong Zhang    Notes:
1753c464158bSHong Zhang 
1754c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1755c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1756c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1757c464158bSHong Zhang    matrices.
1758c464158bSHong Zhang 
1759c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1760c464158bSHong Zhang @*/
1761ca01db9bSBarry Smith int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A)
1762c464158bSHong Zhang {
1763c464158bSHong Zhang   int ierr;
1764c464158bSHong Zhang 
1765c464158bSHong Zhang   PetscFunctionBegin;
1766c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1767c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1768273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
176949b5e25fSSatish Balay   PetscFunctionReturn(0);
177049b5e25fSSatish Balay }
177149b5e25fSSatish Balay 
17724a2ae208SSatish Balay #undef __FUNCT__
17734a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
177449b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
177549b5e25fSSatish Balay {
177649b5e25fSSatish Balay   Mat         C;
177749b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
177849b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
177949b5e25fSSatish Balay 
178049b5e25fSSatish Balay   PetscFunctionBegin;
178129bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
178249b5e25fSSatish Balay 
178349b5e25fSSatish Balay   *B = 0;
1784692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1785692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1786692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1787692f9cbeSHong Zhang 
178849b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1789273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
179049b5e25fSSatish Balay   C->factor         = A->factor;
179149b5e25fSSatish Balay   c->row            = 0;
179249b5e25fSSatish Balay   c->icol           = 0;
179349b5e25fSSatish Balay   c->saved_values   = 0;
179449b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
179549b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
179649b5e25fSSatish Balay 
179749b5e25fSSatish Balay   c->bs         = a->bs;
179849b5e25fSSatish Balay   c->bs2        = a->bs2;
179949b5e25fSSatish Balay   c->mbs        = a->mbs;
180049b5e25fSSatish Balay   c->nbs        = a->nbs;
180149b5e25fSSatish Balay 
1802b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1803b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
180449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
180549b5e25fSSatish Balay     c->imax[i] = a->imax[i];
180649b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
180749b5e25fSSatish Balay   }
180849b5e25fSSatish Balay 
180949b5e25fSSatish Balay   /* allocate the matrix space */
181049b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
181149b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
181282502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
181349b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
181449b5e25fSSatish Balay   c->i = c->j + nz;
181549b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
181649b5e25fSSatish Balay   if (mbs > 0) {
181749b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
181849b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
181949b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
182049b5e25fSSatish Balay     } else {
182149b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
182249b5e25fSSatish Balay     }
182349b5e25fSSatish Balay   }
182449b5e25fSSatish Balay 
1825b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
182649b5e25fSSatish Balay   c->sorted      = a->sorted;
182749b5e25fSSatish Balay   c->roworiented = a->roworiented;
182849b5e25fSSatish Balay   c->nonew       = a->nonew;
182949b5e25fSSatish Balay 
183049b5e25fSSatish Balay   if (a->diag) {
1831b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1832b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
183349b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
183449b5e25fSSatish Balay       c->diag[i] = a->diag[i];
183549b5e25fSSatish Balay     }
183649b5e25fSSatish Balay   } else c->diag        = 0;
183749b5e25fSSatish Balay   c->s_nz               = a->s_nz;
183849b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
183949b5e25fSSatish Balay   c->solve_work         = 0;
18402a1b7f2aSHong Zhang   C->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
184149b5e25fSSatish Balay   c->mult_work          = 0;
184249b5e25fSSatish Balay   *B = C;
1843b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
184449b5e25fSSatish Balay   PetscFunctionReturn(0);
184549b5e25fSSatish Balay }
184649b5e25fSSatish Balay 
18474a2ae208SSatish Balay #undef __FUNCT__
18484a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1849b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
185049b5e25fSSatish Balay {
185149b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
185249b5e25fSSatish Balay   Mat          B;
185349b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18543f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
185549b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
185649b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
185787828ca2SBarry Smith   PetscScalar  *aa;
185849b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
185949b5e25fSSatish Balay 
186049b5e25fSSatish Balay   PetscFunctionBegin;
1861b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
186249b5e25fSSatish Balay   bs2  = bs*bs;
186349b5e25fSSatish Balay 
186449b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
186529bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1866b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
186749b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1868552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
186949b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
187049b5e25fSSatish Balay 
187149b5e25fSSatish Balay   if (header[3] < 0) {
187229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
187349b5e25fSSatish Balay   }
187449b5e25fSSatish Balay 
187529bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
187649b5e25fSSatish Balay 
187749b5e25fSSatish Balay   /*
187849b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
187949b5e25fSSatish Balay     divisible by the blocksize
188049b5e25fSSatish Balay   */
188149b5e25fSSatish Balay   mbs        = M/bs;
188249b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
188349b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
188449b5e25fSSatish Balay   else                  mbs++;
188549b5e25fSSatish Balay   if (extra_rows) {
1886b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
188749b5e25fSSatish Balay   }
188849b5e25fSSatish Balay 
188949b5e25fSSatish Balay   /* read in row lengths */
1890b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
189149b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
189249b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
189349b5e25fSSatish Balay 
189449b5e25fSSatish Balay   /* read in column indices */
1895b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
189649b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
189749b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
189849b5e25fSSatish Balay 
189949b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
190082502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1901a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
190282502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
190349b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
190449b5e25fSSatish Balay   masked   = mask + mbs;
190549b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
190649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
190749b5e25fSSatish Balay     nmask = 0;
190849b5e25fSSatish Balay     for (j=0; j<bs; j++) {
190949b5e25fSSatish Balay       kmax = rowlengths[rowcount];
191049b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19112d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
191203630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
191349b5e25fSSatish Balay       }
191449b5e25fSSatish Balay       rowcount++;
191549b5e25fSSatish Balay     }
1916574b2666SHong Zhang     s_browlengths[i] += nmask;
1917574b2666SHong Zhang 
191849b5e25fSSatish Balay     /* zero out the mask elements we set */
191949b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
192049b5e25fSSatish Balay   }
192149b5e25fSSatish Balay 
192249b5e25fSSatish Balay   /* create our matrix */
19239abb65ffSKris Buschelman   ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr);
19249abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
192578473fd7SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr);
192649b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
192749b5e25fSSatish Balay 
192849b5e25fSSatish Balay   /* set matrix "i" values */
192949b5e25fSSatish Balay   a->i[0] = 0;
193049b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1931574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1932574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
193349b5e25fSSatish Balay   }
19347fe2be48SHong Zhang   a->s_nz = a->i[mbs];
193549b5e25fSSatish Balay 
193649b5e25fSSatish Balay   /* read in nonzero values */
193787828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
193849b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
193949b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
194049b5e25fSSatish Balay 
194149b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
194249b5e25fSSatish Balay   nzcount = 0; jcount = 0;
194349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
194449b5e25fSSatish Balay     nzcountb = nzcount;
194549b5e25fSSatish Balay     nmask    = 0;
194649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
194749b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
194849b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19492d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
195003630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19512d703238SHong Zhang       }
19522d703238SHong Zhang     }
19532d703238SHong Zhang     /* sort the masked values */
19542d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19552d703238SHong Zhang 
19562d703238SHong Zhang     /* set "j" values into matrix */
19572d703238SHong Zhang     maskcount = 1;
19582d703238SHong Zhang     for (j=0; j<nmask; j++) {
195949b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
196049b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
196149b5e25fSSatish Balay     }
1962574b2666SHong Zhang 
196349b5e25fSSatish Balay     /* set "a" values into matrix */
196449b5e25fSSatish Balay     ishift = bs2*a->i[i];
196549b5e25fSSatish Balay     for (j=0; j<bs; j++) {
196649b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
196749b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1968574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1969574b2666SHong Zhang         if (tmp >= i){
197049b5e25fSSatish Balay           block     = mask[tmp] - 1;
197149b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
197249b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1973574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1974574b2666SHong Zhang         }
1975574b2666SHong Zhang         nzcountb++;
197649b5e25fSSatish Balay       }
197749b5e25fSSatish Balay     }
197849b5e25fSSatish Balay     /* zero out the mask elements we set */
197949b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
198049b5e25fSSatish Balay   }
198129bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
198249b5e25fSSatish Balay 
198349b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1984574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
198549b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
198649b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
198749b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
198849b5e25fSSatish Balay 
19899abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19909abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199149b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19929abb65ffSKris Buschelman   *A = B;
199349b5e25fSSatish Balay   PetscFunctionReturn(0);
199449b5e25fSSatish Balay }
1995574b2666SHong Zhang 
1996d06b337dSHong Zhang #undef __FUNCT__
1997d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1998c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1999d06b337dSHong Zhang {
2000d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2001d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
2002d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
2003d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
2004d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
2005d06b337dSHong Zhang 
2006d06b337dSHong Zhang   PetscFunctionBegin;
2007b965ef7fSBarry Smith   its = its*lits;
200891723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
2009d06b337dSHong Zhang 
2010d06b337dSHong Zhang   if (bs > 1)
2011d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2012d06b337dSHong Zhang 
2013d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2014d06b337dSHong Zhang   if (xx != bb) {
2015d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
2016d06b337dSHong Zhang   } else {
2017d06b337dSHong Zhang     b = x;
2018d06b337dSHong Zhang   }
2019d06b337dSHong Zhang 
2020d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
2021d06b337dSHong Zhang 
2022d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
20232798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2024d06b337dSHong Zhang       for (i=0; i<m; i++)
2025d06b337dSHong Zhang         t[i] = b[i];
2026d06b337dSHong Zhang 
2027d06b337dSHong Zhang       for (i=0; i<m; i++){
202844706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2029d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2030d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2031d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2032d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2033d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
203444706875SHong Zhang         PetscLogFlops(2*nz-1);
2035d06b337dSHong Zhang       }
2036d06b337dSHong Zhang     }
2037d06b337dSHong Zhang 
20382798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2039d06b337dSHong Zhang       for (i=0; i<m; i++)
2040d06b337dSHong Zhang         t[i] = b[i];
2041d06b337dSHong Zhang 
2042d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2043d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2044d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2045d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2046d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
204744706875SHong Zhang         PetscLogFlops(2*nz-1);
2048d06b337dSHong Zhang       }
2049d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2050d06b337dSHong Zhang         d  = *(aa + ai[i]);
2051d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2052d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2053d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2054d06b337dSHong Zhang         sum = t[i];
2055d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
205644706875SHong Zhang         PetscLogFlops(2*nz-1);
2057d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2058d06b337dSHong Zhang       }
2059d06b337dSHong Zhang     }
2060d06b337dSHong Zhang     its--;
2061d06b337dSHong Zhang   }
2062d06b337dSHong Zhang 
2063d06b337dSHong Zhang   while (its--) {
206444706875SHong Zhang     /*
206544706875SHong Zhang        forward sweep:
206644706875SHong Zhang        for i=0,...,m-1:
206744706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
206844706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
206944706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2070d06b337dSHong Zhang 
207144706875SHong Zhang     */
20722798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2073d06b337dSHong Zhang       for (i=0; i<m; i++)
2074d06b337dSHong Zhang         t[i] = b[i];
2075d06b337dSHong Zhang 
2076d06b337dSHong Zhang       for (i=0; i<m; i++){
207744706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2078d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2079d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2080d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2081d06b337dSHong Zhang         sum = t[i];
2082d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2083d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2084d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
208544706875SHong Zhang         PetscLogFlops(4*nz-2);
2086d06b337dSHong Zhang       }
2087d06b337dSHong Zhang     }
2088d06b337dSHong Zhang 
20892798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
209044706875SHong Zhang       /*
209144706875SHong Zhang        backward sweep:
209244706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
209344706875SHong Zhang        for i=m-1,...,0:
209444706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
209544706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
209644706875SHong Zhang       */
2097d06b337dSHong Zhang       for (i=0; i<m; i++)
2098d06b337dSHong Zhang         t[i] = b[i];
2099d06b337dSHong Zhang 
2100d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2101d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2102d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2103d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2104d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
210544706875SHong Zhang         PetscLogFlops(2*nz-1);
2106d06b337dSHong Zhang       }
2107d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2108d06b337dSHong Zhang         d  = *(aa + ai[i]);
2109d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2110d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2111d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2112d06b337dSHong Zhang         sum = t[i];
2113d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
211444706875SHong Zhang         PetscLogFlops(2*nz-1);
2115d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2116d06b337dSHong Zhang       }
2117d06b337dSHong Zhang     }
2118d06b337dSHong Zhang   }
2119d06b337dSHong Zhang 
2120d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
2121d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2122d06b337dSHong Zhang   if (bb != xx) {
2123d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2124d06b337dSHong Zhang   }
21252798e883SHong Zhang 
2126d06b337dSHong Zhang   PetscFunctionReturn(0);
2127d06b337dSHong Zhang }
2128d06b337dSHong Zhang 
2129d06b337dSHong Zhang 
2130d06b337dSHong Zhang 
2131d06b337dSHong Zhang 
213249b5e25fSSatish Balay 
213349b5e25fSSatish Balay 
2134