xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 6c6c5352cdaa6aa3a8bd2967d4dd73a4e9ba1251)
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 */
62*6c6c5352SBarry Smith     int nz = a->i[n];
63*6c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
64a1373b80SHong Zhang     for (i=0; i<n+1; i++) a->i[i]++;
65a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
66a1373b80SHong Zhang   } else {
67a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
68a6ece127SHong Zhang   }
6949b5e25fSSatish Balay   PetscFunctionReturn(0);
7049b5e25fSSatish Balay }
7149b5e25fSSatish Balay 
724a2ae208SSatish Balay #undef __FUNCT__
734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
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) {
83*6c6c5352SBarry Smith     int nz = a->i[n]-1;
84*6c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
85a6ece127SHong Zhang     for (i=0; i<n+1; i++) a->i[i]--;
86a6ece127SHong Zhang   }
87a6ece127SHong Zhang   PetscFunctionReturn(0);
8849b5e25fSSatish Balay }
8949b5e25fSSatish Balay 
904a2ae208SSatish Balay #undef __FUNCT__
914a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ"
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)
110*6c6c5352SBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, NZ=%d",A->m,a->nz);
11149b5e25fSSatish Balay #endif
11249b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
11349b5e25fSSatish Balay   if (!a->singlemalloc) {
11449b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
11563c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
11649b5e25fSSatish Balay   }
11749b5e25fSSatish Balay   if (a->row) {
11849b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
11949b5e25fSSatish Balay   }
12049b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
12149b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
12249b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
12349b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
124d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
125d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
126d59c15a7SBarry Smith   if (a->mult_work)   {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
12749b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
128c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
1291a3463dfSHong Zhang 
1301a3463dfSHong Zhang   if (a->inew){
1311a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1321a3463dfSHong Zhang     a->inew = 0;
1331a3463dfSHong Zhang   }
13449b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13549b5e25fSSatish Balay   PetscFunctionReturn(0);
13649b5e25fSSatish Balay }
13749b5e25fSSatish Balay 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
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");
1829a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
1839a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1849a4540c5SBarry Smith   case MAT_HERMITIAN:
1859a4540c5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
18677e54ba9SKris Buschelman   case MAT_SYMMETRIC:
18777e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
1889a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
1899a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1909a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
19177e54ba9SKris Buschelman     break;
1924d9d31abSKris Buschelman   default:
19329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
19449b5e25fSSatish Balay   }
19549b5e25fSSatish Balay   PetscFunctionReturn(0);
19649b5e25fSSatish Balay }
19749b5e25fSSatish Balay 
1984a2ae208SSatish Balay #undef __FUNCT__
1994a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
20087828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
20149b5e25fSSatish Balay {
20249b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
20382502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
20449b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
20587828ca2SBarry Smith   PetscScalar  *v_i;
20649b5e25fSSatish Balay 
20749b5e25fSSatish Balay   PetscFunctionBegin;
20849b5e25fSSatish Balay   bs  = a->bs;
20949b5e25fSSatish Balay   ai  = a->i;
21049b5e25fSSatish Balay   aj  = a->j;
21149b5e25fSSatish Balay   aa  = a->a;
21249b5e25fSSatish Balay   bs2 = a->bs2;
21349b5e25fSSatish Balay 
214a45adfd6SMatthew Knepley   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %d out of range", row);
21549b5e25fSSatish Balay 
21649b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21749b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21849b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
21949b5e25fSSatish Balay   *ncols = bs*M;
22049b5e25fSSatish Balay 
22149b5e25fSSatish Balay   if (v) {
22249b5e25fSSatish Balay     *v = 0;
22349b5e25fSSatish Balay     if (*ncols) {
22487828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
22549b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22649b5e25fSSatish Balay         v_i  = *v + i*bs;
22749b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22849b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
22949b5e25fSSatish Balay       }
23049b5e25fSSatish Balay     }
23149b5e25fSSatish Balay   }
23249b5e25fSSatish Balay 
23349b5e25fSSatish Balay   if (cols) {
23449b5e25fSSatish Balay     *cols = 0;
23549b5e25fSSatish Balay     if (*ncols) {
23682502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
23749b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23849b5e25fSSatish Balay         cols_i = *cols + i*bs;
23949b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24049b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24149b5e25fSSatish Balay       }
24249b5e25fSSatish Balay     }
24349b5e25fSSatish Balay   }
24449b5e25fSSatish Balay 
24549b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2465ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2475ddb2528SHong Zhang #ifdef column_search
24849b5e25fSSatish Balay   v_i    = *v    + M*bs;
24949b5e25fSSatish Balay   cols_i = *cols + M*bs;
25049b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25149b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25249b5e25fSSatish Balay     for (j=0; j<M; j++){
25349b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25449b5e25fSSatish Balay       if (itmp == bn){
25549b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25649b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25749b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25849b5e25fSSatish Balay           *v_i++    = aa_i[k];
25949b5e25fSSatish Balay         }
26049b5e25fSSatish Balay         *ncols += bs;
26149b5e25fSSatish Balay         break;
26249b5e25fSSatish Balay       }
26349b5e25fSSatish Balay     }
26449b5e25fSSatish Balay   }
2655ddb2528SHong Zhang #endif
26649b5e25fSSatish Balay 
26749b5e25fSSatish Balay   PetscFunctionReturn(0);
26849b5e25fSSatish Balay }
26949b5e25fSSatish Balay 
2704a2ae208SSatish Balay #undef __FUNCT__
2714a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
27287828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
27349b5e25fSSatish Balay {
27449b5e25fSSatish Balay   int ierr;
27549b5e25fSSatish Balay 
27649b5e25fSSatish Balay   PetscFunctionBegin;
27749b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27849b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
27949b5e25fSSatish Balay   PetscFunctionReturn(0);
28049b5e25fSSatish Balay }
28149b5e25fSSatish Balay 
2824a2ae208SSatish Balay #undef __FUNCT__
2834a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
28449b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
28549b5e25fSSatish Balay {
2868115998fSBarry Smith   int ierr;
28749b5e25fSSatish Balay   PetscFunctionBegin;
288999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2898115998fSBarry Smith   PetscFunctionReturn(0);
29049b5e25fSSatish Balay }
29149b5e25fSSatish Balay 
2924a2ae208SSatish Balay #undef __FUNCT__
2934a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary"
294b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
29549b5e25fSSatish Balay {
29649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
29749b5e25fSSatish Balay   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
29887828ca2SBarry Smith   PetscScalar  *aa;
29949b5e25fSSatish Balay   FILE         *file;
30049b5e25fSSatish Balay 
30149b5e25fSSatish Balay   PetscFunctionBegin;
302b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
30382502324SSatish Balay   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
304552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
30549b5e25fSSatish Balay 
306c464158bSHong Zhang   col_lens[1] = A->m;
307c464158bSHong Zhang   col_lens[2] = A->m;
308*6c6c5352SBarry Smith   col_lens[3] = a->nz*bs2;
30949b5e25fSSatish Balay 
31049b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
31149b5e25fSSatish Balay   count = 0;
31249b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31449b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
31549b5e25fSSatish Balay     }
31649b5e25fSSatish Balay   }
317c464158bSHong Zhang   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
31849b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
31949b5e25fSSatish Balay 
32049b5e25fSSatish Balay   /* store column indices (zero start index) */
321*6c6c5352SBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
32249b5e25fSSatish Balay   count = 0;
32349b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
32449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
32549b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
32649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
32749b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
32849b5e25fSSatish Balay         }
32949b5e25fSSatish Balay       }
33049b5e25fSSatish Balay     }
33149b5e25fSSatish Balay   }
332*6c6c5352SBarry Smith   ierr = PetscBinaryWrite(fd,jj,bs2*a->nz,PETSC_INT,0);CHKERRQ(ierr);
33349b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
33449b5e25fSSatish Balay 
33549b5e25fSSatish Balay   /* store nonzero values */
336*6c6c5352SBarry Smith   ierr  = PetscMalloc((a->nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
33749b5e25fSSatish Balay   count = 0;
33849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
33949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
34049b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
34149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
34249b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
34349b5e25fSSatish Balay         }
34449b5e25fSSatish Balay       }
34549b5e25fSSatish Balay     }
34649b5e25fSSatish Balay   }
347*6c6c5352SBarry Smith   ierr = PetscBinaryWrite(fd,aa,bs2*a->nz,PETSC_SCALAR,0);CHKERRQ(ierr);
34849b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
34949b5e25fSSatish Balay 
350b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
35149b5e25fSSatish Balay   if (file) {
35249b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
35349b5e25fSSatish Balay   }
35449b5e25fSSatish Balay   PetscFunctionReturn(0);
35549b5e25fSSatish Balay }
35649b5e25fSSatish Balay 
3574a2ae208SSatish Balay #undef __FUNCT__
3584a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
359b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
36049b5e25fSSatish Balay {
36149b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
362fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
363fb9695e5SSatish Balay   char              *name;
364f3ef73ceSBarry Smith   PetscViewerFormat format;
36549b5e25fSSatish Balay 
36649b5e25fSSatish Balay   PetscFunctionBegin;
36780fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
368b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
369456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
370b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
371fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
37229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
373fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
374b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
37549b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
37649b5e25fSSatish Balay       for (j=0; j<bs; j++) {
377b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
37849b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
37949b5e25fSSatish Balay           for (l=0; l<bs; l++) {
38049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
38149b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
38349b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38449b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
38649b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38749b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38862b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38949b5e25fSSatish Balay             }
39049b5e25fSSatish Balay #else
39149b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
39262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
39349b5e25fSSatish Balay             }
39449b5e25fSSatish Balay #endif
39549b5e25fSSatish Balay           }
39649b5e25fSSatish Balay         }
397b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39849b5e25fSSatish Balay       }
39949b5e25fSSatish Balay     }
400b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
40149b5e25fSSatish Balay   } else {
402b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
40349b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
40449b5e25fSSatish Balay       for (j=0; j<bs; j++) {
405b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
40649b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
40749b5e25fSSatish Balay           for (l=0; l<bs; l++) {
40849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
40949b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
41062b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
41149b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41249b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
41362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
41449b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41549b5e25fSSatish Balay             } else {
41662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41749b5e25fSSatish Balay             }
41849b5e25fSSatish Balay #else
419b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
42049b5e25fSSatish Balay #endif
42149b5e25fSSatish Balay           }
42249b5e25fSSatish Balay         }
423b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
42449b5e25fSSatish Balay       }
42549b5e25fSSatish Balay     }
426b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
42749b5e25fSSatish Balay   }
428b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
42949b5e25fSSatish Balay   PetscFunctionReturn(0);
43049b5e25fSSatish Balay }
43149b5e25fSSatish Balay 
4324a2ae208SSatish Balay #undef __FUNCT__
4334a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
434b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
43549b5e25fSSatish Balay {
43649b5e25fSSatish Balay   Mat          A = (Mat) Aa;
43749b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
43849b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
43949b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
44049b5e25fSSatish Balay   MatScalar    *aa;
44149b5e25fSSatish Balay   MPI_Comm     comm;
442b0a32e0cSBarry Smith   PetscViewer  viewer;
44349b5e25fSSatish Balay 
44449b5e25fSSatish Balay   PetscFunctionBegin;
44549b5e25fSSatish Balay   /*
44649b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
44749b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
44849b5e25fSSatish Balay    rest should return immediately.
44949b5e25fSSatish Balay   */
45049b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
45149b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
45249b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
45349b5e25fSSatish Balay 
45449b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
45549b5e25fSSatish Balay 
456b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
457b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
45849b5e25fSSatish Balay 
45949b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
460b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
46149b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46249b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
463c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
46449b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46549b5e25fSSatish Balay       aa = a->a + j*bs2;
46649b5e25fSSatish Balay       for (k=0; k<bs; k++) {
46749b5e25fSSatish Balay         for (l=0; l<bs; l++) {
46849b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
469b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47049b5e25fSSatish Balay         }
47149b5e25fSSatish Balay       }
47249b5e25fSSatish Balay     }
47349b5e25fSSatish Balay   }
474b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
47549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
47649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
477c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
47849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
47949b5e25fSSatish Balay       aa = a->a + j*bs2;
48049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
48149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
48249b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
483b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
48449b5e25fSSatish Balay         }
48549b5e25fSSatish Balay       }
48649b5e25fSSatish Balay     }
48749b5e25fSSatish Balay   }
48849b5e25fSSatish Balay 
489b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
49049b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
49149b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
492c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
49349b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
49449b5e25fSSatish Balay       aa = a->a + j*bs2;
49549b5e25fSSatish Balay       for (k=0; k<bs; k++) {
49649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
49749b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
498b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
49949b5e25fSSatish Balay         }
50049b5e25fSSatish Balay       }
50149b5e25fSSatish Balay     }
50249b5e25fSSatish Balay   }
50349b5e25fSSatish Balay   PetscFunctionReturn(0);
50449b5e25fSSatish Balay }
50549b5e25fSSatish Balay 
5064a2ae208SSatish Balay #undef __FUNCT__
5074a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
508b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
50949b5e25fSSatish Balay {
51049b5e25fSSatish Balay   int          ierr;
51149b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
512b0a32e0cSBarry Smith   PetscDraw    draw;
51349b5e25fSSatish Balay   PetscTruth   isnull;
51449b5e25fSSatish Balay 
51549b5e25fSSatish Balay   PetscFunctionBegin;
51649b5e25fSSatish Balay 
517b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
518b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
51949b5e25fSSatish Balay 
52049b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
521c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
52249b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
523b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
524b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
52549b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
52649b5e25fSSatish Balay   PetscFunctionReturn(0);
52749b5e25fSSatish Balay }
52849b5e25fSSatish Balay 
5294a2ae208SSatish Balay #undef __FUNCT__
5304a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
531b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
53249b5e25fSSatish Balay {
53349b5e25fSSatish Balay   int        ierr;
53449b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
53549b5e25fSSatish Balay 
53649b5e25fSSatish Balay   PetscFunctionBegin;
537b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
538b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
539fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
540fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
54149b5e25fSSatish Balay   if (issocket) {
54229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
54349b5e25fSSatish Balay   } else if (isascii){
54449b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
54549b5e25fSSatish Balay   } else if (isbinary) {
54649b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
54749b5e25fSSatish Balay   } else if (isdraw) {
54849b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
54949b5e25fSSatish Balay   } else {
550b7aaefc3SHong Zhang     SETERRQ1(1,"Viewer type %s not supported by SeqSBAIJ matrices",((PetscObject)viewer)->type_name);
55149b5e25fSSatish Balay   }
55249b5e25fSSatish Balay   PetscFunctionReturn(0);
55349b5e25fSSatish Balay }
55449b5e25fSSatish Balay 
55549b5e25fSSatish Balay 
5564a2ae208SSatish Balay #undef __FUNCT__
5574a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
558f15d580aSBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[])
55949b5e25fSSatish Balay {
560045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
56149b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
56249b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
56349b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
56449b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
56549b5e25fSSatish Balay 
56649b5e25fSSatish Balay   PetscFunctionBegin;
56749b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
56849b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
569590ac198SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
570590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
57149b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
57249b5e25fSSatish Balay     nrow = ailen[brow];
57349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
574590ac198SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
575590ac198SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1);
57649b5e25fSSatish Balay       col  = in[l] ;
57749b5e25fSSatish Balay       bcol = col/bs;
57849b5e25fSSatish Balay       cidx = col%bs;
57949b5e25fSSatish Balay       ridx = row%bs;
58049b5e25fSSatish Balay       high = nrow;
58149b5e25fSSatish Balay       low  = 0; /* assume unsorted */
58249b5e25fSSatish Balay       while (high-low > 5) {
58349b5e25fSSatish Balay         t = (low+high)/2;
58449b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
58549b5e25fSSatish Balay         else             low  = t;
58649b5e25fSSatish Balay       }
58749b5e25fSSatish Balay       for (i=low; i<high; i++) {
58849b5e25fSSatish Balay         if (rp[i] > bcol) break;
58949b5e25fSSatish Balay         if (rp[i] == bcol) {
59049b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
59149b5e25fSSatish Balay           goto finished;
59249b5e25fSSatish Balay         }
59349b5e25fSSatish Balay       }
59449b5e25fSSatish Balay       *v++ = zero;
59549b5e25fSSatish Balay       finished:;
59649b5e25fSSatish Balay     }
59749b5e25fSSatish Balay   }
59849b5e25fSSatish Balay   PetscFunctionReturn(0);
59949b5e25fSSatish Balay }
60049b5e25fSSatish Balay 
60149b5e25fSSatish Balay 
6024a2ae208SSatish Balay #undef __FUNCT__
6034a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
604f15d580aSBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
60549b5e25fSSatish Balay {
6060880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
6070880e062SHong Zhang   int             *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
6080880e062SHong Zhang   int             *imax=a->imax,*ai=a->i,*ailen=a->ilen;
6090880e062SHong Zhang   int             *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
6100880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
611f15d580aSBarry Smith   const MatScalar *value = v;
612f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
6130880e062SHong Zhang 
61449b5e25fSSatish Balay   PetscFunctionBegin;
6150880e062SHong Zhang   if (roworiented) {
6160880e062SHong Zhang     stepval = (n-1)*bs;
6170880e062SHong Zhang   } else {
6180880e062SHong Zhang     stepval = (m-1)*bs;
6190880e062SHong Zhang   }
6200880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
6210880e062SHong Zhang     row  = im[k];
6220880e062SHong Zhang     if (row < 0) continue;
6230880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
624590ac198SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1);
6250880e062SHong Zhang #endif
6260880e062SHong Zhang     rp   = aj + ai[row];
6270880e062SHong Zhang     ap   = aa + bs2*ai[row];
6280880e062SHong Zhang     rmax = imax[row];
6290880e062SHong Zhang     nrow = ailen[row];
6300880e062SHong Zhang     low  = 0;
6310880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
6320880e062SHong Zhang       if (in[l] < 0) continue;
6330880e062SHong Zhang       col = in[l];
634b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g)
635b1823623SSatish Balay       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",col,a->nbs-1);
636b1823623SSatish Balay #endif
6370880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
6380880e062SHong Zhang       if (roworiented) {
6390880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
6400880e062SHong Zhang       } else {
6410880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
6420880e062SHong Zhang       }
6430880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
6440880e062SHong Zhang       while (high-low > 7) {
6450880e062SHong Zhang         t = (low+high)/2;
6460880e062SHong Zhang         if (rp[t] > col) high = t;
6470880e062SHong Zhang         else             low  = t;
6480880e062SHong Zhang       }
6490880e062SHong Zhang       for (i=low; i<high; i++) {
6500880e062SHong Zhang         if (rp[i] > col) break;
6510880e062SHong Zhang         if (rp[i] == col) {
6520880e062SHong Zhang           bap  = ap +  bs2*i;
6530880e062SHong Zhang           if (roworiented) {
6540880e062SHong Zhang             if (is == ADD_VALUES) {
6550880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6560880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6570880e062SHong Zhang                   bap[jj] += *value++;
6580880e062SHong Zhang                 }
6590880e062SHong Zhang               }
6600880e062SHong Zhang             } else {
6610880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6620880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6630880e062SHong Zhang                   bap[jj] = *value++;
6640880e062SHong Zhang                 }
6650880e062SHong Zhang               }
6660880e062SHong Zhang             }
6670880e062SHong Zhang           } else {
6680880e062SHong Zhang             if (is == ADD_VALUES) {
6690880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6700880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6710880e062SHong Zhang                   *bap++ += *value++;
6720880e062SHong Zhang                 }
6730880e062SHong Zhang               }
6740880e062SHong Zhang             } else {
6750880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6760880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6770880e062SHong Zhang                   *bap++  = *value++;
6780880e062SHong Zhang                 }
6790880e062SHong Zhang               }
6800880e062SHong Zhang             }
6810880e062SHong Zhang           }
6820880e062SHong Zhang           goto noinsert2;
6830880e062SHong Zhang         }
6840880e062SHong Zhang       }
6850880e062SHong Zhang       if (nonew == 1) goto noinsert2;
686a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6870880e062SHong Zhang       if (nrow >= rmax) {
6880880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
6890880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
6900880e062SHong Zhang         MatScalar *new_a;
6910880e062SHong Zhang 
692a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6930880e062SHong Zhang 
6940880e062SHong Zhang         /* malloc new storage space */
6950880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
6960880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
6970880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
6980880e062SHong Zhang         new_i   = new_j + new_nz;
6990880e062SHong Zhang 
7000880e062SHong Zhang         /* copy over old data into new slots */
7010880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
7020880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
7030880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
7040880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
7050880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
7060880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
7070880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
7080880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
7090880e062SHong Zhang         /* free up old matrix storage */
7100880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
7110880e062SHong Zhang         if (!a->singlemalloc) {
7120880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
7130880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
7140880e062SHong Zhang         }
7150880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
7160880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
7170880e062SHong Zhang 
7180880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
7190880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
7200880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
721*6c6c5352SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
7220880e062SHong Zhang         a->reallocs++;
723*6c6c5352SBarry Smith         a->nz++;
7240880e062SHong Zhang       }
7250880e062SHong Zhang       N = nrow++ - 1;
7260880e062SHong Zhang       /* shift up all the later entries in this row */
7270880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
7280880e062SHong Zhang         rp[ii+1] = rp[ii];
7290880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
7300880e062SHong Zhang       }
7310880e062SHong Zhang       if (N >= i) {
7320880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
7330880e062SHong Zhang       }
7340880e062SHong Zhang       rp[i] = col;
7350880e062SHong Zhang       bap   = ap +  bs2*i;
7360880e062SHong Zhang       if (roworiented) {
7370880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7380880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
7390880e062SHong Zhang             bap[jj] = *value++;
7400880e062SHong Zhang           }
7410880e062SHong Zhang         }
7420880e062SHong Zhang       } else {
7430880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7440880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
7450880e062SHong Zhang             *bap++  = *value++;
7460880e062SHong Zhang           }
7470880e062SHong Zhang         }
7480880e062SHong Zhang       }
7490880e062SHong Zhang       noinsert2:;
7500880e062SHong Zhang       low = i;
7510880e062SHong Zhang     }
7520880e062SHong Zhang     ailen[row] = nrow;
7530880e062SHong Zhang   }
7540880e062SHong Zhang   PetscFunctionReturn(0);
75549b5e25fSSatish Balay }
75649b5e25fSSatish Balay 
7574a2ae208SSatish Balay #undef __FUNCT__
7584a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
75949b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
76049b5e25fSSatish Balay {
76149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
76249b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
763c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
76449b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
76549b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
76649b5e25fSSatish Balay 
76749b5e25fSSatish Balay   PetscFunctionBegin;
76849b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
76949b5e25fSSatish Balay 
77049b5e25fSSatish Balay   if (m) rmax = ailen[0];
77149b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
77249b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
77349b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
77449b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
77549b5e25fSSatish Balay     if (fshift) {
77649b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
77749b5e25fSSatish Balay       N = ailen[i];
77849b5e25fSSatish Balay       for (j=0; j<N; j++) {
77949b5e25fSSatish Balay         ip[j-fshift] = ip[j];
78049b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
78149b5e25fSSatish Balay       }
78249b5e25fSSatish Balay     }
78349b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
78449b5e25fSSatish Balay   }
78549b5e25fSSatish Balay   if (mbs) {
78649b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
78749b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
78849b5e25fSSatish Balay   }
78949b5e25fSSatish Balay   /* reset ilen and imax for each row */
79049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
79149b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
79249b5e25fSSatish Balay   }
793*6c6c5352SBarry Smith   a->nz = ai[mbs];
79449b5e25fSSatish Balay 
795b424e231SHong Zhang   /* diagonals may have moved, reset it */
796b424e231SHong Zhang   if (a->diag) {
797b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
79849b5e25fSSatish Balay   }
799b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
800*6c6c5352SBarry Smith            m,A->m,a->bs,fshift*bs2,a->nz*bs2);
801b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
80249b5e25fSSatish Balay            a->reallocs);
803b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
80449b5e25fSSatish Balay   a->reallocs          = 0;
80549b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
80649b5e25fSSatish Balay 
80749b5e25fSSatish Balay   PetscFunctionReturn(0);
80849b5e25fSSatish Balay }
80949b5e25fSSatish Balay 
81049b5e25fSSatish Balay /*
81149b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
81249b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
81349b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
81449b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
81549b5e25fSSatish Balay */
8164a2ae208SSatish Balay #undef __FUNCT__
8174a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
818db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
81949b5e25fSSatish Balay {
82049b5e25fSSatish Balay   int        i,j,k,row;
82149b5e25fSSatish Balay   PetscTruth flg;
82249b5e25fSSatish Balay 
82349b5e25fSSatish Balay   PetscFunctionBegin;
82449b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
82549b5e25fSSatish Balay     row = idx[i];
82649b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
82749b5e25fSSatish Balay       sizes[j] = 1;
82849b5e25fSSatish Balay       i++;
82949b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
83049b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
83149b5e25fSSatish Balay       i++;
83249b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
83349b5e25fSSatish Balay       flg = PETSC_TRUE;
83449b5e25fSSatish Balay       for (k=1; k<bs; k++) {
83549b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
83649b5e25fSSatish Balay           flg = PETSC_FALSE;
83749b5e25fSSatish Balay           break;
83849b5e25fSSatish Balay         }
83949b5e25fSSatish Balay       }
84049b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
84149b5e25fSSatish Balay         sizes[j] = bs;
84249b5e25fSSatish Balay         i+= bs;
84349b5e25fSSatish Balay       } else {
84449b5e25fSSatish Balay         sizes[j] = 1;
84549b5e25fSSatish Balay         i++;
84649b5e25fSSatish Balay       }
84749b5e25fSSatish Balay     }
84849b5e25fSSatish Balay   }
84949b5e25fSSatish Balay   *bs_max = j;
85049b5e25fSSatish Balay   PetscFunctionReturn(0);
85149b5e25fSSatish Balay }
85249b5e25fSSatish Balay 
8534a2ae208SSatish Balay #undef __FUNCT__
8544a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
855268466fbSBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag)
85649b5e25fSSatish Balay {
85749b5e25fSSatish Balay   PetscFunctionBegin;
858c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
85949b5e25fSSatish Balay }
86049b5e25fSSatish Balay 
86149b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
86249b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
86349b5e25fSSatish Balay */
86449b5e25fSSatish Balay 
8654a2ae208SSatish Balay #undef __FUNCT__
8664a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
867f15d580aSBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
86849b5e25fSSatish Balay {
86949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
87049b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
87149b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
87249b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
87349b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
87449b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
87549b5e25fSSatish Balay 
87649b5e25fSSatish Balay   PetscFunctionBegin;
87749b5e25fSSatish Balay 
87849b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
87949b5e25fSSatish Balay     row  = im[k];       /* row number */
88049b5e25fSSatish Balay     brow = row/bs;      /* block row number */
88149b5e25fSSatish Balay     if (row < 0) continue;
88249b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
883590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
88449b5e25fSSatish Balay #endif
88549b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
88649b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
88749b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
88849b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
88949b5e25fSSatish Balay     low  = 0;
89049b5e25fSSatish Balay 
89149b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
89249b5e25fSSatish Balay       if (in[l] < 0) continue;
89349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
894590ac198SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m-1);
89549b5e25fSSatish Balay #endif
89649b5e25fSSatish Balay       col = in[l];
89749b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
89849b5e25fSSatish Balay 
89949b5e25fSSatish Balay       if (brow <= bcol){
90049b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
9018549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
90249b5e25fSSatish Balay           /* element value a(k,l) */
90349b5e25fSSatish Balay           if (roworiented) {
90449b5e25fSSatish Balay             value = v[l + k*n];
90549b5e25fSSatish Balay           } else {
90649b5e25fSSatish Balay             value = v[k + l*m];
90749b5e25fSSatish Balay           }
90849b5e25fSSatish Balay 
90949b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
91049b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
91149b5e25fSSatish Balay           while (high-low > 7) {
91249b5e25fSSatish Balay             t = (low+high)/2;
91349b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
91449b5e25fSSatish Balay             else              low  = t;
91549b5e25fSSatish Balay           }
91649b5e25fSSatish Balay           for (i=low; i<high; i++) {
91749b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
91849b5e25fSSatish Balay             if (rp[i] > bcol) break;
91949b5e25fSSatish Balay             if (rp[i] == bcol) {
92049b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
92149b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
92249b5e25fSSatish Balay               else                  *bap  = value;
9238549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9248549e402SHong Zhang               if (brow == bcol && ridx < cidx){
9258549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
9268549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
9278549e402SHong Zhang                 else                  *bap  = value;
9288549e402SHong Zhang               }
92949b5e25fSSatish Balay               goto noinsert1;
93049b5e25fSSatish Balay             }
93149b5e25fSSatish Balay           }
93249b5e25fSSatish Balay 
93349b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
934a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
93549b5e25fSSatish Balay       if (nrow >= rmax) {
93649b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
93749b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
93849b5e25fSSatish Balay         MatScalar *new_a;
93949b5e25fSSatish Balay 
940a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
94149b5e25fSSatish Balay 
94249b5e25fSSatish Balay         /* Malloc new storage space */
94349b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
94482502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
94549b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
94649b5e25fSSatish Balay         new_i = new_j + new_nz;
94749b5e25fSSatish Balay 
94849b5e25fSSatish Balay         /* copy over old data into new slots */
94949b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
95049b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
95149b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
95249b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
95349b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
95449b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
95549b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
95649b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
95749b5e25fSSatish Balay         /* free up old matrix storage */
95849b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
95949b5e25fSSatish Balay         if (!a->singlemalloc) {
96049b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
96149b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
96249b5e25fSSatish Balay         }
96349b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
96449b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
96549b5e25fSSatish Balay 
96649b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
96749b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
968b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
969*6c6c5352SBarry Smith         a->maxnz += bs2*CHUNKSIZE;
97049b5e25fSSatish Balay         a->reallocs++;
971*6c6c5352SBarry Smith         a->nz++;
97249b5e25fSSatish Balay       }
97349b5e25fSSatish Balay 
97449b5e25fSSatish Balay       N = nrow++ - 1;
97549b5e25fSSatish Balay       /* shift up all the later entries in this row */
97649b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
97749b5e25fSSatish Balay         rp[ii+1] = rp[ii];
97849b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
97949b5e25fSSatish Balay       }
98049b5e25fSSatish Balay       if (N>=i) {
98149b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
98249b5e25fSSatish Balay       }
98349b5e25fSSatish Balay       rp[i]                      = bcol;
98449b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
98549b5e25fSSatish Balay       noinsert1:;
98649b5e25fSSatish Balay       low = i;
98749b5e25fSSatish Balay       /* } */
9888549e402SHong Zhang         }
98949b5e25fSSatish Balay       } /* end of if .. if..  */
99049b5e25fSSatish Balay     }                     /* end of loop over added columns */
99149b5e25fSSatish Balay     ailen[brow] = nrow;
99249b5e25fSSatish Balay   }                       /* end of loop over added rows */
99349b5e25fSSatish Balay 
99449b5e25fSSatish Balay   PetscFunctionReturn(0);
99549b5e25fSSatish Balay }
99649b5e25fSSatish Balay 
99715e8a5b3SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
99815e8a5b3SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
99952ebccd6SSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int);
100049b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
1001268466fbSBarry Smith extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
100249b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
100349b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
1004268466fbSBarry Smith extern int MatScale_SeqSBAIJ(const PetscScalar*,Mat);
100549b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
100649b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
100749b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
100849b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
100949b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
101049b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
101113a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
101249b5e25fSSatish Balay 
101349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
101449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
101549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
101649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
101749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
101849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
101949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
102049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
102149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
102249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
102349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
102449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
102549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
102649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
102749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
102849b5e25fSSatish Balay 
1029d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1030d59c15a7SBarry Smith 
103107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
103207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
103307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
103407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
103507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
103607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
103707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
103807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
103907f98182SSatish Balay 
10405f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
10415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
10425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
10435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
10445f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
10455f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
10465f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
104757d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
10483e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
104949b5e25fSSatish Balay 
105049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
105149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
105249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
105349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
105449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
105549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
105649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
105749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
105849b5e25fSSatish Balay 
105949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
106049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
106149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
106249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
106349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
106449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
106549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
106649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
106749b5e25fSSatish Balay 
10684d101231SSatish Balay #ifdef HAVE_MatICCFactor
10691a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
10704a2ae208SSatish Balay #undef __FUNCT__
10714d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
107215e8a5b3SHong Zhang int MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
107349b5e25fSSatish Balay {
10744ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
107549b5e25fSSatish Balay   Mat         outA;
107649b5e25fSSatish Balay   int         ierr;
107749b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
107849b5e25fSSatish Balay 
107949b5e25fSSatish Balay   PetscFunctionBegin;
10801a3463dfSHong Zhang   /*
108129bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
108249b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
108349b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
108449b5e25fSSatish Balay   if (!row_identity || !col_identity) {
108529bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
108649b5e25fSSatish Balay   }
10871a3463dfSHong Zhang   */
108849b5e25fSSatish Balay 
108949b5e25fSSatish Balay   outA          = inA;
10901260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
109149b5e25fSSatish Balay 
109249b5e25fSSatish Balay   if (!a->diag) {
10931a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
109449b5e25fSSatish Balay   }
109549b5e25fSSatish Balay   /*
109649b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
109749b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
109849b5e25fSSatish Balay   */
109949b5e25fSSatish Balay   switch (a->bs) {
110049b5e25fSSatish Balay   case 1:
11010fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
110207f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1103d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
11044d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
110549b5e25fSSatish Balay   case 2:
11061a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
11071a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11081a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11094d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
111049b5e25fSSatish Balay     break;
111149b5e25fSSatish Balay   case 3:
11121a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
11131a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11141a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11154d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
111649b5e25fSSatish Balay     break;
111749b5e25fSSatish Balay   case 4:
11181a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
11191a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11201a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11214d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
112249b5e25fSSatish Balay     break;
112349b5e25fSSatish Balay   case 5:
11241a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
11251a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11261a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11274d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
112849b5e25fSSatish Balay     break;
112949b5e25fSSatish Balay   case 6:
11301a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
11311a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11321a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11334d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
113449b5e25fSSatish Balay     break;
113549b5e25fSSatish Balay   case 7:
11361a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
11371a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11381a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11394d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
114049b5e25fSSatish Balay     break;
114149b5e25fSSatish Balay   default:
114249b5e25fSSatish Balay     a->row        = row;
11431a3463dfSHong Zhang     a->icol       = col;
114449b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
114549b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
114649b5e25fSSatish Balay 
114749b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
114849b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1149b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
115049b5e25fSSatish Balay 
115149b5e25fSSatish Balay     if (!a->solve_work) {
115287828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
115387828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
115449b5e25fSSatish Balay     }
115549b5e25fSSatish Balay   }
115649b5e25fSSatish Balay 
11571a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
115849b5e25fSSatish Balay 
115949b5e25fSSatish Balay   PetscFunctionReturn(0);
116049b5e25fSSatish Balay }
1161950f1e5bSHong Zhang #endif
1162950f1e5bSHong Zhang 
11634a2ae208SSatish Balay #undef __FUNCT__
11644a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
116549b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
116649b5e25fSSatish Balay {
116749b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
116849b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
116949b5e25fSSatish Balay   int               ierr;
117049b5e25fSSatish Balay 
117149b5e25fSSatish Balay   PetscFunctionBegin;
117249b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
117349b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
117449b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
117549b5e25fSSatish Balay   PetscFunctionReturn(0);
117649b5e25fSSatish Balay }
117749b5e25fSSatish Balay 
117849b5e25fSSatish Balay EXTERN_C_BEGIN
11794a2ae208SSatish Balay #undef __FUNCT__
11804a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
118149b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
118249b5e25fSSatish Balay {
1183045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
118449b5e25fSSatish Balay   int         i,nz,n;
118549b5e25fSSatish Balay 
118649b5e25fSSatish Balay   PetscFunctionBegin;
1187*6c6c5352SBarry Smith   nz = baij->maxnz;
1188c464158bSHong Zhang   n  = mat->n;
118949b5e25fSSatish Balay   for (i=0; i<nz; i++) {
119049b5e25fSSatish Balay     baij->j[i] = indices[i];
119149b5e25fSSatish Balay   }
1192*6c6c5352SBarry Smith   baij->nz = nz;
119349b5e25fSSatish Balay   for (i=0; i<n; i++) {
119449b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
119549b5e25fSSatish Balay   }
119649b5e25fSSatish Balay 
119749b5e25fSSatish Balay   PetscFunctionReturn(0);
119849b5e25fSSatish Balay }
119949b5e25fSSatish Balay EXTERN_C_END
120049b5e25fSSatish Balay 
12014a2ae208SSatish Balay #undef __FUNCT__
12024a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
120349b5e25fSSatish Balay /*@
120419585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
120549b5e25fSSatish Balay        in the matrix.
120649b5e25fSSatish Balay 
120749b5e25fSSatish Balay   Input Parameters:
120819585528SSatish Balay +  mat     - the SeqSBAIJ matrix
120949b5e25fSSatish Balay -  indices - the column indices
121049b5e25fSSatish Balay 
121149b5e25fSSatish Balay   Level: advanced
121249b5e25fSSatish Balay 
121349b5e25fSSatish Balay   Notes:
121449b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
121549b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
121649b5e25fSSatish Balay   of the MatSetValues() operation.
121749b5e25fSSatish Balay 
121849b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
121919585528SSatish Balay   MatCreateSeqSBAIJ().
122049b5e25fSSatish Balay 
1221ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
122249b5e25fSSatish Balay 
1223ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
122449b5e25fSSatish Balay @*/
122549b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
122649b5e25fSSatish Balay {
122749b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
122849b5e25fSSatish Balay 
122949b5e25fSSatish Balay   PetscFunctionBegin;
123049b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1231c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
123249b5e25fSSatish Balay   if (f) {
123349b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
123449b5e25fSSatish Balay   } else {
123529bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
123649b5e25fSSatish Balay   }
123749b5e25fSSatish Balay   PetscFunctionReturn(0);
123849b5e25fSSatish Balay }
123949b5e25fSSatish Balay 
12404a2ae208SSatish Balay #undef __FUNCT__
12414a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1242273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1243273d9f13SBarry Smith {
1244273d9f13SBarry Smith   int        ierr;
1245273d9f13SBarry Smith 
1246273d9f13SBarry Smith   PetscFunctionBegin;
1247273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1248273d9f13SBarry Smith   PetscFunctionReturn(0);
1249273d9f13SBarry Smith }
1250273d9f13SBarry Smith 
1251a6ece127SHong Zhang #undef __FUNCT__
1252a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
12534e7234bfSBarry Smith int MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1254a6ece127SHong Zhang {
1255a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1256a6ece127SHong Zhang   PetscFunctionBegin;
1257a6ece127SHong Zhang   *array = a->a;
1258a6ece127SHong Zhang   PetscFunctionReturn(0);
1259a6ece127SHong Zhang }
1260a6ece127SHong Zhang 
1261a6ece127SHong Zhang #undef __FUNCT__
1262a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
12634e7234bfSBarry Smith int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1264a6ece127SHong Zhang {
1265a6ece127SHong Zhang   PetscFunctionBegin;
1266a6ece127SHong Zhang   PetscFunctionReturn(0);
1267a6ece127SHong Zhang }
1268a6ece127SHong Zhang 
126942ee4b1aSHong Zhang #include "petscblaslapack.h"
127042ee4b1aSHong Zhang #undef __FUNCT__
127142ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1272268466fbSBarry Smith int MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
127342ee4b1aSHong Zhang {
127442ee4b1aSHong Zhang   Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1275c4319e64SHong Zhang   int          ierr,one=1,i,bs=y->bs,bs2,j;
127642ee4b1aSHong Zhang 
127742ee4b1aSHong Zhang   PetscFunctionBegin;
127842ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1279*6c6c5352SBarry Smith     BLaxpy_(&x->nz,(PetscScalar*)a,x->a,&one,y->a,&one);
1280c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1281c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1282c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1283c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1284c537a176SHong Zhang     }
1285c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1286c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1287c4319e64SHong Zhang       y->XtoY = X;
1288c537a176SHong Zhang     }
1289c4319e64SHong Zhang     bs2 = bs*bs;
1290*6c6c5352SBarry Smith     for (i=0; i<x->nz; i++) {
1291c4319e64SHong Zhang       j = 0;
1292c4319e64SHong Zhang       while (j < bs2){
12936550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1294c4319e64SHong Zhang         j++;
1295c537a176SHong Zhang       }
1296c4319e64SHong Zhang     }
1297*6c6c5352SBarry Smith     PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %d/%d = %g\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));
129842ee4b1aSHong Zhang   } else {
129942ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
130042ee4b1aSHong Zhang   }
130142ee4b1aSHong Zhang   PetscFunctionReturn(0);
130242ee4b1aSHong Zhang }
130342ee4b1aSHong Zhang 
1304efcf0fc3SBarry Smith #undef __FUNCT__
1305efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1306efcf0fc3SBarry Smith int MatIsSymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1307efcf0fc3SBarry Smith {
1308efcf0fc3SBarry Smith   PetscFunctionBegin;
1309efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1310efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1311efcf0fc3SBarry Smith }
1312efcf0fc3SBarry Smith 
1313efcf0fc3SBarry Smith #undef __FUNCT__
1314efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1315efcf0fc3SBarry Smith int MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1316efcf0fc3SBarry Smith {
1317efcf0fc3SBarry Smith   PetscFunctionBegin;
1318efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1319efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1320efcf0fc3SBarry Smith }
1321efcf0fc3SBarry Smith 
1322efcf0fc3SBarry Smith #undef __FUNCT__
1323efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1324efcf0fc3SBarry Smith int MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg)
1325efcf0fc3SBarry Smith {
1326efcf0fc3SBarry Smith   PetscFunctionBegin;
1327efcf0fc3SBarry Smith   *flg = PETSC_FALSE;
1328efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1329efcf0fc3SBarry Smith }
1330efcf0fc3SBarry Smith 
133149b5e25fSSatish Balay /* -------------------------------------------------------------------*/
133249b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
133349b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
133449b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
133549b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
133697304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
133749b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
133849b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
133949b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
134049b5e25fSSatish Balay        0,
134149b5e25fSSatish Balay        0,
134297304618SKris Buschelman /*10*/ 0,
134349b5e25fSSatish Balay        0,
13445f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1345d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
134649b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
134797304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
134849b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
134949b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
135049b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
135149b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
135297304618SKris Buschelman /*20*/ 0,
135349b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
135449b5e25fSSatish Balay        0,
135549b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
135649b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
135797304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ,
135849b5e25fSSatish Balay        0,
135949b5e25fSSatish Balay        0,
13605f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13615f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
136297304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1363c464158bSHong Zhang        0,
13644d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1365a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1366a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
136797304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
136849b5e25fSSatish Balay        0,
136949b5e25fSSatish Balay        0,
137049b5e25fSSatish Balay        0,
1371950f1e5bSHong Zhang        0,
137297304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
137349b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
137449b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
137549b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
137649b5e25fSSatish Balay        0,
137797304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
137849b5e25fSSatish Balay        MatScale_SeqSBAIJ,
137949b5e25fSSatish Balay        0,
138049b5e25fSSatish Balay        0,
138149b5e25fSSatish Balay        0,
138297304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ,
138349b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
138449b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
138549b5e25fSSatish Balay        0,
138649b5e25fSSatish Balay        0,
138797304618SKris Buschelman /*55*/ 0,
138849b5e25fSSatish Balay        0,
138949b5e25fSSatish Balay        0,
139049b5e25fSSatish Balay        0,
139149b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
139297304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
139349b5e25fSSatish Balay        0,
139449b5e25fSSatish Balay        0,
13958a124369SBarry Smith        MatGetPetscMaps_Petsc,
1396d959ec07SHong Zhang        0,
139797304618SKris Buschelman /*65*/ 0,
1398d959ec07SHong Zhang        0,
1399d959ec07SHong Zhang        0,
1400d959ec07SHong Zhang        0,
1401d959ec07SHong Zhang        0,
140297304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
14033e0d88b5SBarry Smith        0,
14043e0d88b5SBarry Smith        0,
14053e0d88b5SBarry Smith        0,
14063e0d88b5SBarry Smith        0,
140797304618SKris Buschelman /*75*/ 0,
14083e0d88b5SBarry Smith        0,
14093e0d88b5SBarry Smith        0,
14103e0d88b5SBarry Smith        0,
14113e0d88b5SBarry Smith        0,
141297304618SKris Buschelman /*80*/ 0,
14133e0d88b5SBarry Smith        0,
14143e0d88b5SBarry Smith        0,
14153e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
141697304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
14173e0d88b5SBarry Smith #else
141897304618SKris Buschelman        0,
14193e0d88b5SBarry Smith #endif
1420efcf0fc3SBarry Smith /*85*/ MatLoad_SeqSBAIJ,
1421efcf0fc3SBarry Smith        MatIsSymmetric_SeqSBAIJ,
1422efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1423efcf0fc3SBarry Smith        MatIsHermitian_SeqSBAIJ
14243e0d88b5SBarry Smith };
142549b5e25fSSatish Balay 
142649b5e25fSSatish Balay EXTERN_C_BEGIN
14274a2ae208SSatish Balay #undef __FUNCT__
14284a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
142949b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
143049b5e25fSSatish Balay {
14314afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1432c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
143349b5e25fSSatish Balay   int         ierr;
143449b5e25fSSatish Balay 
143549b5e25fSSatish Balay   PetscFunctionBegin;
143649b5e25fSSatish Balay   if (aij->nonew != 1) {
143729bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
143849b5e25fSSatish Balay   }
143949b5e25fSSatish Balay 
144049b5e25fSSatish Balay   /* allocate space for values if not already there */
144149b5e25fSSatish Balay   if (!aij->saved_values) {
144287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
144349b5e25fSSatish Balay   }
144449b5e25fSSatish Balay 
144549b5e25fSSatish Balay   /* copy values over */
144687828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
144749b5e25fSSatish Balay   PetscFunctionReturn(0);
144849b5e25fSSatish Balay }
144949b5e25fSSatish Balay EXTERN_C_END
145049b5e25fSSatish Balay 
145149b5e25fSSatish Balay EXTERN_C_BEGIN
14524a2ae208SSatish Balay #undef __FUNCT__
14534a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
145449b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
145549b5e25fSSatish Balay {
14564afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1457c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
145849b5e25fSSatish Balay 
145949b5e25fSSatish Balay   PetscFunctionBegin;
146049b5e25fSSatish Balay   if (aij->nonew != 1) {
146129bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
146249b5e25fSSatish Balay   }
146349b5e25fSSatish Balay   if (!aij->saved_values) {
146429bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
146549b5e25fSSatish Balay   }
146649b5e25fSSatish Balay 
146749b5e25fSSatish Balay   /* copy values over */
146887828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
146949b5e25fSSatish Balay   PetscFunctionReturn(0);
147049b5e25fSSatish Balay }
147149b5e25fSSatish Balay EXTERN_C_END
147249b5e25fSSatish Balay 
14738549e402SHong Zhang EXTERN_C_BEGIN
14744a2ae208SSatish Balay #undef __FUNCT__
1475a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1476a23d5eceSKris Buschelman int MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz)
147749b5e25fSSatish Balay {
1478c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
14790c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
148049b5e25fSSatish Balay   PetscTruth   flg;
148149b5e25fSSatish Balay 
148249b5e25fSSatish Balay   PetscFunctionBegin;
1483273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1484e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1485c464158bSHong Zhang   mbs  = B->m/bs;
148649b5e25fSSatish Balay   bs2  = bs*bs;
148749b5e25fSSatish Balay 
1488c464158bSHong Zhang   if (mbs*bs != B->m) {
148929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
149049b5e25fSSatish Balay   }
149149b5e25fSSatish Balay 
1492435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1493435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
149449b5e25fSSatish Balay   if (nnz) {
149549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
149629bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
149780fe4e49SBarry 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);
149849b5e25fSSatish Balay     }
149949b5e25fSSatish Balay   }
150049b5e25fSSatish Balay 
1501e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
150249b5e25fSSatish Balay   if (!flg) {
150349b5e25fSSatish Balay     switch (bs) {
150449b5e25fSSatish Balay     case 1:
15054ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
150649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1507d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
150849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
150949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
151049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
151149b5e25fSSatish Balay       break;
151249b5e25fSSatish Balay     case 2:
15134ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
151449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
151549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
151649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
151749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
151849b5e25fSSatish Balay       break;
151949b5e25fSSatish Balay     case 3:
15205f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
152149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
152249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
152349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
152449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
152549b5e25fSSatish Balay       break;
152649b5e25fSSatish Balay     case 4:
15275f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
152849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
152949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
153049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
153149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
153249b5e25fSSatish Balay       break;
153349b5e25fSSatish Balay     case 5:
15345f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
153549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
153649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
153749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
153849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
153949b5e25fSSatish Balay       break;
154049b5e25fSSatish Balay     case 6:
15415f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
154249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
154349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
154449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
154549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
154649b5e25fSSatish Balay       break;
154749b5e25fSSatish Balay     case 7:
1548de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
154949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
155049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1551de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
155249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
155349b5e25fSSatish Balay       break;
155449b5e25fSSatish Balay     }
155549b5e25fSSatish Balay   }
155649b5e25fSSatish Balay 
155749b5e25fSSatish Balay   b->mbs = mbs;
15584afc71dfSHong Zhang   b->nbs = mbs;
1559b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
156049b5e25fSSatish Balay   if (!nnz) {
1561435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
156249b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
156349b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
15648cef66ccSHong Zhang       b->imax[i] = nz;
156549b5e25fSSatish Balay     }
1566153ea458SHong Zhang     nz = nz*mbs; /* total nz */
156749b5e25fSSatish Balay   } else {
156849b5e25fSSatish Balay     nz = 0;
15698cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
157049b5e25fSSatish Balay   }
1571*6c6c5352SBarry Smith   /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
157249b5e25fSSatish Balay 
157349b5e25fSSatish Balay   /* allocate the matrix space */
1574*6c6c5352SBarry Smith   len  = nz*sizeof(int) + nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
157582502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
1576*6c6c5352SBarry Smith   ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1577*6c6c5352SBarry Smith   b->j = (int*)(b->a + nz*bs2);
1578*6c6c5352SBarry Smith   ierr = PetscMemzero(b->j,nz*sizeof(int));CHKERRQ(ierr);
1579*6c6c5352SBarry Smith   b->i = b->j + nz;
158049b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
158149b5e25fSSatish Balay 
158249b5e25fSSatish Balay   /* pointer to beginning of each row */
158349b5e25fSSatish Balay   b->i[0] = 0;
158449b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
158549b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
158649b5e25fSSatish Balay   }
158749b5e25fSSatish Balay 
158849b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
15895033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1590b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
159149b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
159249b5e25fSSatish Balay 
159349b5e25fSSatish Balay   b->bs               = bs;
159449b5e25fSSatish Balay   b->bs2              = bs2;
1595*6c6c5352SBarry Smith   b->nz             = 0;
1596*6c6c5352SBarry Smith   b->maxnz          = nz*bs2;
1597153ea458SHong Zhang 
159816cdd363SHong Zhang   b->inew             = 0;
159916cdd363SHong Zhang   b->jnew             = 0;
160016cdd363SHong Zhang   b->anew             = 0;
160116cdd363SHong Zhang   b->a2anew           = 0;
16021a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1603c464158bSHong Zhang   PetscFunctionReturn(0);
1604c464158bSHong Zhang }
1605a23d5eceSKris Buschelman EXTERN_C_END
1606153ea458SHong Zhang 
16070bad9183SKris Buschelman /*MC
1608fafad747SKris Buschelman    MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
16090bad9183SKris Buschelman    based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
16100bad9183SKris Buschelman 
16110bad9183SKris Buschelman    Options Database Keys:
16120bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
16130bad9183SKris Buschelman 
16140bad9183SKris Buschelman   Level: beginner
16150bad9183SKris Buschelman 
16160bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ
16170bad9183SKris Buschelman M*/
16180bad9183SKris Buschelman 
1619a23d5eceSKris Buschelman EXTERN_C_BEGIN
1620a23d5eceSKris Buschelman #undef __FUNCT__
1621a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1622a23d5eceSKris Buschelman int MatCreate_SeqSBAIJ(Mat B)
1623a23d5eceSKris Buschelman {
1624a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
1625a23d5eceSKris Buschelman   int          ierr,size;
1626a23d5eceSKris Buschelman 
1627a23d5eceSKris Buschelman   PetscFunctionBegin;
1628a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1629a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1630a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1631a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1632a23d5eceSKris Buschelman 
1633a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1634a23d5eceSKris Buschelman   B->data = (void*)b;
1635a23d5eceSKris Buschelman   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1636a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1637a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1638a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1639a23d5eceSKris Buschelman   B->factor           = 0;
1640a23d5eceSKris Buschelman   B->lupivotthreshold = 1.0;
1641a23d5eceSKris Buschelman   B->mapping          = 0;
1642a23d5eceSKris Buschelman   b->row              = 0;
1643a23d5eceSKris Buschelman   b->icol             = 0;
1644a23d5eceSKris Buschelman   b->reallocs         = 0;
1645a23d5eceSKris Buschelman   b->saved_values     = 0;
1646a23d5eceSKris Buschelman 
1647a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1648a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1649a23d5eceSKris Buschelman 
1650a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1651a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1652a23d5eceSKris Buschelman   b->nonew            = 0;
1653a23d5eceSKris Buschelman   b->diag             = 0;
1654a23d5eceSKris Buschelman   b->solve_work       = 0;
1655a23d5eceSKris Buschelman   b->mult_work        = 0;
1656a23d5eceSKris Buschelman   B->spptr            = 0;
1657a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1658a23d5eceSKris Buschelman   b->xtoy             = 0;
1659a23d5eceSKris Buschelman   b->XtoY             = 0;
1660a23d5eceSKris Buschelman 
1661a23d5eceSKris Buschelman   b->inew             = 0;
1662a23d5eceSKris Buschelman   b->jnew             = 0;
1663a23d5eceSKris Buschelman   b->anew             = 0;
1664a23d5eceSKris Buschelman   b->a2anew           = 0;
1665a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1666a23d5eceSKris Buschelman 
1667a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1668a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1669a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1670a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1671a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1672a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1673a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1674a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1675a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1676a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1677a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1678a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
1679a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1680a23d5eceSKris Buschelman }
1681a23d5eceSKris Buschelman EXTERN_C_END
1682a23d5eceSKris Buschelman 
1683a23d5eceSKris Buschelman #undef __FUNCT__
1684a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1685a23d5eceSKris Buschelman /*@C
1686a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1687a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1688a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1689a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1690a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1691a23d5eceSKris Buschelman 
1692a23d5eceSKris Buschelman    Collective on Mat
1693a23d5eceSKris Buschelman 
1694a23d5eceSKris Buschelman    Input Parameters:
1695a23d5eceSKris Buschelman +  A - the symmetric matrix
1696a23d5eceSKris Buschelman .  bs - size of block
1697a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1698a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1699a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1700a23d5eceSKris Buschelman 
1701a23d5eceSKris Buschelman    Options Database Keys:
1702a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1703a23d5eceSKris Buschelman                      block calculations (much slower)
1704a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1705a23d5eceSKris Buschelman 
1706a23d5eceSKris Buschelman    Level: intermediate
1707a23d5eceSKris Buschelman 
1708a23d5eceSKris Buschelman    Notes:
1709a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1710a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1711a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1712a23d5eceSKris Buschelman    matrices.
1713a23d5eceSKris Buschelman 
1714a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1715a23d5eceSKris Buschelman @*/
1716ca01db9bSBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) {
1717ca01db9bSBarry Smith   int ierr,(*f)(Mat,int,int,const int[]);
1718a23d5eceSKris Buschelman 
1719a23d5eceSKris Buschelman   PetscFunctionBegin;
1720a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1721a23d5eceSKris Buschelman   if (f) {
1722a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1723a23d5eceSKris Buschelman   }
1724a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1725a23d5eceSKris Buschelman }
172649b5e25fSSatish Balay 
17274a2ae208SSatish Balay #undef __FUNCT__
17284a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1729c464158bSHong Zhang /*@C
1730c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1731c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1732c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1733c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1734c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
173549b5e25fSSatish Balay 
1736c464158bSHong Zhang    Collective on MPI_Comm
1737c464158bSHong Zhang 
1738c464158bSHong Zhang    Input Parameters:
1739c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1740c464158bSHong Zhang .  bs - size of block
1741c464158bSHong Zhang .  m - number of rows, or number of columns
1742c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1743744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1744744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1745c464158bSHong Zhang 
1746c464158bSHong Zhang    Output Parameter:
1747c464158bSHong Zhang .  A - the symmetric matrix
1748c464158bSHong Zhang 
1749c464158bSHong Zhang    Options Database Keys:
1750c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1751c464158bSHong Zhang                      block calculations (much slower)
1752c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1753c464158bSHong Zhang 
1754c464158bSHong Zhang    Level: intermediate
1755c464158bSHong Zhang 
1756c464158bSHong Zhang    Notes:
1757c464158bSHong Zhang 
1758c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1759c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1760c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1761c464158bSHong Zhang    matrices.
1762c464158bSHong Zhang 
1763c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1764c464158bSHong Zhang @*/
1765ca01db9bSBarry Smith int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A)
1766c464158bSHong Zhang {
1767c464158bSHong Zhang   int ierr;
1768c464158bSHong Zhang 
1769c464158bSHong Zhang   PetscFunctionBegin;
1770c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1771c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1772273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
177349b5e25fSSatish Balay   PetscFunctionReturn(0);
177449b5e25fSSatish Balay }
177549b5e25fSSatish Balay 
17764a2ae208SSatish Balay #undef __FUNCT__
17774a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
177849b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
177949b5e25fSSatish Balay {
178049b5e25fSSatish Balay   Mat         C;
178149b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
1782*6c6c5352SBarry Smith   int         i,len,mbs = a->mbs,nz = a->nz,bs2 =a->bs2,ierr;
178349b5e25fSSatish Balay 
178449b5e25fSSatish Balay   PetscFunctionBegin;
178529bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
178649b5e25fSSatish Balay 
178749b5e25fSSatish Balay   *B = 0;
1788692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1789692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1790692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1791692f9cbeSHong Zhang 
179249b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1793273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
179449b5e25fSSatish Balay   C->factor         = A->factor;
179549b5e25fSSatish Balay   c->row            = 0;
179649b5e25fSSatish Balay   c->icol           = 0;
179749b5e25fSSatish Balay   c->saved_values   = 0;
179849b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
179949b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
180049b5e25fSSatish Balay 
180149b5e25fSSatish Balay   c->bs         = a->bs;
180249b5e25fSSatish Balay   c->bs2        = a->bs2;
180349b5e25fSSatish Balay   c->mbs        = a->mbs;
180449b5e25fSSatish Balay   c->nbs        = a->nbs;
180549b5e25fSSatish Balay 
1806b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1807b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
180849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
180949b5e25fSSatish Balay     c->imax[i] = a->imax[i];
181049b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
181149b5e25fSSatish Balay   }
181249b5e25fSSatish Balay 
181349b5e25fSSatish Balay   /* allocate the matrix space */
181449b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
181549b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
181682502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
181749b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
181849b5e25fSSatish Balay   c->i = c->j + nz;
181949b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
182049b5e25fSSatish Balay   if (mbs > 0) {
182149b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
182249b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
182349b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
182449b5e25fSSatish Balay     } else {
182549b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
182649b5e25fSSatish Balay     }
182749b5e25fSSatish Balay   }
182849b5e25fSSatish Balay 
1829b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
183049b5e25fSSatish Balay   c->sorted      = a->sorted;
183149b5e25fSSatish Balay   c->roworiented = a->roworiented;
183249b5e25fSSatish Balay   c->nonew       = a->nonew;
183349b5e25fSSatish Balay 
183449b5e25fSSatish Balay   if (a->diag) {
1835b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1836b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
183749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
183849b5e25fSSatish Balay       c->diag[i] = a->diag[i];
183949b5e25fSSatish Balay     }
184049b5e25fSSatish Balay   } else c->diag        = 0;
1841*6c6c5352SBarry Smith   c->nz               = a->nz;
1842*6c6c5352SBarry Smith   c->maxnz            = a->maxnz;
184349b5e25fSSatish Balay   c->solve_work         = 0;
18442a1b7f2aSHong Zhang   C->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
184549b5e25fSSatish Balay   c->mult_work          = 0;
184649b5e25fSSatish Balay   *B = C;
1847b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
184849b5e25fSSatish Balay   PetscFunctionReturn(0);
184949b5e25fSSatish Balay }
185049b5e25fSSatish Balay 
18514a2ae208SSatish Balay #undef __FUNCT__
18524a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
18538e9aea5cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A)
185449b5e25fSSatish Balay {
185549b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
185649b5e25fSSatish Balay   Mat          B;
185749b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18583f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
185949b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
186049b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
186187828ca2SBarry Smith   PetscScalar  *aa;
186249b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
186349b5e25fSSatish Balay 
186449b5e25fSSatish Balay   PetscFunctionBegin;
1865b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
186649b5e25fSSatish Balay   bs2  = bs*bs;
186749b5e25fSSatish Balay 
186849b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
186929bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1870b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
187149b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1872552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
187349b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
187449b5e25fSSatish Balay 
187549b5e25fSSatish Balay   if (header[3] < 0) {
187629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
187749b5e25fSSatish Balay   }
187849b5e25fSSatish Balay 
187929bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
188049b5e25fSSatish Balay 
188149b5e25fSSatish Balay   /*
188249b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
188349b5e25fSSatish Balay     divisible by the blocksize
188449b5e25fSSatish Balay   */
188549b5e25fSSatish Balay   mbs        = M/bs;
188649b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
188749b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
188849b5e25fSSatish Balay   else                  mbs++;
188949b5e25fSSatish Balay   if (extra_rows) {
1890b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
189149b5e25fSSatish Balay   }
189249b5e25fSSatish Balay 
189349b5e25fSSatish Balay   /* read in row lengths */
1894b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
189549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
189649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
189749b5e25fSSatish Balay 
189849b5e25fSSatish Balay   /* read in column indices */
1899b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
190049b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
190149b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
190249b5e25fSSatish Balay 
190349b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
190482502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1905a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
190682502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
190749b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
190849b5e25fSSatish Balay   masked   = mask + mbs;
190949b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
191049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
191149b5e25fSSatish Balay     nmask = 0;
191249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
191349b5e25fSSatish Balay       kmax = rowlengths[rowcount];
191449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19152d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
191603630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
191749b5e25fSSatish Balay       }
191849b5e25fSSatish Balay       rowcount++;
191949b5e25fSSatish Balay     }
1920574b2666SHong Zhang     s_browlengths[i] += nmask;
1921574b2666SHong Zhang 
192249b5e25fSSatish Balay     /* zero out the mask elements we set */
192349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
192449b5e25fSSatish Balay   }
192549b5e25fSSatish Balay 
192649b5e25fSSatish Balay   /* create our matrix */
19279abb65ffSKris Buschelman   ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr);
19289abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
192978473fd7SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr);
193049b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
193149b5e25fSSatish Balay 
193249b5e25fSSatish Balay   /* set matrix "i" values */
193349b5e25fSSatish Balay   a->i[0] = 0;
193449b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1935574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1936574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
193749b5e25fSSatish Balay   }
1938*6c6c5352SBarry Smith   a->nz = a->i[mbs];
193949b5e25fSSatish Balay 
194049b5e25fSSatish Balay   /* read in nonzero values */
194187828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
194249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
194349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
194449b5e25fSSatish Balay 
194549b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
194649b5e25fSSatish Balay   nzcount = 0; jcount = 0;
194749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
194849b5e25fSSatish Balay     nzcountb = nzcount;
194949b5e25fSSatish Balay     nmask    = 0;
195049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
195149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
195249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19532d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
195403630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19552d703238SHong Zhang       }
19562d703238SHong Zhang     }
19572d703238SHong Zhang     /* sort the masked values */
19582d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19592d703238SHong Zhang 
19602d703238SHong Zhang     /* set "j" values into matrix */
19612d703238SHong Zhang     maskcount = 1;
19622d703238SHong Zhang     for (j=0; j<nmask; j++) {
196349b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
196449b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
196549b5e25fSSatish Balay     }
1966574b2666SHong Zhang 
196749b5e25fSSatish Balay     /* set "a" values into matrix */
196849b5e25fSSatish Balay     ishift = bs2*a->i[i];
196949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
197049b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
197149b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1972574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1973574b2666SHong Zhang         if (tmp >= i){
197449b5e25fSSatish Balay           block     = mask[tmp] - 1;
197549b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
197649b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1977574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1978574b2666SHong Zhang         }
1979574b2666SHong Zhang         nzcountb++;
198049b5e25fSSatish Balay       }
198149b5e25fSSatish Balay     }
198249b5e25fSSatish Balay     /* zero out the mask elements we set */
198349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
198449b5e25fSSatish Balay   }
1985*6c6c5352SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
198649b5e25fSSatish Balay 
198749b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1988574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
198949b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
199049b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
199149b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
199249b5e25fSSatish Balay 
19939abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19949abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199549b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19969abb65ffSKris Buschelman   *A = B;
199749b5e25fSSatish Balay   PetscFunctionReturn(0);
199849b5e25fSSatish Balay }
1999574b2666SHong Zhang 
2000d06b337dSHong Zhang #undef __FUNCT__
2001d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
2002c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
2003d06b337dSHong Zhang {
2004d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2005d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
2006d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
2007d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
2008d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
2009d06b337dSHong Zhang 
2010d06b337dSHong Zhang   PetscFunctionBegin;
2011b965ef7fSBarry Smith   its = its*lits;
201291723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
2013d06b337dSHong Zhang 
2014d06b337dSHong Zhang   if (bs > 1)
2015d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2016d06b337dSHong Zhang 
2017d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2018d06b337dSHong Zhang   if (xx != bb) {
2019d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
2020d06b337dSHong Zhang   } else {
2021d06b337dSHong Zhang     b = x;
2022d06b337dSHong Zhang   }
2023d06b337dSHong Zhang 
2024d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
2025d06b337dSHong Zhang 
2026d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
20272798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2028d06b337dSHong Zhang       for (i=0; i<m; i++)
2029d06b337dSHong Zhang         t[i] = b[i];
2030d06b337dSHong Zhang 
2031d06b337dSHong Zhang       for (i=0; i<m; i++){
203244706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2033d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2034d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2035d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2036d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2037d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
203844706875SHong Zhang         PetscLogFlops(2*nz-1);
2039d06b337dSHong Zhang       }
2040d06b337dSHong Zhang     }
2041d06b337dSHong Zhang 
20422798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2043d06b337dSHong Zhang       for (i=0; i<m; i++)
2044d06b337dSHong Zhang         t[i] = b[i];
2045d06b337dSHong Zhang 
2046d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2047d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2048d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2049d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2050d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
205144706875SHong Zhang         PetscLogFlops(2*nz-1);
2052d06b337dSHong Zhang       }
2053d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2054d06b337dSHong Zhang         d  = *(aa + ai[i]);
2055d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2056d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2057d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2058d06b337dSHong Zhang         sum = t[i];
2059d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
206044706875SHong Zhang         PetscLogFlops(2*nz-1);
2061d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2062d06b337dSHong Zhang       }
2063d06b337dSHong Zhang     }
2064d06b337dSHong Zhang     its--;
2065d06b337dSHong Zhang   }
2066d06b337dSHong Zhang 
2067d06b337dSHong Zhang   while (its--) {
206844706875SHong Zhang     /*
206944706875SHong Zhang        forward sweep:
207044706875SHong Zhang        for i=0,...,m-1:
207144706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
207244706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
207344706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2074d06b337dSHong Zhang 
207544706875SHong Zhang     */
20762798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2077d06b337dSHong Zhang       for (i=0; i<m; i++)
2078d06b337dSHong Zhang         t[i] = b[i];
2079d06b337dSHong Zhang 
2080d06b337dSHong Zhang       for (i=0; i<m; i++){
208144706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2082d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2083d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2084d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2085d06b337dSHong Zhang         sum = t[i];
2086d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2087d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2088d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
208944706875SHong Zhang         PetscLogFlops(4*nz-2);
2090d06b337dSHong Zhang       }
2091d06b337dSHong Zhang     }
2092d06b337dSHong Zhang 
20932798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
209444706875SHong Zhang       /*
209544706875SHong Zhang        backward sweep:
209644706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
209744706875SHong Zhang        for i=m-1,...,0:
209844706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
209944706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
210044706875SHong Zhang       */
2101d06b337dSHong Zhang       for (i=0; i<m; i++)
2102d06b337dSHong Zhang         t[i] = b[i];
2103d06b337dSHong Zhang 
2104d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2105d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2106d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2107d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2108d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
210944706875SHong Zhang         PetscLogFlops(2*nz-1);
2110d06b337dSHong Zhang       }
2111d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2112d06b337dSHong Zhang         d  = *(aa + ai[i]);
2113d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2114d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2115d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2116d06b337dSHong Zhang         sum = t[i];
2117d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
211844706875SHong Zhang         PetscLogFlops(2*nz-1);
2119d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2120d06b337dSHong Zhang       }
2121d06b337dSHong Zhang     }
2122d06b337dSHong Zhang   }
2123d06b337dSHong Zhang 
2124d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
2125d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2126d06b337dSHong Zhang   if (bb != xx) {
2127d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2128d06b337dSHong Zhang   }
21292798e883SHong Zhang 
2130d06b337dSHong Zhang   PetscFunctionReturn(0);
2131d06b337dSHong Zhang }
2132d06b337dSHong Zhang 
2133d06b337dSHong Zhang 
2134d06b337dSHong Zhang 
2135d06b337dSHong Zhang 
213649b5e25fSSatish Balay 
213749b5e25fSSatish Balay 
2138