xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision d59c15a76c23fdab77d951e0cfe3025ce3c8d81c)
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"
5149b5e25fSSatish Balay 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"
74045c9aa0SHong Zhang 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);}
124*d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
125*d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
126*d59c15a7SBarry 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);}
1281a3463dfSHong Zhang 
1291a3463dfSHong Zhang   if (a->inew){
1301a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1311a3463dfSHong Zhang     a->inew = 0;
1321a3463dfSHong Zhang   }
13349b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13449b5e25fSSatish Balay   PetscFunctionReturn(0);
13549b5e25fSSatish Balay }
13649b5e25fSSatish Balay 
1374a2ae208SSatish Balay #undef __FUNCT__
1384a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
13949b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14049b5e25fSSatish Balay {
141045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14249b5e25fSSatish Balay 
14349b5e25fSSatish Balay   PetscFunctionBegin;
1444d9d31abSKris Buschelman   switch (op) {
1454d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1464d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1474d9d31abSKris Buschelman     break;
1484d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1494d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1504d9d31abSKris Buschelman     break;
1514d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1524d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1534d9d31abSKris Buschelman     break;
1544d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1554d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1564d9d31abSKris Buschelman     break;
1574d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1584d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1594d9d31abSKris Buschelman     break;
1604d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1614d9d31abSKris Buschelman     a->nonew          = 1;
1624d9d31abSKris Buschelman     break;
1634d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1644d9d31abSKris Buschelman     a->nonew          = -1;
1654d9d31abSKris Buschelman     break;
1664d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1674d9d31abSKris Buschelman     a->nonew          = -2;
1684d9d31abSKris Buschelman     break;
1694d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1704d9d31abSKris Buschelman     a->nonew          = 0;
1714d9d31abSKris Buschelman     break;
1724d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1734d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1744d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1754d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1764d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
177d03495bdSKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
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 
204c464158bSHong Zhang   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
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"
54887828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,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;
55929bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
560c464158bSHong Zhang     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
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 */
56429bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
565c464158bSHong Zhang       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
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"
59487828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,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;
6010880e062SHong Zhang   MatScalar   *value = v,*ap,*aa = a->a,*bap;
6020880e062SHong Zhang 
60349b5e25fSSatish Balay   PetscFunctionBegin;
6040880e062SHong Zhang   if (roworiented) {
6050880e062SHong Zhang     stepval = (n-1)*bs;
6060880e062SHong Zhang   } else {
6070880e062SHong Zhang     stepval = (m-1)*bs;
6080880e062SHong Zhang   }
6090880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
6100880e062SHong Zhang     row  = im[k];
6110880e062SHong Zhang     if (row < 0) continue;
6120880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
6130880e062SHong Zhang     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
6140880e062SHong Zhang #endif
6150880e062SHong Zhang     rp   = aj + ai[row];
6160880e062SHong Zhang     ap   = aa + bs2*ai[row];
6170880e062SHong Zhang     rmax = imax[row];
6180880e062SHong Zhang     nrow = ailen[row];
6190880e062SHong Zhang     low  = 0;
6200880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
6210880e062SHong Zhang       if (in[l] < 0) continue;
6220880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
6230880e062SHong Zhang       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
6240880e062SHong Zhang #endif
6250880e062SHong Zhang       col = in[l];
6260880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
6270880e062SHong Zhang       if (roworiented) {
6280880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
6290880e062SHong Zhang       } else {
6300880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
6310880e062SHong Zhang       }
6320880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
6330880e062SHong Zhang       while (high-low > 7) {
6340880e062SHong Zhang         t = (low+high)/2;
6350880e062SHong Zhang         if (rp[t] > col) high = t;
6360880e062SHong Zhang         else             low  = t;
6370880e062SHong Zhang       }
6380880e062SHong Zhang       for (i=low; i<high; i++) {
6390880e062SHong Zhang         if (rp[i] > col) break;
6400880e062SHong Zhang         if (rp[i] == col) {
6410880e062SHong Zhang           bap  = ap +  bs2*i;
6420880e062SHong Zhang           if (roworiented) {
6430880e062SHong Zhang             if (is == ADD_VALUES) {
6440880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6450880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6460880e062SHong Zhang                   bap[jj] += *value++;
6470880e062SHong Zhang                 }
6480880e062SHong Zhang               }
6490880e062SHong Zhang             } else {
6500880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6510880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6520880e062SHong Zhang                   bap[jj] = *value++;
6530880e062SHong Zhang                 }
6540880e062SHong Zhang               }
6550880e062SHong Zhang             }
6560880e062SHong Zhang           } else {
6570880e062SHong Zhang             if (is == ADD_VALUES) {
6580880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6590880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6600880e062SHong Zhang                   *bap++ += *value++;
6610880e062SHong Zhang                 }
6620880e062SHong Zhang               }
6630880e062SHong Zhang             } else {
6640880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6650880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6660880e062SHong Zhang                   *bap++  = *value++;
6670880e062SHong Zhang                 }
6680880e062SHong Zhang               }
6690880e062SHong Zhang             }
6700880e062SHong Zhang           }
6710880e062SHong Zhang           goto noinsert2;
6720880e062SHong Zhang         }
6730880e062SHong Zhang       }
6740880e062SHong Zhang       if (nonew == 1) goto noinsert2;
6750880e062SHong Zhang       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
6760880e062SHong Zhang       if (nrow >= rmax) {
6770880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
6780880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
6790880e062SHong Zhang         MatScalar *new_a;
6800880e062SHong Zhang 
6810880e062SHong Zhang         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
6820880e062SHong Zhang 
6830880e062SHong Zhang         /* malloc new storage space */
6840880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
6850880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
6860880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
6870880e062SHong Zhang         new_i   = new_j + new_nz;
6880880e062SHong Zhang 
6890880e062SHong Zhang         /* copy over old data into new slots */
6900880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
6910880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
6920880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
6930880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
6940880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
6950880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
6960880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
6970880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
6980880e062SHong Zhang         /* free up old matrix storage */
6990880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
7000880e062SHong Zhang         if (!a->singlemalloc) {
7010880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
7020880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
7030880e062SHong Zhang         }
7040880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
7050880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
7060880e062SHong Zhang 
7070880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
7080880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
7090880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
7100880e062SHong Zhang         a->s_maxnz += bs2*CHUNKSIZE;
7110880e062SHong Zhang         a->reallocs++;
7120880e062SHong Zhang         a->s_nz++;
7130880e062SHong Zhang       }
7140880e062SHong Zhang       N = nrow++ - 1;
7150880e062SHong Zhang       /* shift up all the later entries in this row */
7160880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
7170880e062SHong Zhang         rp[ii+1] = rp[ii];
7180880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
7190880e062SHong Zhang       }
7200880e062SHong Zhang       if (N >= i) {
7210880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
7220880e062SHong Zhang       }
7230880e062SHong Zhang       rp[i] = col;
7240880e062SHong Zhang       bap   = ap +  bs2*i;
7250880e062SHong Zhang       if (roworiented) {
7260880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7270880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
7280880e062SHong Zhang             bap[jj] = *value++;
7290880e062SHong Zhang           }
7300880e062SHong Zhang         }
7310880e062SHong Zhang       } else {
7320880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
7330880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
7340880e062SHong Zhang             *bap++  = *value++;
7350880e062SHong Zhang           }
7360880e062SHong Zhang         }
7370880e062SHong Zhang       }
7380880e062SHong Zhang       noinsert2:;
7390880e062SHong Zhang       low = i;
7400880e062SHong Zhang     }
7410880e062SHong Zhang     ailen[row] = nrow;
7420880e062SHong Zhang   }
7430880e062SHong Zhang   PetscFunctionReturn(0);
74449b5e25fSSatish Balay }
74549b5e25fSSatish Balay 
7464a2ae208SSatish Balay #undef __FUNCT__
7474a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
74849b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
74949b5e25fSSatish Balay {
75049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
75149b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
752c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
75349b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
75449b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
755cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES)
756cf4441caSHong Zhang   PetscTruth   flag;
757cf4441caSHong Zhang #endif
75849b5e25fSSatish Balay 
75949b5e25fSSatish Balay   PetscFunctionBegin;
76049b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
76149b5e25fSSatish Balay 
76249b5e25fSSatish Balay   if (m) rmax = ailen[0];
76349b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
76449b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
76549b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
76649b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
76749b5e25fSSatish Balay     if (fshift) {
76849b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
76949b5e25fSSatish Balay       N = ailen[i];
77049b5e25fSSatish Balay       for (j=0; j<N; j++) {
77149b5e25fSSatish Balay         ip[j-fshift] = ip[j];
77249b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
77349b5e25fSSatish Balay       }
77449b5e25fSSatish Balay     }
77549b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
77649b5e25fSSatish Balay   }
77749b5e25fSSatish Balay   if (mbs) {
77849b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
77949b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
78049b5e25fSSatish Balay   }
78149b5e25fSSatish Balay   /* reset ilen and imax for each row */
78249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
78349b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
78449b5e25fSSatish Balay   }
78549b5e25fSSatish Balay   a->s_nz = ai[mbs];
78649b5e25fSSatish Balay 
787b424e231SHong Zhang   /* diagonals may have moved, reset it */
788b424e231SHong Zhang   if (a->diag) {
789b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
79049b5e25fSSatish Balay   }
791b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
792c464158bSHong Zhang            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
793b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
79449b5e25fSSatish Balay            a->reallocs);
795b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
79649b5e25fSSatish Balay   a->reallocs          = 0;
79749b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
79849b5e25fSSatish Balay 
799cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES)
8002c535e4dSHong Zhang   ierr = PetscOptionsHasName(A->prefix,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr);
801cf4441caSHong Zhang   if (flag) { ierr = MatUseSpooles_SeqSBAIJ(A);CHKERRQ(ierr); }
802cf4441caSHong Zhang #endif
803cf4441caSHong Zhang 
80449b5e25fSSatish Balay   PetscFunctionReturn(0);
80549b5e25fSSatish Balay }
80649b5e25fSSatish Balay 
80749b5e25fSSatish Balay /*
80849b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
80949b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
81049b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
81149b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
81249b5e25fSSatish Balay */
8134a2ae208SSatish Balay #undef __FUNCT__
8144a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
815db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
81649b5e25fSSatish Balay {
81749b5e25fSSatish Balay   int        i,j,k,row;
81849b5e25fSSatish Balay   PetscTruth flg;
81949b5e25fSSatish Balay 
82049b5e25fSSatish Balay   PetscFunctionBegin;
82149b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
82249b5e25fSSatish Balay     row = idx[i];
82349b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
82449b5e25fSSatish Balay       sizes[j] = 1;
82549b5e25fSSatish Balay       i++;
82649b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
82749b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
82849b5e25fSSatish Balay       i++;
82949b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
83049b5e25fSSatish Balay       flg = PETSC_TRUE;
83149b5e25fSSatish Balay       for (k=1; k<bs; k++) {
83249b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
83349b5e25fSSatish Balay           flg = PETSC_FALSE;
83449b5e25fSSatish Balay           break;
83549b5e25fSSatish Balay         }
83649b5e25fSSatish Balay       }
83749b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
83849b5e25fSSatish Balay         sizes[j] = bs;
83949b5e25fSSatish Balay         i+= bs;
84049b5e25fSSatish Balay       } else {
84149b5e25fSSatish Balay         sizes[j] = 1;
84249b5e25fSSatish Balay         i++;
84349b5e25fSSatish Balay       }
84449b5e25fSSatish Balay     }
84549b5e25fSSatish Balay   }
84649b5e25fSSatish Balay   *bs_max = j;
84749b5e25fSSatish Balay   PetscFunctionReturn(0);
84849b5e25fSSatish Balay }
84949b5e25fSSatish Balay 
8504a2ae208SSatish Balay #undef __FUNCT__
8514a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
85287828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag)
85349b5e25fSSatish Balay {
85449b5e25fSSatish Balay   PetscFunctionBegin;
855c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
85649b5e25fSSatish Balay }
85749b5e25fSSatish Balay 
85849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
85949b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
86049b5e25fSSatish Balay */
86149b5e25fSSatish Balay 
8624a2ae208SSatish Balay #undef __FUNCT__
8634a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
86487828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
86549b5e25fSSatish Balay {
86649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
86749b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
86849b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
86949b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
87049b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
87149b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
87249b5e25fSSatish Balay 
87349b5e25fSSatish Balay   PetscFunctionBegin;
87449b5e25fSSatish Balay 
87549b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
87649b5e25fSSatish Balay     row  = im[k];       /* row number */
87749b5e25fSSatish Balay     brow = row/bs;      /* block row number */
87849b5e25fSSatish Balay     if (row < 0) continue;
87949b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
880c464158bSHong Zhang     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
88149b5e25fSSatish Balay #endif
88249b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
88349b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
88449b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
88549b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
88649b5e25fSSatish Balay     low  = 0;
88749b5e25fSSatish Balay 
88849b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
88949b5e25fSSatish Balay       if (in[l] < 0) continue;
89049b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
891c464158bSHong Zhang       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m);
89249b5e25fSSatish Balay #endif
89349b5e25fSSatish Balay       col = in[l];
89449b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
89549b5e25fSSatish Balay 
89649b5e25fSSatish Balay       if (brow <= bcol){
89749b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8988549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
89949b5e25fSSatish Balay           /* element value a(k,l) */
90049b5e25fSSatish Balay           if (roworiented) {
90149b5e25fSSatish Balay             value = v[l + k*n];
90249b5e25fSSatish Balay           } else {
90349b5e25fSSatish Balay             value = v[k + l*m];
90449b5e25fSSatish Balay           }
90549b5e25fSSatish Balay 
90649b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
90749b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
90849b5e25fSSatish Balay           while (high-low > 7) {
90949b5e25fSSatish Balay             t = (low+high)/2;
91049b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
91149b5e25fSSatish Balay             else              low  = t;
91249b5e25fSSatish Balay           }
91349b5e25fSSatish Balay           for (i=low; i<high; i++) {
91449b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
91549b5e25fSSatish Balay             if (rp[i] > bcol) break;
91649b5e25fSSatish Balay             if (rp[i] == bcol) {
91749b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
91849b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
91949b5e25fSSatish Balay               else                  *bap  = value;
9208549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9218549e402SHong Zhang               if (brow == bcol && ridx < cidx){
9228549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
9238549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
9248549e402SHong Zhang                 else                  *bap  = value;
9258549e402SHong Zhang               }
92649b5e25fSSatish Balay               goto noinsert1;
92749b5e25fSSatish Balay             }
92849b5e25fSSatish Balay           }
92949b5e25fSSatish Balay 
93049b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
93129bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
93249b5e25fSSatish Balay       if (nrow >= rmax) {
93349b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
93449b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
93549b5e25fSSatish Balay         MatScalar *new_a;
93649b5e25fSSatish Balay 
93729bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
93849b5e25fSSatish Balay 
93949b5e25fSSatish Balay         /* Malloc new storage space */
94049b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
94182502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
94249b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
94349b5e25fSSatish Balay         new_i = new_j + new_nz;
94449b5e25fSSatish Balay 
94549b5e25fSSatish Balay         /* copy over old data into new slots */
94649b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
94749b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
94849b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
94949b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
95049b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
95149b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
95249b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
95349b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
95449b5e25fSSatish Balay         /* free up old matrix storage */
95549b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
95649b5e25fSSatish Balay         if (!a->singlemalloc) {
95749b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
95849b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
95949b5e25fSSatish Balay         }
96049b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
96149b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
96249b5e25fSSatish Balay 
96349b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
96449b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
965b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
96649b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
96749b5e25fSSatish Balay         a->reallocs++;
96849b5e25fSSatish Balay         a->s_nz++;
96949b5e25fSSatish Balay       }
97049b5e25fSSatish Balay 
97149b5e25fSSatish Balay       N = nrow++ - 1;
97249b5e25fSSatish Balay       /* shift up all the later entries in this row */
97349b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
97449b5e25fSSatish Balay         rp[ii+1] = rp[ii];
97549b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
97649b5e25fSSatish Balay       }
97749b5e25fSSatish Balay       if (N>=i) {
97849b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
97949b5e25fSSatish Balay       }
98049b5e25fSSatish Balay       rp[i]                      = bcol;
98149b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
98249b5e25fSSatish Balay       noinsert1:;
98349b5e25fSSatish Balay       low = i;
98449b5e25fSSatish Balay       /* } */
9858549e402SHong Zhang         }
98649b5e25fSSatish Balay       } /* end of if .. if..  */
98749b5e25fSSatish Balay     }                     /* end of loop over added columns */
98849b5e25fSSatish Balay     ailen[brow] = nrow;
98949b5e25fSSatish Balay   }                       /* end of loop over added rows */
99049b5e25fSSatish Balay 
99149b5e25fSSatish Balay   PetscFunctionReturn(0);
99249b5e25fSSatish Balay }
99349b5e25fSSatish Balay 
994915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
995915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
99649b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
99749b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
99849b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
99949b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
100049b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
100187828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat);
100249b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
100349b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
100449b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
100549b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
100649b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
100749b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
100813a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
100949b5e25fSSatish Balay 
101049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
101149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
101249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
101349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
101449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
101549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
101649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
101749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
101849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
101949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
102049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
102149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
102249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
102349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
102449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
102549b5e25fSSatish Balay 
1026*d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1027*d59c15a7SBarry Smith 
102807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
102907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
103007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
103107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
103207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
103307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
103407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
103507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
103607f98182SSatish Balay 
10375f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
10385f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
10395f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
10405f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
10415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
10425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
10435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
104457d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
10453e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
104649b5e25fSSatish Balay 
104749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
104849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
104949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
105049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
105149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
105249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
105349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
105449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
105549b5e25fSSatish Balay 
105649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
105749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
105849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
105949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
106049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
106149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
106249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
106349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
106449b5e25fSSatish Balay 
10654d101231SSatish Balay #ifdef HAVE_MatICCFactor
10661a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
10674a2ae208SSatish Balay #undef __FUNCT__
10684d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
10694d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level)
107049b5e25fSSatish Balay {
10714ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
107249b5e25fSSatish Balay   Mat         outA;
107349b5e25fSSatish Balay   int         ierr;
107449b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
107549b5e25fSSatish Balay 
107649b5e25fSSatish Balay   PetscFunctionBegin;
10771a3463dfSHong Zhang   /*
107829bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
107949b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
108049b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
108149b5e25fSSatish Balay   if (!row_identity || !col_identity) {
108229bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
108349b5e25fSSatish Balay   }
10841a3463dfSHong Zhang   */
108549b5e25fSSatish Balay 
108649b5e25fSSatish Balay   outA          = inA;
10871260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
108849b5e25fSSatish Balay 
108949b5e25fSSatish Balay   if (!a->diag) {
10901a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
109149b5e25fSSatish Balay   }
109249b5e25fSSatish Balay   /*
109349b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
109449b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
109549b5e25fSSatish Balay   */
109649b5e25fSSatish Balay   switch (a->bs) {
109749b5e25fSSatish Balay   case 1:
10980fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
109907f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1100*d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
11014d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
110249b5e25fSSatish Balay   case 2:
11031a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
11041a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11051a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11064d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
110749b5e25fSSatish Balay     break;
110849b5e25fSSatish Balay   case 3:
11091a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
11101a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11111a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11124d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
111349b5e25fSSatish Balay     break;
111449b5e25fSSatish Balay   case 4:
11151a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
11161a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11171a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11184d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
111949b5e25fSSatish Balay     break;
112049b5e25fSSatish Balay   case 5:
11211a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
11221a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11231a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11244d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
112549b5e25fSSatish Balay     break;
112649b5e25fSSatish Balay   case 6:
11271a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
11281a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11291a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11304d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
113149b5e25fSSatish Balay     break;
113249b5e25fSSatish Balay   case 7:
11331a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
11341a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11351a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11364d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
113749b5e25fSSatish Balay     break;
113849b5e25fSSatish Balay   default:
113949b5e25fSSatish Balay     a->row        = row;
11401a3463dfSHong Zhang     a->icol       = col;
114149b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
114249b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
114349b5e25fSSatish Balay 
114449b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
114549b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1146b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
114749b5e25fSSatish Balay 
114849b5e25fSSatish Balay     if (!a->solve_work) {
114987828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
115087828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
115149b5e25fSSatish Balay     }
115249b5e25fSSatish Balay   }
115349b5e25fSSatish Balay 
11541a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay   PetscFunctionReturn(0);
115749b5e25fSSatish Balay }
1158950f1e5bSHong Zhang #endif
1159950f1e5bSHong Zhang 
11604a2ae208SSatish Balay #undef __FUNCT__
11614a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
116249b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
116349b5e25fSSatish Balay {
116449b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
116549b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
116649b5e25fSSatish Balay   int               ierr;
116749b5e25fSSatish Balay 
116849b5e25fSSatish Balay   PetscFunctionBegin;
116949b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
117049b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
117149b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
117249b5e25fSSatish Balay   PetscFunctionReturn(0);
117349b5e25fSSatish Balay }
117449b5e25fSSatish Balay 
117549b5e25fSSatish Balay EXTERN_C_BEGIN
11764a2ae208SSatish Balay #undef __FUNCT__
11774a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
117849b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
117949b5e25fSSatish Balay {
1180045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
118149b5e25fSSatish Balay   int         i,nz,n;
118249b5e25fSSatish Balay 
118349b5e25fSSatish Balay   PetscFunctionBegin;
1184045c9aa0SHong Zhang   nz = baij->s_maxnz;
1185c464158bSHong Zhang   n  = mat->n;
118649b5e25fSSatish Balay   for (i=0; i<nz; i++) {
118749b5e25fSSatish Balay     baij->j[i] = indices[i];
118849b5e25fSSatish Balay   }
1189045c9aa0SHong Zhang   baij->s_nz = nz;
119049b5e25fSSatish Balay   for (i=0; i<n; i++) {
119149b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
119249b5e25fSSatish Balay   }
119349b5e25fSSatish Balay 
119449b5e25fSSatish Balay   PetscFunctionReturn(0);
119549b5e25fSSatish Balay }
119649b5e25fSSatish Balay EXTERN_C_END
119749b5e25fSSatish Balay 
11984a2ae208SSatish Balay #undef __FUNCT__
11994a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
120049b5e25fSSatish Balay /*@
120119585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
120249b5e25fSSatish Balay        in the matrix.
120349b5e25fSSatish Balay 
120449b5e25fSSatish Balay   Input Parameters:
120519585528SSatish Balay +  mat     - the SeqSBAIJ matrix
120649b5e25fSSatish Balay -  indices - the column indices
120749b5e25fSSatish Balay 
120849b5e25fSSatish Balay   Level: advanced
120949b5e25fSSatish Balay 
121049b5e25fSSatish Balay   Notes:
121149b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
121249b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
121349b5e25fSSatish Balay   of the MatSetValues() operation.
121449b5e25fSSatish Balay 
121549b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
121619585528SSatish Balay   MatCreateSeqSBAIJ().
121749b5e25fSSatish Balay 
1218ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
121949b5e25fSSatish Balay 
1220ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
122149b5e25fSSatish Balay @*/
122249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
122349b5e25fSSatish Balay {
122449b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
122549b5e25fSSatish Balay 
122649b5e25fSSatish Balay   PetscFunctionBegin;
122749b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1228c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
122949b5e25fSSatish Balay   if (f) {
123049b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
123149b5e25fSSatish Balay   } else {
123229bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
123349b5e25fSSatish Balay   }
123449b5e25fSSatish Balay   PetscFunctionReturn(0);
123549b5e25fSSatish Balay }
123649b5e25fSSatish Balay 
12374a2ae208SSatish Balay #undef __FUNCT__
12384a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1239273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1240273d9f13SBarry Smith {
1241273d9f13SBarry Smith   int        ierr;
1242273d9f13SBarry Smith 
1243273d9f13SBarry Smith   PetscFunctionBegin;
1244273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1245273d9f13SBarry Smith   PetscFunctionReturn(0);
1246273d9f13SBarry Smith }
1247273d9f13SBarry Smith 
1248a6ece127SHong Zhang #undef __FUNCT__
1249a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1250a6ece127SHong Zhang int MatGetArray_SeqSBAIJ(Mat A,PetscScalar **array)
1251a6ece127SHong Zhang {
1252a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1253a6ece127SHong Zhang   PetscFunctionBegin;
1254a6ece127SHong Zhang   *array = a->a;
1255a6ece127SHong Zhang   PetscFunctionReturn(0);
1256a6ece127SHong Zhang }
1257a6ece127SHong Zhang 
1258a6ece127SHong Zhang #undef __FUNCT__
1259a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1260a6ece127SHong Zhang int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar **array)
1261a6ece127SHong Zhang {
1262a6ece127SHong Zhang   PetscFunctionBegin;
1263a6ece127SHong Zhang   PetscFunctionReturn(0);
1264a6ece127SHong Zhang }
1265a6ece127SHong Zhang 
126642ee4b1aSHong Zhang #include "petscblaslapack.h"
126742ee4b1aSHong Zhang #undef __FUNCT__
126842ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
126942ee4b1aSHong Zhang int MatAXPY_SeqSBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
127042ee4b1aSHong Zhang {
127142ee4b1aSHong Zhang   int          ierr,one=1;
127242ee4b1aSHong Zhang   Mat_SeqSBAIJ *x  = (Mat_SeqSBAIJ *)X->data,*y = (Mat_SeqSBAIJ *)Y->data;
127342ee4b1aSHong Zhang 
127442ee4b1aSHong Zhang   PetscFunctionBegin;
127542ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
127642ee4b1aSHong Zhang     BLaxpy_(&x->s_nz,a,x->a,&one,y->a,&one);
127742ee4b1aSHong Zhang   } else {
127842ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
127942ee4b1aSHong Zhang   }
128042ee4b1aSHong Zhang   PetscFunctionReturn(0);
128142ee4b1aSHong Zhang }
128242ee4b1aSHong Zhang 
128349b5e25fSSatish Balay /* -------------------------------------------------------------------*/
128449b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
128549b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
128649b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
128749b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
128849b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
128949b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
129049b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
129149b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
129249b5e25fSSatish Balay        0,
129349b5e25fSSatish Balay        0,
129449b5e25fSSatish Balay        0,
129549b5e25fSSatish Balay        0,
12965f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1297d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
129849b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
129949b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
130049b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
130149b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
130249b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
130349b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
130449b5e25fSSatish Balay        0,
130549b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
130649b5e25fSSatish Balay        0,
130749b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
130849b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
130949b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
131049b5e25fSSatish Balay        0,
131149b5e25fSSatish Balay        0,
13125f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13135f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
1314273d9f13SBarry Smith        MatSetUpPreallocation_SeqSBAIJ,
1315c464158bSHong Zhang        0,
13164d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1317a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1318a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
131949b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
132049b5e25fSSatish Balay        0,
132149b5e25fSSatish Balay        0,
132249b5e25fSSatish Balay        0,
1323950f1e5bSHong Zhang        0,
132442ee4b1aSHong Zhang        MatAXPY_SeqSBAIJ,
132549b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
132649b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
132749b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
132849b5e25fSSatish Balay        0,
132949b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
133049b5e25fSSatish Balay        MatScale_SeqSBAIJ,
133149b5e25fSSatish Balay        0,
133249b5e25fSSatish Balay        0,
133349b5e25fSSatish Balay        0,
133449b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
133549b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
133649b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
133749b5e25fSSatish Balay        0,
133849b5e25fSSatish Balay        0,
133949b5e25fSSatish Balay        0,
134049b5e25fSSatish Balay        0,
134149b5e25fSSatish Balay        0,
134249b5e25fSSatish Balay        0,
134349b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
134449b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
134549b5e25fSSatish Balay        0,
134649b5e25fSSatish Balay        0,
13478a124369SBarry Smith        MatGetPetscMaps_Petsc,
1348d959ec07SHong Zhang        0,
1349d959ec07SHong Zhang        0,
1350d959ec07SHong Zhang        0,
1351d959ec07SHong Zhang        0,
1352d959ec07SHong Zhang        0,
1353d959ec07SHong Zhang        0,
13543e0d88b5SBarry Smith        MatGetRowMax_SeqSBAIJ,
13553e0d88b5SBarry Smith        0,
13563e0d88b5SBarry Smith        0,
13573e0d88b5SBarry Smith        0,
13583e0d88b5SBarry Smith        0,
13593e0d88b5SBarry Smith        0,
13603e0d88b5SBarry Smith        0,
13613e0d88b5SBarry Smith        0,
13623e0d88b5SBarry Smith        0,
13633e0d88b5SBarry Smith        0,
13643e0d88b5SBarry Smith        0,
13653e0d88b5SBarry Smith        0,
13663e0d88b5SBarry Smith        0,
13673e0d88b5SBarry Smith        0,
13683e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
13693e0d88b5SBarry Smith        MatGetInertia_SeqSBAIJ
13703e0d88b5SBarry Smith #else
13713e0d88b5SBarry Smith        0
13723e0d88b5SBarry Smith #endif
13733e0d88b5SBarry Smith };
137449b5e25fSSatish Balay 
137549b5e25fSSatish Balay EXTERN_C_BEGIN
13764a2ae208SSatish Balay #undef __FUNCT__
13774a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
137849b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
137949b5e25fSSatish Balay {
13804afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1381c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
138249b5e25fSSatish Balay   int         ierr;
138349b5e25fSSatish Balay 
138449b5e25fSSatish Balay   PetscFunctionBegin;
138549b5e25fSSatish Balay   if (aij->nonew != 1) {
138629bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
138749b5e25fSSatish Balay   }
138849b5e25fSSatish Balay 
138949b5e25fSSatish Balay   /* allocate space for values if not already there */
139049b5e25fSSatish Balay   if (!aij->saved_values) {
139187828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
139249b5e25fSSatish Balay   }
139349b5e25fSSatish Balay 
139449b5e25fSSatish Balay   /* copy values over */
139587828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
139649b5e25fSSatish Balay   PetscFunctionReturn(0);
139749b5e25fSSatish Balay }
139849b5e25fSSatish Balay EXTERN_C_END
139949b5e25fSSatish Balay 
140049b5e25fSSatish Balay EXTERN_C_BEGIN
14014a2ae208SSatish Balay #undef __FUNCT__
14024a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
140349b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
140449b5e25fSSatish Balay {
14054afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1406c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
140749b5e25fSSatish Balay 
140849b5e25fSSatish Balay   PetscFunctionBegin;
140949b5e25fSSatish Balay   if (aij->nonew != 1) {
141029bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
141149b5e25fSSatish Balay   }
141249b5e25fSSatish Balay   if (!aij->saved_values) {
141329bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
141449b5e25fSSatish Balay   }
141549b5e25fSSatish Balay 
141649b5e25fSSatish Balay   /* copy values over */
141787828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
141849b5e25fSSatish Balay   PetscFunctionReturn(0);
141949b5e25fSSatish Balay }
142049b5e25fSSatish Balay EXTERN_C_END
142149b5e25fSSatish Balay 
14228549e402SHong Zhang EXTERN_C_BEGIN
14234a2ae208SSatish Balay #undef __FUNCT__
14244a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ"
1425c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B)
1426c464158bSHong Zhang {
1427c464158bSHong Zhang   Mat_SeqSBAIJ *b;
14280c955e93SHong Zhang   int          ierr,size;
1429c464158bSHong Zhang 
1430c464158bSHong Zhang   PetscFunctionBegin;
1431c464158bSHong Zhang   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1432c464158bSHong Zhang   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1433273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1434273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1435c464158bSHong Zhang 
1436b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1437b0a32e0cSBarry Smith   B->data = (void*)b;
1438c464158bSHong Zhang   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1439c464158bSHong Zhang   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1440c464158bSHong Zhang   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1441c464158bSHong Zhang   B->ops->view        = MatView_SeqSBAIJ;
1442c464158bSHong Zhang   B->factor           = 0;
1443c464158bSHong Zhang   B->lupivotthreshold = 1.0;
1444c464158bSHong Zhang   B->mapping          = 0;
1445c464158bSHong Zhang   b->row              = 0;
1446c464158bSHong Zhang   b->icol             = 0;
1447c464158bSHong Zhang   b->reallocs         = 0;
1448c464158bSHong Zhang   b->saved_values     = 0;
1449c464158bSHong Zhang 
14508a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
14518a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1452c464158bSHong Zhang 
1453c464158bSHong Zhang   b->sorted           = PETSC_FALSE;
1454c464158bSHong Zhang   b->roworiented      = PETSC_TRUE;
1455c464158bSHong Zhang   b->nonew            = 0;
1456c464158bSHong Zhang   b->diag             = 0;
1457c464158bSHong Zhang   b->solve_work       = 0;
1458c464158bSHong Zhang   b->mult_work        = 0;
14592a1b7f2aSHong Zhang   B->spptr            = 0;
1460c464158bSHong Zhang   b->keepzeroedrows   = PETSC_FALSE;
1461c464158bSHong Zhang 
1462c464158bSHong Zhang   b->inew             = 0;
1463c464158bSHong Zhang   b->jnew             = 0;
1464c464158bSHong Zhang   b->anew             = 0;
1465c464158bSHong Zhang   b->a2anew           = 0;
1466c464158bSHong Zhang   b->permute          = PETSC_FALSE;
1467c464158bSHong Zhang 
1468c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1469c464158bSHong Zhang                                      "MatStoreValues_SeqSBAIJ",
1470c464158bSHong Zhang                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1471c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1472c464158bSHong Zhang                                      "MatRetrieveValues_SeqSBAIJ",
1473c464158bSHong Zhang                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1474c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1475c464158bSHong Zhang                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1476c464158bSHong Zhang                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1477c464158bSHong Zhang   PetscFunctionReturn(0);
1478c464158bSHong Zhang }
14798549e402SHong Zhang EXTERN_C_END
1480c464158bSHong Zhang 
14814a2ae208SSatish Balay #undef __FUNCT__
14824a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
148349b5e25fSSatish Balay /*@C
1484273d9f13SBarry Smith    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
148549b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
148649b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
148749b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
148849b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
148949b5e25fSSatish Balay 
1490c464158bSHong Zhang    Collective on Mat
149149b5e25fSSatish Balay 
149249b5e25fSSatish Balay    Input Parameters:
1493c464158bSHong Zhang +  A - the symmetric matrix
149449b5e25fSSatish Balay .  bs - size of block
149549b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
1496744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1497744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
149849b5e25fSSatish Balay 
149949b5e25fSSatish Balay    Options Database Keys:
150049b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
150149b5e25fSSatish Balay                      block calculations (much slower)
150249b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
150349b5e25fSSatish Balay 
150449b5e25fSSatish Balay    Level: intermediate
150549b5e25fSSatish Balay 
150649b5e25fSSatish Balay    Notes:
150749b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
150849b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
150949b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
151049b5e25fSSatish Balay    matrices.
151149b5e25fSSatish Balay 
1512a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
151349b5e25fSSatish Balay @*/
1514273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
151549b5e25fSSatish Balay {
1516c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
15170c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
151849b5e25fSSatish Balay   PetscTruth   flg;
15194afc71dfSHong Zhang   int          s_nz;
152049b5e25fSSatish Balay 
152149b5e25fSSatish Balay   PetscFunctionBegin;
1522273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1523e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1524c464158bSHong Zhang   mbs  = B->m/bs;
152549b5e25fSSatish Balay   bs2  = bs*bs;
152649b5e25fSSatish Balay 
1527c464158bSHong Zhang   if (mbs*bs != B->m) {
152829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
152949b5e25fSSatish Balay   }
153049b5e25fSSatish Balay 
1531435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1532435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
153349b5e25fSSatish Balay   if (nnz) {
153449b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
153529bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
153680fe4e49SBarry 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);
153749b5e25fSSatish Balay     }
153849b5e25fSSatish Balay   }
153949b5e25fSSatish Balay 
1540e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
154149b5e25fSSatish Balay   if (!flg) {
154249b5e25fSSatish Balay     switch (bs) {
154349b5e25fSSatish Balay     case 1:
15444ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
154549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1546*d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
154749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
154849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
154949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
155049b5e25fSSatish Balay       break;
155149b5e25fSSatish Balay     case 2:
15524ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
155349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
155449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
155549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
155649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
155749b5e25fSSatish Balay       break;
155849b5e25fSSatish Balay     case 3:
15595f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
156049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
156149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
156249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
156349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
156449b5e25fSSatish Balay       break;
156549b5e25fSSatish Balay     case 4:
15665f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
156749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
156849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
156949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
157049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
157149b5e25fSSatish Balay       break;
157249b5e25fSSatish Balay     case 5:
15735f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
157449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
157549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
157649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
157749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
157849b5e25fSSatish Balay       break;
157949b5e25fSSatish Balay     case 6:
15805f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
158149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
158249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
158349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
158449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
158549b5e25fSSatish Balay       break;
158649b5e25fSSatish Balay     case 7:
1587de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
158849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
158949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1590de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
159149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
159249b5e25fSSatish Balay       break;
159349b5e25fSSatish Balay     }
159449b5e25fSSatish Balay   }
159549b5e25fSSatish Balay 
159649b5e25fSSatish Balay   b->mbs = mbs;
15974afc71dfSHong Zhang   b->nbs = mbs;
1598b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
159949b5e25fSSatish Balay   if (!nnz) {
1600435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
160149b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
160249b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
16038cef66ccSHong Zhang       b->imax[i] = nz;
160449b5e25fSSatish Balay     }
1605153ea458SHong Zhang     nz = nz*mbs; /* total nz */
160649b5e25fSSatish Balay   } else {
160749b5e25fSSatish Balay     nz = 0;
16088cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
160949b5e25fSSatish Balay   }
16108cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
16118cef66ccSHong Zhang   s_nz = nz;
161249b5e25fSSatish Balay 
161349b5e25fSSatish Balay   /* allocate the matrix space */
1614c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
161582502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
161649b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
161749b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
161849b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
161949b5e25fSSatish Balay   b->i = b->j + s_nz;
162049b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
162149b5e25fSSatish Balay 
162249b5e25fSSatish Balay   /* pointer to beginning of each row */
162349b5e25fSSatish Balay   b->i[0] = 0;
162449b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
162549b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
162649b5e25fSSatish Balay   }
162749b5e25fSSatish Balay 
162849b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
16295033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1630b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
163149b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
163249b5e25fSSatish Balay 
163349b5e25fSSatish Balay   b->bs               = bs;
163449b5e25fSSatish Balay   b->bs2              = bs2;
163549b5e25fSSatish Balay   b->s_nz             = 0;
163649b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1637153ea458SHong Zhang 
163816cdd363SHong Zhang   b->inew             = 0;
163916cdd363SHong Zhang   b->jnew             = 0;
164016cdd363SHong Zhang   b->anew             = 0;
164116cdd363SHong Zhang   b->a2anew           = 0;
16421a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1643c464158bSHong Zhang   PetscFunctionReturn(0);
1644c464158bSHong Zhang }
1645153ea458SHong Zhang 
164649b5e25fSSatish Balay 
16474a2ae208SSatish Balay #undef __FUNCT__
16484a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1649c464158bSHong Zhang /*@C
1650c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1651c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1652c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1653c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1654c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
165549b5e25fSSatish Balay 
1656c464158bSHong Zhang    Collective on MPI_Comm
1657c464158bSHong Zhang 
1658c464158bSHong Zhang    Input Parameters:
1659c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1660c464158bSHong Zhang .  bs - size of block
1661c464158bSHong Zhang .  m - number of rows, or number of columns
1662c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1663744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1664744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1665c464158bSHong Zhang 
1666c464158bSHong Zhang    Output Parameter:
1667c464158bSHong Zhang .  A - the symmetric matrix
1668c464158bSHong Zhang 
1669c464158bSHong Zhang    Options Database Keys:
1670c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1671c464158bSHong Zhang                      block calculations (much slower)
1672c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1673c464158bSHong Zhang 
1674c464158bSHong Zhang    Level: intermediate
1675c464158bSHong Zhang 
1676c464158bSHong Zhang    Notes:
1677c464158bSHong Zhang 
1678c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1679c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1680c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1681c464158bSHong Zhang    matrices.
1682c464158bSHong Zhang 
1683c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1684c464158bSHong Zhang @*/
1685c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1686c464158bSHong Zhang {
1687c464158bSHong Zhang   int ierr;
1688c464158bSHong Zhang 
1689c464158bSHong Zhang   PetscFunctionBegin;
1690c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1691c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1692273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
169349b5e25fSSatish Balay   PetscFunctionReturn(0);
169449b5e25fSSatish Balay }
169549b5e25fSSatish Balay 
16964a2ae208SSatish Balay #undef __FUNCT__
16974a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
169849b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
169949b5e25fSSatish Balay {
170049b5e25fSSatish Balay   Mat         C;
170149b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
170249b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
170349b5e25fSSatish Balay 
170449b5e25fSSatish Balay   PetscFunctionBegin;
170529bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
170649b5e25fSSatish Balay 
170749b5e25fSSatish Balay   *B = 0;
1708692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1709692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1710692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1711692f9cbeSHong Zhang 
171249b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1713273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
171449b5e25fSSatish Balay   C->factor         = A->factor;
171549b5e25fSSatish Balay   c->row            = 0;
171649b5e25fSSatish Balay   c->icol           = 0;
171749b5e25fSSatish Balay   c->saved_values   = 0;
171849b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
171949b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
172049b5e25fSSatish Balay 
172149b5e25fSSatish Balay   c->bs         = a->bs;
172249b5e25fSSatish Balay   c->bs2        = a->bs2;
172349b5e25fSSatish Balay   c->mbs        = a->mbs;
172449b5e25fSSatish Balay   c->nbs        = a->nbs;
172549b5e25fSSatish Balay 
1726b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1727b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
172849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
172949b5e25fSSatish Balay     c->imax[i] = a->imax[i];
173049b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
173149b5e25fSSatish Balay   }
173249b5e25fSSatish Balay 
173349b5e25fSSatish Balay   /* allocate the matrix space */
173449b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
173549b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
173682502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
173749b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
173849b5e25fSSatish Balay   c->i = c->j + nz;
173949b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
174049b5e25fSSatish Balay   if (mbs > 0) {
174149b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
174249b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
174349b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
174449b5e25fSSatish Balay     } else {
174549b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
174649b5e25fSSatish Balay     }
174749b5e25fSSatish Balay   }
174849b5e25fSSatish Balay 
1749b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
175049b5e25fSSatish Balay   c->sorted      = a->sorted;
175149b5e25fSSatish Balay   c->roworiented = a->roworiented;
175249b5e25fSSatish Balay   c->nonew       = a->nonew;
175349b5e25fSSatish Balay 
175449b5e25fSSatish Balay   if (a->diag) {
1755b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1756b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
175749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
175849b5e25fSSatish Balay       c->diag[i] = a->diag[i];
175949b5e25fSSatish Balay     }
176049b5e25fSSatish Balay   } else c->diag        = 0;
176149b5e25fSSatish Balay   c->s_nz               = a->s_nz;
176249b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
176349b5e25fSSatish Balay   c->solve_work         = 0;
17642a1b7f2aSHong Zhang   C->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
176549b5e25fSSatish Balay   c->mult_work          = 0;
176649b5e25fSSatish Balay   *B = C;
1767b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
176849b5e25fSSatish Balay   PetscFunctionReturn(0);
176949b5e25fSSatish Balay }
177049b5e25fSSatish Balay 
17718549e402SHong Zhang EXTERN_C_BEGIN
17724a2ae208SSatish Balay #undef __FUNCT__
17734a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1774b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
177549b5e25fSSatish Balay {
177649b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
177749b5e25fSSatish Balay   Mat          B;
177849b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
17793f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
178049b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
178149b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
178287828ca2SBarry Smith   PetscScalar  *aa;
178349b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
178449b5e25fSSatish Balay 
178549b5e25fSSatish Balay   PetscFunctionBegin;
1786b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
178749b5e25fSSatish Balay   bs2  = bs*bs;
178849b5e25fSSatish Balay 
178949b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
179029bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1791b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
179249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1793552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
179449b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
179549b5e25fSSatish Balay 
179649b5e25fSSatish Balay   if (header[3] < 0) {
179729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
179849b5e25fSSatish Balay   }
179949b5e25fSSatish Balay 
180029bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
180149b5e25fSSatish Balay 
180249b5e25fSSatish Balay   /*
180349b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
180449b5e25fSSatish Balay     divisible by the blocksize
180549b5e25fSSatish Balay   */
180649b5e25fSSatish Balay   mbs        = M/bs;
180749b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
180849b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
180949b5e25fSSatish Balay   else                  mbs++;
181049b5e25fSSatish Balay   if (extra_rows) {
1811b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
181249b5e25fSSatish Balay   }
181349b5e25fSSatish Balay 
181449b5e25fSSatish Balay   /* read in row lengths */
1815b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
181649b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
181749b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
181849b5e25fSSatish Balay 
181949b5e25fSSatish Balay   /* read in column indices */
1820b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
182149b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
182249b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
182349b5e25fSSatish Balay 
182449b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
182582502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1826a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
182782502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
182849b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
182949b5e25fSSatish Balay   masked   = mask + mbs;
183049b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
183149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
183249b5e25fSSatish Balay     nmask = 0;
183349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
183449b5e25fSSatish Balay       kmax = rowlengths[rowcount];
183549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18362d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
183703630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
183849b5e25fSSatish Balay       }
183949b5e25fSSatish Balay       rowcount++;
184049b5e25fSSatish Balay     }
1841574b2666SHong Zhang     s_browlengths[i] += nmask;
1842574b2666SHong Zhang 
184349b5e25fSSatish Balay     /* zero out the mask elements we set */
184449b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
184549b5e25fSSatish Balay   }
184649b5e25fSSatish Balay 
184749b5e25fSSatish Balay   /* create our matrix */
18487fe2be48SHong Zhang   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr);
184949b5e25fSSatish Balay   B = *A;
185049b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
185149b5e25fSSatish Balay 
185249b5e25fSSatish Balay   /* set matrix "i" values */
185349b5e25fSSatish Balay   a->i[0] = 0;
185449b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1855574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1856574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
185749b5e25fSSatish Balay   }
18587fe2be48SHong Zhang   a->s_nz = a->i[mbs];
185949b5e25fSSatish Balay 
186049b5e25fSSatish Balay   /* read in nonzero values */
186187828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
186249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
186349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
186449b5e25fSSatish Balay 
186549b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
186649b5e25fSSatish Balay   nzcount = 0; jcount = 0;
186749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
186849b5e25fSSatish Balay     nzcountb = nzcount;
186949b5e25fSSatish Balay     nmask    = 0;
187049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
187149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
187249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18732d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
187403630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
18752d703238SHong Zhang       }
18762d703238SHong Zhang     }
18772d703238SHong Zhang     /* sort the masked values */
18782d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
18792d703238SHong Zhang 
18802d703238SHong Zhang     /* set "j" values into matrix */
18812d703238SHong Zhang     maskcount = 1;
18822d703238SHong Zhang     for (j=0; j<nmask; j++) {
188349b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
188449b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
188549b5e25fSSatish Balay     }
1886574b2666SHong Zhang 
188749b5e25fSSatish Balay     /* set "a" values into matrix */
188849b5e25fSSatish Balay     ishift = bs2*a->i[i];
188949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
189049b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
189149b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1892574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1893574b2666SHong Zhang         if (tmp >= i){
189449b5e25fSSatish Balay           block     = mask[tmp] - 1;
189549b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
189649b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1897574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1898574b2666SHong Zhang         }
1899574b2666SHong Zhang         nzcountb++;
190049b5e25fSSatish Balay       }
190149b5e25fSSatish Balay     }
190249b5e25fSSatish Balay     /* zero out the mask elements we set */
190349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
190449b5e25fSSatish Balay   }
190529bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
190649b5e25fSSatish Balay 
190749b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1908574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
190949b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
191049b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
191149b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
191249b5e25fSSatish Balay 
191349b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
191449b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
191549b5e25fSSatish Balay   PetscFunctionReturn(0);
191649b5e25fSSatish Balay }
19178549e402SHong Zhang EXTERN_C_END
1918574b2666SHong Zhang 
1919d06b337dSHong Zhang #undef __FUNCT__
1920d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1921c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1922d06b337dSHong Zhang {
1923d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1924d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1925d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
1926d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
1927d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1928d06b337dSHong Zhang 
1929d06b337dSHong Zhang   PetscFunctionBegin;
1930b965ef7fSBarry Smith   its = its*lits;
193191723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1932d06b337dSHong Zhang 
1933d06b337dSHong Zhang   if (bs > 1)
1934d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1935d06b337dSHong Zhang 
1936d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1937d06b337dSHong Zhang   if (xx != bb) {
1938d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1939d06b337dSHong Zhang   } else {
1940d06b337dSHong Zhang     b = x;
1941d06b337dSHong Zhang   }
1942d06b337dSHong Zhang 
1943d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1944d06b337dSHong Zhang 
1945d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19462798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1947d06b337dSHong Zhang       for (i=0; i<m; i++)
1948d06b337dSHong Zhang         t[i] = b[i];
1949d06b337dSHong Zhang 
1950d06b337dSHong Zhang       for (i=0; i<m; i++){
195144706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1952d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1953d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1954d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1955d06b337dSHong Zhang         x[i] = omega*t[i]/d;
1956d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
195744706875SHong Zhang         PetscLogFlops(2*nz-1);
1958d06b337dSHong Zhang       }
1959d06b337dSHong Zhang     }
1960d06b337dSHong Zhang 
19612798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1962d06b337dSHong Zhang       for (i=0; i<m; i++)
1963d06b337dSHong Zhang         t[i] = b[i];
1964d06b337dSHong Zhang 
1965d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
1966d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1967d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1968d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1969d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
197044706875SHong Zhang         PetscLogFlops(2*nz-1);
1971d06b337dSHong Zhang       }
1972d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
1973d06b337dSHong Zhang         d  = *(aa + ai[i]);
1974d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1975d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1976d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1977d06b337dSHong Zhang         sum = t[i];
1978d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
197944706875SHong Zhang         PetscLogFlops(2*nz-1);
1980d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
1981d06b337dSHong Zhang       }
1982d06b337dSHong Zhang     }
1983d06b337dSHong Zhang     its--;
1984d06b337dSHong Zhang   }
1985d06b337dSHong Zhang 
1986d06b337dSHong Zhang   while (its--) {
198744706875SHong Zhang     /*
198844706875SHong Zhang        forward sweep:
198944706875SHong Zhang        for i=0,...,m-1:
199044706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
199144706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
199244706875SHong Zhang          b      = b - x[i]*U^T(i,:);
1993d06b337dSHong Zhang 
199444706875SHong Zhang     */
19952798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1996d06b337dSHong Zhang       for (i=0; i<m; i++)
1997d06b337dSHong Zhang         t[i] = b[i];
1998d06b337dSHong Zhang 
1999d06b337dSHong Zhang       for (i=0; i<m; i++){
200044706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2001d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2002d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2003d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2004d06b337dSHong Zhang         sum = t[i];
2005d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2006d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2007d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
200844706875SHong Zhang         PetscLogFlops(4*nz-2);
2009d06b337dSHong Zhang       }
2010d06b337dSHong Zhang     }
2011d06b337dSHong Zhang 
20122798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
201344706875SHong Zhang       /*
201444706875SHong Zhang        backward sweep:
201544706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
201644706875SHong Zhang        for i=m-1,...,0:
201744706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
201844706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
201944706875SHong Zhang       */
2020d06b337dSHong Zhang       for (i=0; i<m; i++)
2021d06b337dSHong Zhang         t[i] = b[i];
2022d06b337dSHong Zhang 
2023d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2024d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2025d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2026d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2027d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
202844706875SHong Zhang         PetscLogFlops(2*nz-1);
2029d06b337dSHong Zhang       }
2030d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2031d06b337dSHong Zhang         d  = *(aa + ai[i]);
2032d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2033d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2034d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2035d06b337dSHong Zhang         sum = t[i];
2036d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
203744706875SHong Zhang         PetscLogFlops(2*nz-1);
2038d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2039d06b337dSHong Zhang       }
2040d06b337dSHong Zhang     }
2041d06b337dSHong Zhang   }
2042d06b337dSHong Zhang 
2043d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
2044d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2045d06b337dSHong Zhang   if (bb != xx) {
2046d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2047d06b337dSHong Zhang   }
20482798e883SHong Zhang 
2049d06b337dSHong Zhang   PetscFunctionReturn(0);
2050d06b337dSHong Zhang }
2051d06b337dSHong Zhang 
2052d06b337dSHong Zhang 
2053d06b337dSHong Zhang 
2054d06b337dSHong Zhang 
205549b5e25fSSatish Balay 
205649b5e25fSSatish Balay 
2057