xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision a45adfd634f538af0df75033dc42d977b3d42dfd)
173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
4a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h"         /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/vec/vecimpl.h"
949b5e25fSSatish Balay #include "src/inline/spops.h"
10aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1149b5e25fSSatish Balay 
1249b5e25fSSatish Balay #define CHUNKSIZE  10
1349b5e25fSSatish Balay 
1449b5e25fSSatish Balay /*
1549b5e25fSSatish Balay      Checks for missing diagonals
1649b5e25fSSatish Balay */
174a2ae208SSatish Balay #undef __FUNCT__
184a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
1949b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A)
2049b5e25fSSatish Balay {
21045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2249b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2349b5e25fSSatish Balay 
2449b5e25fSSatish Balay   PetscFunctionBegin;
25045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2649b5e25fSSatish Balay   diag = a->diag;
2749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
28e11b0e14SHong Zhang     if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i);
2949b5e25fSSatish Balay   }
3049b5e25fSSatish Balay   PetscFunctionReturn(0);
3149b5e25fSSatish Balay }
3249b5e25fSSatish Balay 
334a2ae208SSatish Balay #undef __FUNCT__
344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
3549b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A)
3649b5e25fSSatish Balay {
37045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
38b424e231SHong Zhang   int          i,mbs = a->mbs,ierr;
3949b5e25fSSatish Balay 
4049b5e25fSSatish Balay   PetscFunctionBegin;
4149b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4249b5e25fSSatish Balay 
43b424e231SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr);
44b424e231SHong Zhang   PetscLogObjectMemory(A,(mbs+1)*sizeof(int));
45b424e231SHong Zhang   for (i=0; i<mbs; i++) a->diag[i] = a->i[i];
4649b5e25fSSatish Balay   PetscFunctionReturn(0);
4749b5e25fSSatish Balay }
4849b5e25fSSatish Balay 
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
514e7234bfSBarry Smith static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done)
5249b5e25fSSatish Balay {
53a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
549e37e433SSatish Balay   int         n = a->mbs,i;
5549b5e25fSSatish Balay 
5649b5e25fSSatish Balay   PetscFunctionBegin;
57d3e5a4abSHong Zhang   *nn = n;
58a1373b80SHong Zhang   if (!ia) PetscFunctionReturn(0);
59a1373b80SHong Zhang 
60a6ece127SHong Zhang   if (oshift == 1) {
61a6ece127SHong Zhang     /* temporarily add 1 to i and j indices */
62a6ece127SHong Zhang     int s_nz = a->i[n];
63a6ece127SHong Zhang     for (i=0; i<s_nz; i++) a->j[i]++;
64a1373b80SHong Zhang     for (i=0; i<n+1; i++) a->i[i]++;
65a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
66a1373b80SHong Zhang   } else {
67a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
68a6ece127SHong Zhang   }
6949b5e25fSSatish Balay   PetscFunctionReturn(0);
7049b5e25fSSatish Balay }
7149b5e25fSSatish Balay 
724a2ae208SSatish Balay #undef __FUNCT__
734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
744e7234bfSBarry Smith static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int *ia[],int *ja[],PetscTruth *done)
7549b5e25fSSatish Balay {
76b7aaefc3SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
77a15ac5e4SHong Zhang   int         i,n = a->mbs;
78a6ece127SHong Zhang 
7949b5e25fSSatish Balay   PetscFunctionBegin;
8049b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
81a6ece127SHong Zhang 
82a6ece127SHong Zhang   if (oshift == 1) {
83a6ece127SHong Zhang     int s_nz = a->i[n]-1;
84a6ece127SHong Zhang     for (i=0; i<s_nz; i++) a->j[i]--;
85a6ece127SHong Zhang     for (i=0; i<n+1; i++) a->i[i]--;
86a6ece127SHong Zhang   }
87a6ece127SHong Zhang   PetscFunctionReturn(0);
8849b5e25fSSatish Balay }
8949b5e25fSSatish Balay 
904a2ae208SSatish Balay #undef __FUNCT__
914a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ"
9249b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
9349b5e25fSSatish Balay {
9449b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
9549b5e25fSSatish Balay 
9649b5e25fSSatish Balay   PetscFunctionBegin;
9749b5e25fSSatish Balay   *bs = sbaij->bs;
9849b5e25fSSatish Balay   PetscFunctionReturn(0);
9949b5e25fSSatish Balay }
10049b5e25fSSatish Balay 
1014a2ae208SSatish Balay #undef __FUNCT__
1024a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
10349b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
10449b5e25fSSatish Balay {
10549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10649b5e25fSSatish Balay   int         ierr;
10749b5e25fSSatish Balay 
10849b5e25fSSatish Balay   PetscFunctionBegin;
10949b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
110b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz);
11149b5e25fSSatish Balay #endif
11249b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
11349b5e25fSSatish Balay   if (!a->singlemalloc) {
11449b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
11563c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
11649b5e25fSSatish Balay   }
11749b5e25fSSatish Balay   if (a->row) {
11849b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
11949b5e25fSSatish Balay   }
12049b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
12149b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
12249b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
12349b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
124d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
125d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
126d59c15a7SBarry Smith   if (a->mult_work)   {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
12749b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
128c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
1291a3463dfSHong Zhang 
1301a3463dfSHong Zhang   if (a->inew){
1311a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1321a3463dfSHong Zhang     a->inew = 0;
1331a3463dfSHong Zhang   }
13449b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13549b5e25fSSatish Balay   PetscFunctionReturn(0);
13649b5e25fSSatish Balay }
13749b5e25fSSatish Balay 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
14049b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14149b5e25fSSatish Balay {
142045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14349b5e25fSSatish Balay 
14449b5e25fSSatish Balay   PetscFunctionBegin;
1454d9d31abSKris Buschelman   switch (op) {
1464d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1474d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1484d9d31abSKris Buschelman     break;
1494d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1504d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1514d9d31abSKris Buschelman     break;
1524d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1534d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1564d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1594d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1624d9d31abSKris Buschelman     a->nonew          = 1;
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1654d9d31abSKris Buschelman     a->nonew          = -1;
1664d9d31abSKris Buschelman     break;
1674d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1684d9d31abSKris Buschelman     a->nonew          = -2;
1694d9d31abSKris Buschelman     break;
1704d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1714d9d31abSKris Buschelman     a->nonew          = 0;
1724d9d31abSKris Buschelman     break;
1734d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1744d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1754d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1764d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1774d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
178b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1794d9d31abSKris Buschelman     break;
1804d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
18129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1824d9d31abSKris Buschelman   default:
18329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
18449b5e25fSSatish Balay   }
18549b5e25fSSatish Balay   PetscFunctionReturn(0);
18649b5e25fSSatish Balay }
18749b5e25fSSatish Balay 
1884a2ae208SSatish Balay #undef __FUNCT__
1894a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19087828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
19149b5e25fSSatish Balay {
19249b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
19382502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
19449b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
19587828ca2SBarry Smith   PetscScalar  *v_i;
19649b5e25fSSatish Balay 
19749b5e25fSSatish Balay   PetscFunctionBegin;
19849b5e25fSSatish Balay   bs  = a->bs;
19949b5e25fSSatish Balay   ai  = a->i;
20049b5e25fSSatish Balay   aj  = a->j;
20149b5e25fSSatish Balay   aa  = a->a;
20249b5e25fSSatish Balay   bs2 = a->bs2;
20349b5e25fSSatish Balay 
204*a45adfd6SMatthew Knepley   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %d out of range", row);
20549b5e25fSSatish Balay 
20649b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
20749b5e25fSSatish Balay   bp  = row % bs; /* Block position */
20849b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
20949b5e25fSSatish Balay   *ncols = bs*M;
21049b5e25fSSatish Balay 
21149b5e25fSSatish Balay   if (v) {
21249b5e25fSSatish Balay     *v = 0;
21349b5e25fSSatish Balay     if (*ncols) {
21487828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
21549b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
21649b5e25fSSatish Balay         v_i  = *v + i*bs;
21749b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
21849b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
21949b5e25fSSatish Balay       }
22049b5e25fSSatish Balay     }
22149b5e25fSSatish Balay   }
22249b5e25fSSatish Balay 
22349b5e25fSSatish Balay   if (cols) {
22449b5e25fSSatish Balay     *cols = 0;
22549b5e25fSSatish Balay     if (*ncols) {
22682502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
22749b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22849b5e25fSSatish Balay         cols_i = *cols + i*bs;
22949b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
23049b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
23149b5e25fSSatish Balay       }
23249b5e25fSSatish Balay     }
23349b5e25fSSatish Balay   }
23449b5e25fSSatish Balay 
23549b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2365ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2375ddb2528SHong Zhang #ifdef column_search
23849b5e25fSSatish Balay   v_i    = *v    + M*bs;
23949b5e25fSSatish Balay   cols_i = *cols + M*bs;
24049b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
24149b5e25fSSatish Balay     M = ai[i+1] - ai[i];
24249b5e25fSSatish Balay     for (j=0; j<M; j++){
24349b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
24449b5e25fSSatish Balay       if (itmp == bn){
24549b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
24649b5e25fSSatish Balay         for (k=0; k<bs; k++) {
24749b5e25fSSatish Balay           *cols_i++ = i*bs+k;
24849b5e25fSSatish Balay           *v_i++    = aa_i[k];
24949b5e25fSSatish Balay         }
25049b5e25fSSatish Balay         *ncols += bs;
25149b5e25fSSatish Balay         break;
25249b5e25fSSatish Balay       }
25349b5e25fSSatish Balay     }
25449b5e25fSSatish Balay   }
2555ddb2528SHong Zhang #endif
25649b5e25fSSatish Balay 
25749b5e25fSSatish Balay   PetscFunctionReturn(0);
25849b5e25fSSatish Balay }
25949b5e25fSSatish Balay 
2604a2ae208SSatish Balay #undef __FUNCT__
2614a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
26287828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
26349b5e25fSSatish Balay {
26449b5e25fSSatish Balay   int ierr;
26549b5e25fSSatish Balay 
26649b5e25fSSatish Balay   PetscFunctionBegin;
26749b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
26849b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
26949b5e25fSSatish Balay   PetscFunctionReturn(0);
27049b5e25fSSatish Balay }
27149b5e25fSSatish Balay 
2724a2ae208SSatish Balay #undef __FUNCT__
2734a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
27449b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
27549b5e25fSSatish Balay {
2768115998fSBarry Smith   int ierr;
27749b5e25fSSatish Balay   PetscFunctionBegin;
278999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2798115998fSBarry Smith   PetscFunctionReturn(0);
28049b5e25fSSatish Balay }
28149b5e25fSSatish Balay 
2824a2ae208SSatish Balay #undef __FUNCT__
2834a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary"
284b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
28549b5e25fSSatish Balay {
28649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
28749b5e25fSSatish Balay   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
28887828ca2SBarry Smith   PetscScalar  *aa;
28949b5e25fSSatish Balay   FILE         *file;
29049b5e25fSSatish Balay 
29149b5e25fSSatish Balay   PetscFunctionBegin;
292b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
29382502324SSatish Balay   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
294552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
29549b5e25fSSatish Balay 
296c464158bSHong Zhang   col_lens[1] = A->m;
297c464158bSHong Zhang   col_lens[2] = A->m;
29849b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
29949b5e25fSSatish Balay 
30049b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
30149b5e25fSSatish Balay   count = 0;
30249b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
30349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
30449b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
30549b5e25fSSatish Balay     }
30649b5e25fSSatish Balay   }
307c464158bSHong Zhang   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
30849b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
30949b5e25fSSatish Balay 
31049b5e25fSSatish Balay   /* store column indices (zero start index) */
31182502324SSatish Balay   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
31249b5e25fSSatish Balay   count = 0;
31349b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31549b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
31649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
31749b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
31849b5e25fSSatish Balay         }
31949b5e25fSSatish Balay       }
32049b5e25fSSatish Balay     }
32149b5e25fSSatish Balay   }
32249b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
32349b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
32449b5e25fSSatish Balay 
32549b5e25fSSatish Balay   /* store nonzero values */
32687828ca2SBarry Smith   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
32749b5e25fSSatish Balay   count = 0;
32849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
32949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
33049b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
33149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
33249b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
33349b5e25fSSatish Balay         }
33449b5e25fSSatish Balay       }
33549b5e25fSSatish Balay     }
33649b5e25fSSatish Balay   }
33749b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
33849b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
33949b5e25fSSatish Balay 
340b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
34149b5e25fSSatish Balay   if (file) {
34249b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
34349b5e25fSSatish Balay   }
34449b5e25fSSatish Balay   PetscFunctionReturn(0);
34549b5e25fSSatish Balay }
34649b5e25fSSatish Balay 
3474a2ae208SSatish Balay #undef __FUNCT__
3484a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
349b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
35049b5e25fSSatish Balay {
35149b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
352fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
353fb9695e5SSatish Balay   char              *name;
354f3ef73ceSBarry Smith   PetscViewerFormat format;
35549b5e25fSSatish Balay 
35649b5e25fSSatish Balay   PetscFunctionBegin;
35780fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
358b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
359456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
360b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
361fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
36229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
363fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
364b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
36549b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
36649b5e25fSSatish Balay       for (j=0; j<bs; j++) {
367b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
36849b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
36949b5e25fSSatish Balay           for (l=0; l<bs; l++) {
37049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
37149b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
37349b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37449b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37562b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
37649b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37749b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37862b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37949b5e25fSSatish Balay             }
38049b5e25fSSatish Balay #else
38149b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
38262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
38349b5e25fSSatish Balay             }
38449b5e25fSSatish Balay #endif
38549b5e25fSSatish Balay           }
38649b5e25fSSatish Balay         }
387b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
38849b5e25fSSatish Balay       }
38949b5e25fSSatish Balay     }
390b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
39149b5e25fSSatish Balay   } else {
392b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
39349b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
39449b5e25fSSatish Balay       for (j=0; j<bs; j++) {
395b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
39649b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
39749b5e25fSSatish Balay           for (l=0; l<bs; l++) {
39849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
39949b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
40062b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
40149b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40249b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
40362b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
40449b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40549b5e25fSSatish Balay             } else {
40662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40749b5e25fSSatish Balay             }
40849b5e25fSSatish Balay #else
409b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
41049b5e25fSSatish Balay #endif
41149b5e25fSSatish Balay           }
41249b5e25fSSatish Balay         }
413b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
41449b5e25fSSatish Balay       }
41549b5e25fSSatish Balay     }
416b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
41749b5e25fSSatish Balay   }
418b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
41949b5e25fSSatish Balay   PetscFunctionReturn(0);
42049b5e25fSSatish Balay }
42149b5e25fSSatish Balay 
4224a2ae208SSatish Balay #undef __FUNCT__
4234a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
424b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
42549b5e25fSSatish Balay {
42649b5e25fSSatish Balay   Mat          A = (Mat) Aa;
42749b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
42849b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
42949b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
43049b5e25fSSatish Balay   MatScalar    *aa;
43149b5e25fSSatish Balay   MPI_Comm     comm;
432b0a32e0cSBarry Smith   PetscViewer  viewer;
43349b5e25fSSatish Balay 
43449b5e25fSSatish Balay   PetscFunctionBegin;
43549b5e25fSSatish Balay   /*
43649b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
43749b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
43849b5e25fSSatish Balay    rest should return immediately.
43949b5e25fSSatish Balay   */
44049b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
44149b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
44249b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
44349b5e25fSSatish Balay 
44449b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
44549b5e25fSSatish Balay 
446b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
447b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
44849b5e25fSSatish Balay 
44949b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
450b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
45149b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45249b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
453c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
45449b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
45549b5e25fSSatish Balay       aa = a->a + j*bs2;
45649b5e25fSSatish Balay       for (k=0; k<bs; k++) {
45749b5e25fSSatish Balay         for (l=0; l<bs; l++) {
45849b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
459b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
46049b5e25fSSatish Balay         }
46149b5e25fSSatish Balay       }
46249b5e25fSSatish Balay     }
46349b5e25fSSatish Balay   }
464b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
46549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
467c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
46849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46949b5e25fSSatish Balay       aa = a->a + j*bs2;
47049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47249b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
473b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47449b5e25fSSatish Balay         }
47549b5e25fSSatish Balay       }
47649b5e25fSSatish Balay     }
47749b5e25fSSatish Balay   }
47849b5e25fSSatish Balay 
479b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
48049b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
48149b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
482c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
48349b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
48449b5e25fSSatish Balay       aa = a->a + j*bs2;
48549b5e25fSSatish Balay       for (k=0; k<bs; k++) {
48649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
48749b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
488b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
48949b5e25fSSatish Balay         }
49049b5e25fSSatish Balay       }
49149b5e25fSSatish Balay     }
49249b5e25fSSatish Balay   }
49349b5e25fSSatish Balay   PetscFunctionReturn(0);
49449b5e25fSSatish Balay }
49549b5e25fSSatish Balay 
4964a2ae208SSatish Balay #undef __FUNCT__
4974a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
498b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
49949b5e25fSSatish Balay {
50049b5e25fSSatish Balay   int          ierr;
50149b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
502b0a32e0cSBarry Smith   PetscDraw    draw;
50349b5e25fSSatish Balay   PetscTruth   isnull;
50449b5e25fSSatish Balay 
50549b5e25fSSatish Balay   PetscFunctionBegin;
50649b5e25fSSatish Balay 
507b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
508b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
50949b5e25fSSatish Balay 
51049b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
511c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
51249b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
513b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
514b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
51549b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
51649b5e25fSSatish Balay   PetscFunctionReturn(0);
51749b5e25fSSatish Balay }
51849b5e25fSSatish Balay 
5194a2ae208SSatish Balay #undef __FUNCT__
5204a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
521b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
52249b5e25fSSatish Balay {
52349b5e25fSSatish Balay   int        ierr;
52449b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
52549b5e25fSSatish Balay 
52649b5e25fSSatish Balay   PetscFunctionBegin;
527b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
528b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
529fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
530fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
53149b5e25fSSatish Balay   if (issocket) {
53229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
53349b5e25fSSatish Balay   } else if (isascii){
53449b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
53549b5e25fSSatish Balay   } else if (isbinary) {
53649b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
53749b5e25fSSatish Balay   } else if (isdraw) {
53849b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
53949b5e25fSSatish Balay   } else {
540b7aaefc3SHong Zhang     SETERRQ1(1,"Viewer type %s not supported by SeqSBAIJ matrices",((PetscObject)viewer)->type_name);
54149b5e25fSSatish Balay   }
54249b5e25fSSatish Balay   PetscFunctionReturn(0);
54349b5e25fSSatish Balay }
54449b5e25fSSatish Balay 
54549b5e25fSSatish Balay 
5464a2ae208SSatish Balay #undef __FUNCT__
5474a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
548f15d580aSBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],PetscScalar v[])
54949b5e25fSSatish Balay {
550045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
55149b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
55249b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
55349b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
55449b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
55549b5e25fSSatish Balay 
55649b5e25fSSatish Balay   PetscFunctionBegin;
55749b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
55849b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
559590ac198SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
560590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
56149b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
56249b5e25fSSatish Balay     nrow = ailen[brow];
56349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
564590ac198SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %d",in[l]);
565590ac198SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->n-1);
56649b5e25fSSatish Balay       col  = in[l] ;
56749b5e25fSSatish Balay       bcol = col/bs;
56849b5e25fSSatish Balay       cidx = col%bs;
56949b5e25fSSatish Balay       ridx = row%bs;
57049b5e25fSSatish Balay       high = nrow;
57149b5e25fSSatish Balay       low  = 0; /* assume unsorted */
57249b5e25fSSatish Balay       while (high-low > 5) {
57349b5e25fSSatish Balay         t = (low+high)/2;
57449b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
57549b5e25fSSatish Balay         else             low  = t;
57649b5e25fSSatish Balay       }
57749b5e25fSSatish Balay       for (i=low; i<high; i++) {
57849b5e25fSSatish Balay         if (rp[i] > bcol) break;
57949b5e25fSSatish Balay         if (rp[i] == bcol) {
58049b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
58149b5e25fSSatish Balay           goto finished;
58249b5e25fSSatish Balay         }
58349b5e25fSSatish Balay       }
58449b5e25fSSatish Balay       *v++ = zero;
58549b5e25fSSatish Balay       finished:;
58649b5e25fSSatish Balay     }
58749b5e25fSSatish Balay   }
58849b5e25fSSatish Balay   PetscFunctionReturn(0);
58949b5e25fSSatish Balay }
59049b5e25fSSatish Balay 
59149b5e25fSSatish Balay 
5924a2ae208SSatish Balay #undef __FUNCT__
5934a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
594f15d580aSBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
59549b5e25fSSatish Balay {
5960880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
5970880e062SHong Zhang   int             *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
5980880e062SHong Zhang   int             *imax=a->imax,*ai=a->i,*ailen=a->ilen;
5990880e062SHong Zhang   int             *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
6000880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
601f15d580aSBarry Smith   const MatScalar *value = v;
602f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
6030880e062SHong Zhang 
60449b5e25fSSatish Balay   PetscFunctionBegin;
6050880e062SHong Zhang   if (roworiented) {
6060880e062SHong Zhang     stepval = (n-1)*bs;
6070880e062SHong Zhang   } else {
6080880e062SHong Zhang     stepval = (m-1)*bs;
6090880e062SHong Zhang   }
6100880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
6110880e062SHong Zhang     row  = im[k];
6120880e062SHong Zhang     if (row < 0) continue;
6130880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
614590ac198SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->mbs-1);
6150880e062SHong Zhang #endif
6160880e062SHong Zhang     rp   = aj + ai[row];
6170880e062SHong Zhang     ap   = aa + bs2*ai[row];
6180880e062SHong Zhang     rmax = imax[row];
6190880e062SHong Zhang     nrow = ailen[row];
6200880e062SHong Zhang     low  = 0;
6210880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
6220880e062SHong Zhang       if (in[l] < 0) continue;
6230880e062SHong Zhang       col = in[l];
624b1823623SSatish Balay #if defined(PETSC_USE_BOPT_g)
625b1823623SSatish Balay       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",col,a->nbs-1);
626b1823623SSatish Balay #endif
6270880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
6280880e062SHong Zhang       if (roworiented) {
6290880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
6300880e062SHong Zhang       } else {
6310880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
6320880e062SHong Zhang       }
6330880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
6340880e062SHong Zhang       while (high-low > 7) {
6350880e062SHong Zhang         t = (low+high)/2;
6360880e062SHong Zhang         if (rp[t] > col) high = t;
6370880e062SHong Zhang         else             low  = t;
6380880e062SHong Zhang       }
6390880e062SHong Zhang       for (i=low; i<high; i++) {
6400880e062SHong Zhang         if (rp[i] > col) break;
6410880e062SHong Zhang         if (rp[i] == col) {
6420880e062SHong Zhang           bap  = ap +  bs2*i;
6430880e062SHong Zhang           if (roworiented) {
6440880e062SHong Zhang             if (is == ADD_VALUES) {
6450880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6460880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6470880e062SHong Zhang                   bap[jj] += *value++;
6480880e062SHong Zhang                 }
6490880e062SHong Zhang               }
6500880e062SHong Zhang             } else {
6510880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6520880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6530880e062SHong Zhang                   bap[jj] = *value++;
6540880e062SHong Zhang                 }
6550880e062SHong Zhang               }
6560880e062SHong Zhang             }
6570880e062SHong Zhang           } else {
6580880e062SHong Zhang             if (is == ADD_VALUES) {
6590880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6600880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6610880e062SHong Zhang                   *bap++ += *value++;
6620880e062SHong Zhang                 }
6630880e062SHong Zhang               }
6640880e062SHong Zhang             } else {
6650880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6660880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6670880e062SHong Zhang                   *bap++  = *value++;
6680880e062SHong Zhang                 }
6690880e062SHong Zhang               }
6700880e062SHong Zhang             }
6710880e062SHong Zhang           }
6720880e062SHong Zhang           goto noinsert2;
6730880e062SHong Zhang         }
6740880e062SHong Zhang       }
6750880e062SHong Zhang       if (nonew == 1) goto noinsert2;
676*a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6770880e062SHong Zhang       if (nrow >= rmax) {
6780880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
6790880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
6800880e062SHong Zhang         MatScalar *new_a;
6810880e062SHong Zhang 
682*a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
6830880e062SHong Zhang 
6840880e062SHong Zhang         /* malloc new storage space */
6850880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
6860880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
6870880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
6880880e062SHong Zhang         new_i   = new_j + new_nz;
6890880e062SHong Zhang 
6900880e062SHong Zhang         /* copy over old data into new slots */
6910880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
6920880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
6930880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
6940880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
6950880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
6960880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
6970880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
6980880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
6990880e062SHong Zhang         /* free up old matrix storage */
7000880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
7010880e062SHong Zhang         if (!a->singlemalloc) {
7020880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
7030880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
7040880e062SHong Zhang         }
7050880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
7060880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
7070880e062SHong Zhang 
7080880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
7090880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
7100880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
7110880e062SHong Zhang         a->s_maxnz += bs2*CHUNKSIZE;
7120880e062SHong Zhang         a->reallocs++;
7130880e062SHong Zhang         a->s_nz++;
7140880e062SHong Zhang       }
7150880e062SHong Zhang       N = nrow++ - 1;
7160880e062SHong Zhang       /* shift up all the later entries in this row */
7170880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
7180880e062SHong Zhang         rp[ii+1] = rp[ii];
7190880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
7200880e062SHong Zhang       }
7210880e062SHong Zhang       if (N >= i) {
7220880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
7230880e062SHong Zhang       }
7240880e062SHong Zhang       rp[i] = col;
7250880e062SHong Zhang       bap   = ap +  bs2*i;
7260880e062SHong Zhang       if (roworiented) {
7270880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7280880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
7290880e062SHong Zhang             bap[jj] = *value++;
7300880e062SHong Zhang           }
7310880e062SHong Zhang         }
7320880e062SHong Zhang       } else {
7330880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7340880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
7350880e062SHong Zhang             *bap++  = *value++;
7360880e062SHong Zhang           }
7370880e062SHong Zhang         }
7380880e062SHong Zhang       }
7390880e062SHong Zhang       noinsert2:;
7400880e062SHong Zhang       low = i;
7410880e062SHong Zhang     }
7420880e062SHong Zhang     ailen[row] = nrow;
7430880e062SHong Zhang   }
7440880e062SHong Zhang   PetscFunctionReturn(0);
74549b5e25fSSatish Balay }
74649b5e25fSSatish Balay 
7474a2ae208SSatish Balay #undef __FUNCT__
7484a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
74949b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
75049b5e25fSSatish Balay {
75149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
75249b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
753c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
75449b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
75549b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
75649b5e25fSSatish Balay 
75749b5e25fSSatish Balay   PetscFunctionBegin;
75849b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
75949b5e25fSSatish Balay 
76049b5e25fSSatish Balay   if (m) rmax = ailen[0];
76149b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
76249b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
76349b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
76449b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
76549b5e25fSSatish Balay     if (fshift) {
76649b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
76749b5e25fSSatish Balay       N = ailen[i];
76849b5e25fSSatish Balay       for (j=0; j<N; j++) {
76949b5e25fSSatish Balay         ip[j-fshift] = ip[j];
77049b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
77149b5e25fSSatish Balay       }
77249b5e25fSSatish Balay     }
77349b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
77449b5e25fSSatish Balay   }
77549b5e25fSSatish Balay   if (mbs) {
77649b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
77749b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
77849b5e25fSSatish Balay   }
77949b5e25fSSatish Balay   /* reset ilen and imax for each row */
78049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
78149b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
78249b5e25fSSatish Balay   }
78349b5e25fSSatish Balay   a->s_nz = ai[mbs];
78449b5e25fSSatish Balay 
785b424e231SHong Zhang   /* diagonals may have moved, reset it */
786b424e231SHong Zhang   if (a->diag) {
787b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
78849b5e25fSSatish Balay   }
789b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
790c464158bSHong Zhang            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
791b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
79249b5e25fSSatish Balay            a->reallocs);
793b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
79449b5e25fSSatish Balay   a->reallocs          = 0;
79549b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
79649b5e25fSSatish Balay 
79749b5e25fSSatish Balay   PetscFunctionReturn(0);
79849b5e25fSSatish Balay }
79949b5e25fSSatish Balay 
80049b5e25fSSatish Balay /*
80149b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
80249b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
80349b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
80449b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
80549b5e25fSSatish Balay */
8064a2ae208SSatish Balay #undef __FUNCT__
8074a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
808db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
80949b5e25fSSatish Balay {
81049b5e25fSSatish Balay   int        i,j,k,row;
81149b5e25fSSatish Balay   PetscTruth flg;
81249b5e25fSSatish Balay 
81349b5e25fSSatish Balay   PetscFunctionBegin;
81449b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
81549b5e25fSSatish Balay     row = idx[i];
81649b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
81749b5e25fSSatish Balay       sizes[j] = 1;
81849b5e25fSSatish Balay       i++;
81949b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
82049b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
82149b5e25fSSatish Balay       i++;
82249b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
82349b5e25fSSatish Balay       flg = PETSC_TRUE;
82449b5e25fSSatish Balay       for (k=1; k<bs; k++) {
82549b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
82649b5e25fSSatish Balay           flg = PETSC_FALSE;
82749b5e25fSSatish Balay           break;
82849b5e25fSSatish Balay         }
82949b5e25fSSatish Balay       }
83049b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
83149b5e25fSSatish Balay         sizes[j] = bs;
83249b5e25fSSatish Balay         i+= bs;
83349b5e25fSSatish Balay       } else {
83449b5e25fSSatish Balay         sizes[j] = 1;
83549b5e25fSSatish Balay         i++;
83649b5e25fSSatish Balay       }
83749b5e25fSSatish Balay     }
83849b5e25fSSatish Balay   }
83949b5e25fSSatish Balay   *bs_max = j;
84049b5e25fSSatish Balay   PetscFunctionReturn(0);
84149b5e25fSSatish Balay }
84249b5e25fSSatish Balay 
8434a2ae208SSatish Balay #undef __FUNCT__
8444a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
845268466fbSBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,const PetscScalar *diag)
84649b5e25fSSatish Balay {
84749b5e25fSSatish Balay   PetscFunctionBegin;
848c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
84949b5e25fSSatish Balay }
85049b5e25fSSatish Balay 
85149b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
85249b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
85349b5e25fSSatish Balay */
85449b5e25fSSatish Balay 
8554a2ae208SSatish Balay #undef __FUNCT__
8564a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
857f15d580aSBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode is)
85849b5e25fSSatish Balay {
85949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
86049b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
86149b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
86249b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
86349b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
86449b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
86549b5e25fSSatish Balay 
86649b5e25fSSatish Balay   PetscFunctionBegin;
86749b5e25fSSatish Balay 
86849b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
86949b5e25fSSatish Balay     row  = im[k];       /* row number */
87049b5e25fSSatish Balay     brow = row/bs;      /* block row number */
87149b5e25fSSatish Balay     if (row < 0) continue;
87249b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
873590ac198SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m-1);
87449b5e25fSSatish Balay #endif
87549b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
87649b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
87749b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
87849b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
87949b5e25fSSatish Balay     low  = 0;
88049b5e25fSSatish Balay 
88149b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
88249b5e25fSSatish Balay       if (in[l] < 0) continue;
88349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
884590ac198SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m-1);
88549b5e25fSSatish Balay #endif
88649b5e25fSSatish Balay       col = in[l];
88749b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
88849b5e25fSSatish Balay 
88949b5e25fSSatish Balay       if (brow <= bcol){
89049b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8918549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
89249b5e25fSSatish Balay           /* element value a(k,l) */
89349b5e25fSSatish Balay           if (roworiented) {
89449b5e25fSSatish Balay             value = v[l + k*n];
89549b5e25fSSatish Balay           } else {
89649b5e25fSSatish Balay             value = v[k + l*m];
89749b5e25fSSatish Balay           }
89849b5e25fSSatish Balay 
89949b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
90049b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
90149b5e25fSSatish Balay           while (high-low > 7) {
90249b5e25fSSatish Balay             t = (low+high)/2;
90349b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
90449b5e25fSSatish Balay             else              low  = t;
90549b5e25fSSatish Balay           }
90649b5e25fSSatish Balay           for (i=low; i<high; i++) {
90749b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
90849b5e25fSSatish Balay             if (rp[i] > bcol) break;
90949b5e25fSSatish Balay             if (rp[i] == bcol) {
91049b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
91149b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
91249b5e25fSSatish Balay               else                  *bap  = value;
9138549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9148549e402SHong Zhang               if (brow == bcol && ridx < cidx){
9158549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
9168549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
9178549e402SHong Zhang                 else                  *bap  = value;
9188549e402SHong Zhang               }
91949b5e25fSSatish Balay               goto noinsert1;
92049b5e25fSSatish Balay             }
92149b5e25fSSatish Balay           }
92249b5e25fSSatish Balay 
92349b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
924*a45adfd6SMatthew Knepley       else if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
92549b5e25fSSatish Balay       if (nrow >= rmax) {
92649b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
92749b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
92849b5e25fSSatish Balay         MatScalar *new_a;
92949b5e25fSSatish Balay 
930*a45adfd6SMatthew Knepley         if (nonew == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%d, %d) in the matrix", row, col);
93149b5e25fSSatish Balay 
93249b5e25fSSatish Balay         /* Malloc new storage space */
93349b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
93482502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
93549b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
93649b5e25fSSatish Balay         new_i = new_j + new_nz;
93749b5e25fSSatish Balay 
93849b5e25fSSatish Balay         /* copy over old data into new slots */
93949b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
94049b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
94149b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
94249b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
94349b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
94449b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
94549b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
94649b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
94749b5e25fSSatish Balay         /* free up old matrix storage */
94849b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
94949b5e25fSSatish Balay         if (!a->singlemalloc) {
95049b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
95149b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
95249b5e25fSSatish Balay         }
95349b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
95449b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
95549b5e25fSSatish Balay 
95649b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
95749b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
958b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
95949b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
96049b5e25fSSatish Balay         a->reallocs++;
96149b5e25fSSatish Balay         a->s_nz++;
96249b5e25fSSatish Balay       }
96349b5e25fSSatish Balay 
96449b5e25fSSatish Balay       N = nrow++ - 1;
96549b5e25fSSatish Balay       /* shift up all the later entries in this row */
96649b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
96749b5e25fSSatish Balay         rp[ii+1] = rp[ii];
96849b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
96949b5e25fSSatish Balay       }
97049b5e25fSSatish Balay       if (N>=i) {
97149b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
97249b5e25fSSatish Balay       }
97349b5e25fSSatish Balay       rp[i]                      = bcol;
97449b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
97549b5e25fSSatish Balay       noinsert1:;
97649b5e25fSSatish Balay       low = i;
97749b5e25fSSatish Balay       /* } */
9788549e402SHong Zhang         }
97949b5e25fSSatish Balay       } /* end of if .. if..  */
98049b5e25fSSatish Balay     }                     /* end of loop over added columns */
98149b5e25fSSatish Balay     ailen[brow] = nrow;
98249b5e25fSSatish Balay   }                       /* end of loop over added rows */
98349b5e25fSSatish Balay 
98449b5e25fSSatish Balay   PetscFunctionReturn(0);
98549b5e25fSSatish Balay }
98649b5e25fSSatish Balay 
98715e8a5b3SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
98815e8a5b3SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
98952ebccd6SSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS[],int);
99049b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
991268466fbSBarry Smith extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,const IS[],const IS[],MatReuse,Mat*[]);
99249b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
99349b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
994268466fbSBarry Smith extern int MatScale_SeqSBAIJ(const PetscScalar*,Mat);
99549b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
99649b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
99749b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
99849b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
99949b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
100049b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
100113a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
100249b5e25fSSatish Balay 
100349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
100449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
100549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
100649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
100749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
100849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
100949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
101049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
101149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
101249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
101349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
101449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
101549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
101649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
101749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
101849b5e25fSSatish Balay 
1019d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1020d59c15a7SBarry Smith 
102107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
102207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
102307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
102407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
102507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
102607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
102707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
102807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
102907f98182SSatish Balay 
10305f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
10315f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
10325f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
10335f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
10345f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
10355f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
10365f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
103757d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
10383e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
103949b5e25fSSatish Balay 
104049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
104149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
104249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
104349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
104449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
104549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
104649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
104749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
104849b5e25fSSatish Balay 
104949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
105049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
105149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
105249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
105349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
105449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
105549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
105649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
105749b5e25fSSatish Balay 
10584d101231SSatish Balay #ifdef HAVE_MatICCFactor
10591a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
10604a2ae208SSatish Balay #undef __FUNCT__
10614d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
106215e8a5b3SHong Zhang int MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
106349b5e25fSSatish Balay {
10644ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
106549b5e25fSSatish Balay   Mat         outA;
106649b5e25fSSatish Balay   int         ierr;
106749b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
106849b5e25fSSatish Balay 
106949b5e25fSSatish Balay   PetscFunctionBegin;
10701a3463dfSHong Zhang   /*
107129bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
107249b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
107349b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
107449b5e25fSSatish Balay   if (!row_identity || !col_identity) {
107529bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
107649b5e25fSSatish Balay   }
10771a3463dfSHong Zhang   */
107849b5e25fSSatish Balay 
107949b5e25fSSatish Balay   outA          = inA;
10801260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
108149b5e25fSSatish Balay 
108249b5e25fSSatish Balay   if (!a->diag) {
10831a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
108449b5e25fSSatish Balay   }
108549b5e25fSSatish Balay   /*
108649b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
108749b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
108849b5e25fSSatish Balay   */
108949b5e25fSSatish Balay   switch (a->bs) {
109049b5e25fSSatish Balay   case 1:
10910fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
109207f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1093d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
10944d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
109549b5e25fSSatish Balay   case 2:
10961a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
10971a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10981a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10994d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
110049b5e25fSSatish Balay     break;
110149b5e25fSSatish Balay   case 3:
11021a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
11031a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11041a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11054d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
110649b5e25fSSatish Balay     break;
110749b5e25fSSatish Balay   case 4:
11081a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
11091a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11101a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11114d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
111249b5e25fSSatish Balay     break;
111349b5e25fSSatish Balay   case 5:
11141a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
11151a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11161a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11174d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
111849b5e25fSSatish Balay     break;
111949b5e25fSSatish Balay   case 6:
11201a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
11211a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11221a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11234d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
112449b5e25fSSatish Balay     break;
112549b5e25fSSatish Balay   case 7:
11261a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
11271a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11281a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11294d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
113049b5e25fSSatish Balay     break;
113149b5e25fSSatish Balay   default:
113249b5e25fSSatish Balay     a->row        = row;
11331a3463dfSHong Zhang     a->icol       = col;
113449b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
113549b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
113649b5e25fSSatish Balay 
113749b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
113849b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1139b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
114049b5e25fSSatish Balay 
114149b5e25fSSatish Balay     if (!a->solve_work) {
114287828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
114387828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
114449b5e25fSSatish Balay     }
114549b5e25fSSatish Balay   }
114649b5e25fSSatish Balay 
11471a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
114849b5e25fSSatish Balay 
114949b5e25fSSatish Balay   PetscFunctionReturn(0);
115049b5e25fSSatish Balay }
1151950f1e5bSHong Zhang #endif
1152950f1e5bSHong Zhang 
11534a2ae208SSatish Balay #undef __FUNCT__
11544a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
115549b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
115649b5e25fSSatish Balay {
115749b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
115849b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
115949b5e25fSSatish Balay   int               ierr;
116049b5e25fSSatish Balay 
116149b5e25fSSatish Balay   PetscFunctionBegin;
116249b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
116349b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
116449b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
116549b5e25fSSatish Balay   PetscFunctionReturn(0);
116649b5e25fSSatish Balay }
116749b5e25fSSatish Balay 
116849b5e25fSSatish Balay EXTERN_C_BEGIN
11694a2ae208SSatish Balay #undef __FUNCT__
11704a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
117149b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
117249b5e25fSSatish Balay {
1173045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
117449b5e25fSSatish Balay   int         i,nz,n;
117549b5e25fSSatish Balay 
117649b5e25fSSatish Balay   PetscFunctionBegin;
1177045c9aa0SHong Zhang   nz = baij->s_maxnz;
1178c464158bSHong Zhang   n  = mat->n;
117949b5e25fSSatish Balay   for (i=0; i<nz; i++) {
118049b5e25fSSatish Balay     baij->j[i] = indices[i];
118149b5e25fSSatish Balay   }
1182045c9aa0SHong Zhang   baij->s_nz = nz;
118349b5e25fSSatish Balay   for (i=0; i<n; i++) {
118449b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
118549b5e25fSSatish Balay   }
118649b5e25fSSatish Balay 
118749b5e25fSSatish Balay   PetscFunctionReturn(0);
118849b5e25fSSatish Balay }
118949b5e25fSSatish Balay EXTERN_C_END
119049b5e25fSSatish Balay 
11914a2ae208SSatish Balay #undef __FUNCT__
11924a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
119349b5e25fSSatish Balay /*@
119419585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
119549b5e25fSSatish Balay        in the matrix.
119649b5e25fSSatish Balay 
119749b5e25fSSatish Balay   Input Parameters:
119819585528SSatish Balay +  mat     - the SeqSBAIJ matrix
119949b5e25fSSatish Balay -  indices - the column indices
120049b5e25fSSatish Balay 
120149b5e25fSSatish Balay   Level: advanced
120249b5e25fSSatish Balay 
120349b5e25fSSatish Balay   Notes:
120449b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
120549b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
120649b5e25fSSatish Balay   of the MatSetValues() operation.
120749b5e25fSSatish Balay 
120849b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
120919585528SSatish Balay   MatCreateSeqSBAIJ().
121049b5e25fSSatish Balay 
1211ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
121249b5e25fSSatish Balay 
1213ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
121449b5e25fSSatish Balay @*/
121549b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
121649b5e25fSSatish Balay {
121749b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
121849b5e25fSSatish Balay 
121949b5e25fSSatish Balay   PetscFunctionBegin;
122049b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1221c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
122249b5e25fSSatish Balay   if (f) {
122349b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
122449b5e25fSSatish Balay   } else {
122529bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
122649b5e25fSSatish Balay   }
122749b5e25fSSatish Balay   PetscFunctionReturn(0);
122849b5e25fSSatish Balay }
122949b5e25fSSatish Balay 
12304a2ae208SSatish Balay #undef __FUNCT__
12314a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1232273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1233273d9f13SBarry Smith {
1234273d9f13SBarry Smith   int        ierr;
1235273d9f13SBarry Smith 
1236273d9f13SBarry Smith   PetscFunctionBegin;
1237273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1238273d9f13SBarry Smith   PetscFunctionReturn(0);
1239273d9f13SBarry Smith }
1240273d9f13SBarry Smith 
1241a6ece127SHong Zhang #undef __FUNCT__
1242a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
12434e7234bfSBarry Smith int MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1244a6ece127SHong Zhang {
1245a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1246a6ece127SHong Zhang   PetscFunctionBegin;
1247a6ece127SHong Zhang   *array = a->a;
1248a6ece127SHong Zhang   PetscFunctionReturn(0);
1249a6ece127SHong Zhang }
1250a6ece127SHong Zhang 
1251a6ece127SHong Zhang #undef __FUNCT__
1252a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
12534e7234bfSBarry Smith int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1254a6ece127SHong Zhang {
1255a6ece127SHong Zhang   PetscFunctionBegin;
1256a6ece127SHong Zhang   PetscFunctionReturn(0);
1257a6ece127SHong Zhang }
1258a6ece127SHong Zhang 
125942ee4b1aSHong Zhang #include "petscblaslapack.h"
126042ee4b1aSHong Zhang #undef __FUNCT__
126142ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1262268466fbSBarry Smith int MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
126342ee4b1aSHong Zhang {
126442ee4b1aSHong Zhang   Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1265c4319e64SHong Zhang   int          ierr,one=1,i,bs=y->bs,bs2,j;
126642ee4b1aSHong Zhang 
126742ee4b1aSHong Zhang   PetscFunctionBegin;
126842ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1269268466fbSBarry Smith     BLaxpy_(&x->s_nz,(PetscScalar*)a,x->a,&one,y->a,&one);
1270c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1271c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1272c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1273c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1274c537a176SHong Zhang     }
1275c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1276c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1277c4319e64SHong Zhang       y->XtoY = X;
1278c537a176SHong Zhang     }
1279c4319e64SHong Zhang     bs2 = bs*bs;
1280c537a176SHong Zhang     for (i=0; i<x->s_nz; i++) {
1281c4319e64SHong Zhang       j = 0;
1282c4319e64SHong Zhang       while (j < bs2){
12836550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1284c4319e64SHong Zhang         j++;
1285c537a176SHong Zhang       }
1286c4319e64SHong Zhang     }
1287c4319e64SHong Zhang     PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %d/%d = %g\n",bs2*x->s_nz,bs2*y->s_nz,(PetscReal)(bs2*x->s_nz)/(bs2*y->s_nz));
128842ee4b1aSHong Zhang   } else {
128942ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
129042ee4b1aSHong Zhang   }
129142ee4b1aSHong Zhang   PetscFunctionReturn(0);
129242ee4b1aSHong Zhang }
129342ee4b1aSHong Zhang 
129449b5e25fSSatish Balay /* -------------------------------------------------------------------*/
129549b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
129649b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
129749b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
129849b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
129997304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
130049b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
130149b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
130249b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
130349b5e25fSSatish Balay        0,
130449b5e25fSSatish Balay        0,
130597304618SKris Buschelman /*10*/ 0,
130649b5e25fSSatish Balay        0,
13075f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1308d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
130949b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
131097304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
131149b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
131249b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
131349b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
131449b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
131597304618SKris Buschelman /*20*/ 0,
131649b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
131749b5e25fSSatish Balay        0,
131849b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
131949b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
132097304618SKris Buschelman /*25*/ MatZeroRows_SeqSBAIJ,
132149b5e25fSSatish Balay        0,
132249b5e25fSSatish Balay        0,
13235f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13245f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
132597304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1326c464158bSHong Zhang        0,
13274d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1328a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1329a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
133097304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
133149b5e25fSSatish Balay        0,
133249b5e25fSSatish Balay        0,
133349b5e25fSSatish Balay        0,
1334950f1e5bSHong Zhang        0,
133597304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
133649b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
133749b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
133849b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
133949b5e25fSSatish Balay        0,
134097304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
134149b5e25fSSatish Balay        MatScale_SeqSBAIJ,
134249b5e25fSSatish Balay        0,
134349b5e25fSSatish Balay        0,
134449b5e25fSSatish Balay        0,
134597304618SKris Buschelman /*50*/ MatGetBlockSize_SeqSBAIJ,
134649b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
134749b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
134849b5e25fSSatish Balay        0,
134949b5e25fSSatish Balay        0,
135097304618SKris Buschelman /*55*/ 0,
135149b5e25fSSatish Balay        0,
135249b5e25fSSatish Balay        0,
135349b5e25fSSatish Balay        0,
135449b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
135597304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
135649b5e25fSSatish Balay        0,
135749b5e25fSSatish Balay        0,
13588a124369SBarry Smith        MatGetPetscMaps_Petsc,
1359d959ec07SHong Zhang        0,
136097304618SKris Buschelman /*65*/ 0,
1361d959ec07SHong Zhang        0,
1362d959ec07SHong Zhang        0,
1363d959ec07SHong Zhang        0,
1364d959ec07SHong Zhang        0,
136597304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
13663e0d88b5SBarry Smith        0,
13673e0d88b5SBarry Smith        0,
13683e0d88b5SBarry Smith        0,
13693e0d88b5SBarry Smith        0,
137097304618SKris Buschelman /*75*/ 0,
13713e0d88b5SBarry Smith        0,
13723e0d88b5SBarry Smith        0,
13733e0d88b5SBarry Smith        0,
13743e0d88b5SBarry Smith        0,
137597304618SKris Buschelman /*80*/ 0,
13763e0d88b5SBarry Smith        0,
13773e0d88b5SBarry Smith        0,
13783e0d88b5SBarry Smith        0,
13793e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
138097304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
13813e0d88b5SBarry Smith #else
138297304618SKris Buschelman        0,
13833e0d88b5SBarry Smith #endif
138497304618SKris Buschelman /*85*/ MatLoad_SeqSBAIJ
13853e0d88b5SBarry Smith };
138649b5e25fSSatish Balay 
138749b5e25fSSatish Balay EXTERN_C_BEGIN
13884a2ae208SSatish Balay #undef __FUNCT__
13894a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
139049b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
139149b5e25fSSatish Balay {
13924afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1393c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
139449b5e25fSSatish Balay   int         ierr;
139549b5e25fSSatish Balay 
139649b5e25fSSatish Balay   PetscFunctionBegin;
139749b5e25fSSatish Balay   if (aij->nonew != 1) {
139829bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
139949b5e25fSSatish Balay   }
140049b5e25fSSatish Balay 
140149b5e25fSSatish Balay   /* allocate space for values if not already there */
140249b5e25fSSatish Balay   if (!aij->saved_values) {
140387828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
140449b5e25fSSatish Balay   }
140549b5e25fSSatish Balay 
140649b5e25fSSatish Balay   /* copy values over */
140787828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
140849b5e25fSSatish Balay   PetscFunctionReturn(0);
140949b5e25fSSatish Balay }
141049b5e25fSSatish Balay EXTERN_C_END
141149b5e25fSSatish Balay 
141249b5e25fSSatish Balay EXTERN_C_BEGIN
14134a2ae208SSatish Balay #undef __FUNCT__
14144a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
141549b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
141649b5e25fSSatish Balay {
14174afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1418c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
141949b5e25fSSatish Balay 
142049b5e25fSSatish Balay   PetscFunctionBegin;
142149b5e25fSSatish Balay   if (aij->nonew != 1) {
142229bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
142349b5e25fSSatish Balay   }
142449b5e25fSSatish Balay   if (!aij->saved_values) {
142529bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
142649b5e25fSSatish Balay   }
142749b5e25fSSatish Balay 
142849b5e25fSSatish Balay   /* copy values over */
142987828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
143049b5e25fSSatish Balay   PetscFunctionReturn(0);
143149b5e25fSSatish Balay }
143249b5e25fSSatish Balay EXTERN_C_END
143349b5e25fSSatish Balay 
14348549e402SHong Zhang EXTERN_C_BEGIN
14354a2ae208SSatish Balay #undef __FUNCT__
1436a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1437a23d5eceSKris Buschelman int MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz)
143849b5e25fSSatish Balay {
1439c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
14400c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
144149b5e25fSSatish Balay   PetscTruth   flg;
14424afc71dfSHong Zhang   int          s_nz;
144349b5e25fSSatish Balay 
144449b5e25fSSatish Balay   PetscFunctionBegin;
1445273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1446e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1447c464158bSHong Zhang   mbs  = B->m/bs;
144849b5e25fSSatish Balay   bs2  = bs*bs;
144949b5e25fSSatish Balay 
1450c464158bSHong Zhang   if (mbs*bs != B->m) {
145129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
145249b5e25fSSatish Balay   }
145349b5e25fSSatish Balay 
1454435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1455435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
145649b5e25fSSatish Balay   if (nnz) {
145749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
145829bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
145980fe4e49SBarry 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);
146049b5e25fSSatish Balay     }
146149b5e25fSSatish Balay   }
146249b5e25fSSatish Balay 
1463e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
146449b5e25fSSatish Balay   if (!flg) {
146549b5e25fSSatish Balay     switch (bs) {
146649b5e25fSSatish Balay     case 1:
14674ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
146849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1469d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
147049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
147149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
147249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
147349b5e25fSSatish Balay       break;
147449b5e25fSSatish Balay     case 2:
14754ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
147649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
147749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
147849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
147949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
148049b5e25fSSatish Balay       break;
148149b5e25fSSatish Balay     case 3:
14825f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
148349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
148449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
148549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
148649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
148749b5e25fSSatish Balay       break;
148849b5e25fSSatish Balay     case 4:
14895f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
149049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
149149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
149249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
149349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
149449b5e25fSSatish Balay       break;
149549b5e25fSSatish Balay     case 5:
14965f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
149749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
149849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
149949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
150049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
150149b5e25fSSatish Balay       break;
150249b5e25fSSatish Balay     case 6:
15035f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
150449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
150549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
150649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
150749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
150849b5e25fSSatish Balay       break;
150949b5e25fSSatish Balay     case 7:
1510de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
151149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
151249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1513de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
151449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
151549b5e25fSSatish Balay       break;
151649b5e25fSSatish Balay     }
151749b5e25fSSatish Balay   }
151849b5e25fSSatish Balay 
151949b5e25fSSatish Balay   b->mbs = mbs;
15204afc71dfSHong Zhang   b->nbs = mbs;
1521b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
152249b5e25fSSatish Balay   if (!nnz) {
1523435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
152449b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
152549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
15268cef66ccSHong Zhang       b->imax[i] = nz;
152749b5e25fSSatish Balay     }
1528153ea458SHong Zhang     nz = nz*mbs; /* total nz */
152949b5e25fSSatish Balay   } else {
153049b5e25fSSatish Balay     nz = 0;
15318cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
153249b5e25fSSatish Balay   }
15338cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
15348cef66ccSHong Zhang   s_nz = nz;
153549b5e25fSSatish Balay 
153649b5e25fSSatish Balay   /* allocate the matrix space */
1537c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
153882502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
153949b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
154049b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
154149b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
154249b5e25fSSatish Balay   b->i = b->j + s_nz;
154349b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
154449b5e25fSSatish Balay 
154549b5e25fSSatish Balay   /* pointer to beginning of each row */
154649b5e25fSSatish Balay   b->i[0] = 0;
154749b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
154849b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
154949b5e25fSSatish Balay   }
155049b5e25fSSatish Balay 
155149b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
15525033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1553b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
155449b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
155549b5e25fSSatish Balay 
155649b5e25fSSatish Balay   b->bs               = bs;
155749b5e25fSSatish Balay   b->bs2              = bs2;
155849b5e25fSSatish Balay   b->s_nz             = 0;
155949b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1560153ea458SHong Zhang 
156116cdd363SHong Zhang   b->inew             = 0;
156216cdd363SHong Zhang   b->jnew             = 0;
156316cdd363SHong Zhang   b->anew             = 0;
156416cdd363SHong Zhang   b->a2anew           = 0;
15651a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1566c464158bSHong Zhang   PetscFunctionReturn(0);
1567c464158bSHong Zhang }
1568a23d5eceSKris Buschelman EXTERN_C_END
1569153ea458SHong Zhang 
15700bad9183SKris Buschelman /*MC
1571fafad747SKris Buschelman    MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
15720bad9183SKris Buschelman    based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
15730bad9183SKris Buschelman 
15740bad9183SKris Buschelman    Options Database Keys:
15750bad9183SKris Buschelman . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
15760bad9183SKris Buschelman 
15770bad9183SKris Buschelman   Level: beginner
15780bad9183SKris Buschelman 
15790bad9183SKris Buschelman .seealso: MatCreateSeqSBAIJ
15800bad9183SKris Buschelman M*/
15810bad9183SKris Buschelman 
1582a23d5eceSKris Buschelman EXTERN_C_BEGIN
1583a23d5eceSKris Buschelman #undef __FUNCT__
1584a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1585a23d5eceSKris Buschelman int MatCreate_SeqSBAIJ(Mat B)
1586a23d5eceSKris Buschelman {
1587a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
1588a23d5eceSKris Buschelman   int          ierr,size;
1589a23d5eceSKris Buschelman 
1590a23d5eceSKris Buschelman   PetscFunctionBegin;
1591a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1592a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1593a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1594a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1595a23d5eceSKris Buschelman 
1596a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1597a23d5eceSKris Buschelman   B->data = (void*)b;
1598a23d5eceSKris Buschelman   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1599a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1600a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1601a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1602a23d5eceSKris Buschelman   B->factor           = 0;
1603a23d5eceSKris Buschelman   B->lupivotthreshold = 1.0;
1604a23d5eceSKris Buschelman   B->mapping          = 0;
1605a23d5eceSKris Buschelman   b->row              = 0;
1606a23d5eceSKris Buschelman   b->icol             = 0;
1607a23d5eceSKris Buschelman   b->reallocs         = 0;
1608a23d5eceSKris Buschelman   b->saved_values     = 0;
1609a23d5eceSKris Buschelman 
1610a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1611a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1612a23d5eceSKris Buschelman 
1613a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1614a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1615a23d5eceSKris Buschelman   b->nonew            = 0;
1616a23d5eceSKris Buschelman   b->diag             = 0;
1617a23d5eceSKris Buschelman   b->solve_work       = 0;
1618a23d5eceSKris Buschelman   b->mult_work        = 0;
1619a23d5eceSKris Buschelman   B->spptr            = 0;
1620a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1621a23d5eceSKris Buschelman   b->xtoy             = 0;
1622a23d5eceSKris Buschelman   b->XtoY             = 0;
1623a23d5eceSKris Buschelman 
1624a23d5eceSKris Buschelman   b->inew             = 0;
1625a23d5eceSKris Buschelman   b->jnew             = 0;
1626a23d5eceSKris Buschelman   b->anew             = 0;
1627a23d5eceSKris Buschelman   b->a2anew           = 0;
1628a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1629a23d5eceSKris Buschelman 
1630a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1631a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1632a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1633a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1634a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1635a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1636a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1637a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1638a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1639a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1640a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1641a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
1642a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1643a23d5eceSKris Buschelman }
1644a23d5eceSKris Buschelman EXTERN_C_END
1645a23d5eceSKris Buschelman 
1646a23d5eceSKris Buschelman #undef __FUNCT__
1647a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1648a23d5eceSKris Buschelman /*@C
1649a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1650a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1651a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1652a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1653a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1654a23d5eceSKris Buschelman 
1655a23d5eceSKris Buschelman    Collective on Mat
1656a23d5eceSKris Buschelman 
1657a23d5eceSKris Buschelman    Input Parameters:
1658a23d5eceSKris Buschelman +  A - the symmetric matrix
1659a23d5eceSKris Buschelman .  bs - size of block
1660a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1661a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1662a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1663a23d5eceSKris Buschelman 
1664a23d5eceSKris Buschelman    Options Database Keys:
1665a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1666a23d5eceSKris Buschelman                      block calculations (much slower)
1667a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1668a23d5eceSKris Buschelman 
1669a23d5eceSKris Buschelman    Level: intermediate
1670a23d5eceSKris Buschelman 
1671a23d5eceSKris Buschelman    Notes:
1672a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1673a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1674a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1675a23d5eceSKris Buschelman    matrices.
1676a23d5eceSKris Buschelman 
1677a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1678a23d5eceSKris Buschelman @*/
1679ca01db9bSBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,const int nnz[]) {
1680ca01db9bSBarry Smith   int ierr,(*f)(Mat,int,int,const int[]);
1681a23d5eceSKris Buschelman 
1682a23d5eceSKris Buschelman   PetscFunctionBegin;
1683a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1684a23d5eceSKris Buschelman   if (f) {
1685a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1686a23d5eceSKris Buschelman   }
1687a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1688a23d5eceSKris Buschelman }
168949b5e25fSSatish Balay 
16904a2ae208SSatish Balay #undef __FUNCT__
16914a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1692c464158bSHong Zhang /*@C
1693c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1694c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1695c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1696c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1697c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
169849b5e25fSSatish Balay 
1699c464158bSHong Zhang    Collective on MPI_Comm
1700c464158bSHong Zhang 
1701c464158bSHong Zhang    Input Parameters:
1702c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1703c464158bSHong Zhang .  bs - size of block
1704c464158bSHong Zhang .  m - number of rows, or number of columns
1705c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1706744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1707744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1708c464158bSHong Zhang 
1709c464158bSHong Zhang    Output Parameter:
1710c464158bSHong Zhang .  A - the symmetric matrix
1711c464158bSHong Zhang 
1712c464158bSHong Zhang    Options Database Keys:
1713c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1714c464158bSHong Zhang                      block calculations (much slower)
1715c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1716c464158bSHong Zhang 
1717c464158bSHong Zhang    Level: intermediate
1718c464158bSHong Zhang 
1719c464158bSHong Zhang    Notes:
1720c464158bSHong Zhang 
1721c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1722c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1723c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1724c464158bSHong Zhang    matrices.
1725c464158bSHong Zhang 
1726c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1727c464158bSHong Zhang @*/
1728ca01db9bSBarry Smith int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,const int nnz[],Mat *A)
1729c464158bSHong Zhang {
1730c464158bSHong Zhang   int ierr;
1731c464158bSHong Zhang 
1732c464158bSHong Zhang   PetscFunctionBegin;
1733c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1734c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1735273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
173649b5e25fSSatish Balay   PetscFunctionReturn(0);
173749b5e25fSSatish Balay }
173849b5e25fSSatish Balay 
17394a2ae208SSatish Balay #undef __FUNCT__
17404a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
174149b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
174249b5e25fSSatish Balay {
174349b5e25fSSatish Balay   Mat         C;
174449b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
174549b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
174649b5e25fSSatish Balay 
174749b5e25fSSatish Balay   PetscFunctionBegin;
174829bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
174949b5e25fSSatish Balay 
175049b5e25fSSatish Balay   *B = 0;
1751692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1752692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1753692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1754692f9cbeSHong Zhang 
175549b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1756273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
175749b5e25fSSatish Balay   C->factor         = A->factor;
175849b5e25fSSatish Balay   c->row            = 0;
175949b5e25fSSatish Balay   c->icol           = 0;
176049b5e25fSSatish Balay   c->saved_values   = 0;
176149b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
176249b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
176349b5e25fSSatish Balay 
176449b5e25fSSatish Balay   c->bs         = a->bs;
176549b5e25fSSatish Balay   c->bs2        = a->bs2;
176649b5e25fSSatish Balay   c->mbs        = a->mbs;
176749b5e25fSSatish Balay   c->nbs        = a->nbs;
176849b5e25fSSatish Balay 
1769b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1770b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
177149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
177249b5e25fSSatish Balay     c->imax[i] = a->imax[i];
177349b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
177449b5e25fSSatish Balay   }
177549b5e25fSSatish Balay 
177649b5e25fSSatish Balay   /* allocate the matrix space */
177749b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
177849b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
177982502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
178049b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
178149b5e25fSSatish Balay   c->i = c->j + nz;
178249b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
178349b5e25fSSatish Balay   if (mbs > 0) {
178449b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
178549b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
178649b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178749b5e25fSSatish Balay     } else {
178849b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178949b5e25fSSatish Balay     }
179049b5e25fSSatish Balay   }
179149b5e25fSSatish Balay 
1792b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
179349b5e25fSSatish Balay   c->sorted      = a->sorted;
179449b5e25fSSatish Balay   c->roworiented = a->roworiented;
179549b5e25fSSatish Balay   c->nonew       = a->nonew;
179649b5e25fSSatish Balay 
179749b5e25fSSatish Balay   if (a->diag) {
1798b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1799b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
180049b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
180149b5e25fSSatish Balay       c->diag[i] = a->diag[i];
180249b5e25fSSatish Balay     }
180349b5e25fSSatish Balay   } else c->diag        = 0;
180449b5e25fSSatish Balay   c->s_nz               = a->s_nz;
180549b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
180649b5e25fSSatish Balay   c->solve_work         = 0;
18072a1b7f2aSHong Zhang   C->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
180849b5e25fSSatish Balay   c->mult_work          = 0;
180949b5e25fSSatish Balay   *B = C;
1810b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
181149b5e25fSSatish Balay   PetscFunctionReturn(0);
181249b5e25fSSatish Balay }
181349b5e25fSSatish Balay 
18144a2ae208SSatish Balay #undef __FUNCT__
18154a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1816b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
181749b5e25fSSatish Balay {
181849b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
181949b5e25fSSatish Balay   Mat          B;
182049b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18213f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
182249b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
182349b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
182487828ca2SBarry Smith   PetscScalar  *aa;
182549b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
182649b5e25fSSatish Balay 
182749b5e25fSSatish Balay   PetscFunctionBegin;
1828b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
182949b5e25fSSatish Balay   bs2  = bs*bs;
183049b5e25fSSatish Balay 
183149b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
183229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1833b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
183449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1835552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
183649b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
183749b5e25fSSatish Balay 
183849b5e25fSSatish Balay   if (header[3] < 0) {
183929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
184049b5e25fSSatish Balay   }
184149b5e25fSSatish Balay 
184229bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
184349b5e25fSSatish Balay 
184449b5e25fSSatish Balay   /*
184549b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
184649b5e25fSSatish Balay     divisible by the blocksize
184749b5e25fSSatish Balay   */
184849b5e25fSSatish Balay   mbs        = M/bs;
184949b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
185049b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
185149b5e25fSSatish Balay   else                  mbs++;
185249b5e25fSSatish Balay   if (extra_rows) {
1853b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
185449b5e25fSSatish Balay   }
185549b5e25fSSatish Balay 
185649b5e25fSSatish Balay   /* read in row lengths */
1857b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
185849b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
185949b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
186049b5e25fSSatish Balay 
186149b5e25fSSatish Balay   /* read in column indices */
1862b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
186349b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
186449b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
186549b5e25fSSatish Balay 
186649b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
186782502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1868a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
186982502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
187049b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
187149b5e25fSSatish Balay   masked   = mask + mbs;
187249b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
187349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
187449b5e25fSSatish Balay     nmask = 0;
187549b5e25fSSatish Balay     for (j=0; j<bs; j++) {
187649b5e25fSSatish Balay       kmax = rowlengths[rowcount];
187749b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18782d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
187903630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
188049b5e25fSSatish Balay       }
188149b5e25fSSatish Balay       rowcount++;
188249b5e25fSSatish Balay     }
1883574b2666SHong Zhang     s_browlengths[i] += nmask;
1884574b2666SHong Zhang 
188549b5e25fSSatish Balay     /* zero out the mask elements we set */
188649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
188749b5e25fSSatish Balay   }
188849b5e25fSSatish Balay 
188949b5e25fSSatish Balay   /* create our matrix */
18909abb65ffSKris Buschelman   ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr);
18919abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
189278473fd7SKris Buschelman   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,s_browlengths);CHKERRQ(ierr);
189349b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
189449b5e25fSSatish Balay 
189549b5e25fSSatish Balay   /* set matrix "i" values */
189649b5e25fSSatish Balay   a->i[0] = 0;
189749b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1898574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1899574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
190049b5e25fSSatish Balay   }
19017fe2be48SHong Zhang   a->s_nz = a->i[mbs];
190249b5e25fSSatish Balay 
190349b5e25fSSatish Balay   /* read in nonzero values */
190487828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
190549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
190649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
190749b5e25fSSatish Balay 
190849b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
190949b5e25fSSatish Balay   nzcount = 0; jcount = 0;
191049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
191149b5e25fSSatish Balay     nzcountb = nzcount;
191249b5e25fSSatish Balay     nmask    = 0;
191349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
191449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
191549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19162d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
191703630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19182d703238SHong Zhang       }
19192d703238SHong Zhang     }
19202d703238SHong Zhang     /* sort the masked values */
19212d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19222d703238SHong Zhang 
19232d703238SHong Zhang     /* set "j" values into matrix */
19242d703238SHong Zhang     maskcount = 1;
19252d703238SHong Zhang     for (j=0; j<nmask; j++) {
192649b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
192749b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
192849b5e25fSSatish Balay     }
1929574b2666SHong Zhang 
193049b5e25fSSatish Balay     /* set "a" values into matrix */
193149b5e25fSSatish Balay     ishift = bs2*a->i[i];
193249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
193349b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
193449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1935574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1936574b2666SHong Zhang         if (tmp >= i){
193749b5e25fSSatish Balay           block     = mask[tmp] - 1;
193849b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
193949b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1940574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1941574b2666SHong Zhang         }
1942574b2666SHong Zhang         nzcountb++;
194349b5e25fSSatish Balay       }
194449b5e25fSSatish Balay     }
194549b5e25fSSatish Balay     /* zero out the mask elements we set */
194649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
194749b5e25fSSatish Balay   }
194829bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
194949b5e25fSSatish Balay 
195049b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1951574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
195249b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
195349b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
195449b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
195549b5e25fSSatish Balay 
19569abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19579abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
195849b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19599abb65ffSKris Buschelman   *A = B;
196049b5e25fSSatish Balay   PetscFunctionReturn(0);
196149b5e25fSSatish Balay }
1962574b2666SHong Zhang 
1963d06b337dSHong Zhang #undef __FUNCT__
1964d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1965c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1966d06b337dSHong Zhang {
1967d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1968d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1969d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
1970d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
1971d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1972d06b337dSHong Zhang 
1973d06b337dSHong Zhang   PetscFunctionBegin;
1974b965ef7fSBarry Smith   its = its*lits;
197591723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1976d06b337dSHong Zhang 
1977d06b337dSHong Zhang   if (bs > 1)
1978d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1979d06b337dSHong Zhang 
1980d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1981d06b337dSHong Zhang   if (xx != bb) {
1982d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1983d06b337dSHong Zhang   } else {
1984d06b337dSHong Zhang     b = x;
1985d06b337dSHong Zhang   }
1986d06b337dSHong Zhang 
1987d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1988d06b337dSHong Zhang 
1989d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19902798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1991d06b337dSHong Zhang       for (i=0; i<m; i++)
1992d06b337dSHong Zhang         t[i] = b[i];
1993d06b337dSHong Zhang 
1994d06b337dSHong Zhang       for (i=0; i<m; i++){
199544706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1996d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1997d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1998d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1999d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2000d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
200144706875SHong Zhang         PetscLogFlops(2*nz-1);
2002d06b337dSHong Zhang       }
2003d06b337dSHong Zhang     }
2004d06b337dSHong Zhang 
20052798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2006d06b337dSHong Zhang       for (i=0; i<m; i++)
2007d06b337dSHong Zhang         t[i] = b[i];
2008d06b337dSHong Zhang 
2009d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2010d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2011d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2012d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2013d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
201444706875SHong Zhang         PetscLogFlops(2*nz-1);
2015d06b337dSHong Zhang       }
2016d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2017d06b337dSHong Zhang         d  = *(aa + ai[i]);
2018d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2019d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2020d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2021d06b337dSHong Zhang         sum = t[i];
2022d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
202344706875SHong Zhang         PetscLogFlops(2*nz-1);
2024d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2025d06b337dSHong Zhang       }
2026d06b337dSHong Zhang     }
2027d06b337dSHong Zhang     its--;
2028d06b337dSHong Zhang   }
2029d06b337dSHong Zhang 
2030d06b337dSHong Zhang   while (its--) {
203144706875SHong Zhang     /*
203244706875SHong Zhang        forward sweep:
203344706875SHong Zhang        for i=0,...,m-1:
203444706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
203544706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
203644706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2037d06b337dSHong Zhang 
203844706875SHong Zhang     */
20392798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2040d06b337dSHong Zhang       for (i=0; i<m; i++)
2041d06b337dSHong Zhang         t[i] = b[i];
2042d06b337dSHong Zhang 
2043d06b337dSHong Zhang       for (i=0; i<m; i++){
204444706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2045d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2046d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2047d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2048d06b337dSHong Zhang         sum = t[i];
2049d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2050d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2051d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
205244706875SHong Zhang         PetscLogFlops(4*nz-2);
2053d06b337dSHong Zhang       }
2054d06b337dSHong Zhang     }
2055d06b337dSHong Zhang 
20562798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
205744706875SHong Zhang       /*
205844706875SHong Zhang        backward sweep:
205944706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
206044706875SHong Zhang        for i=m-1,...,0:
206144706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
206244706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
206344706875SHong Zhang       */
2064d06b337dSHong Zhang       for (i=0; i<m; i++)
2065d06b337dSHong Zhang         t[i] = b[i];
2066d06b337dSHong Zhang 
2067d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2068d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2069d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2070d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2071d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
207244706875SHong Zhang         PetscLogFlops(2*nz-1);
2073d06b337dSHong Zhang       }
2074d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2075d06b337dSHong Zhang         d  = *(aa + ai[i]);
2076d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2077d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2078d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2079d06b337dSHong Zhang         sum = t[i];
2080d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
208144706875SHong Zhang         PetscLogFlops(2*nz-1);
2082d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2083d06b337dSHong Zhang       }
2084d06b337dSHong Zhang     }
2085d06b337dSHong Zhang   }
2086d06b337dSHong Zhang 
2087d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
2088d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2089d06b337dSHong Zhang   if (bb != xx) {
2090d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2091d06b337dSHong Zhang   }
20922798e883SHong Zhang 
2093d06b337dSHong Zhang   PetscFunctionReturn(0);
2094d06b337dSHong Zhang }
2095d06b337dSHong Zhang 
2096d06b337dSHong Zhang 
2097d06b337dSHong Zhang 
2098d06b337dSHong Zhang 
209949b5e25fSSatish Balay 
210049b5e25fSSatish Balay 
2101