xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision a23d5ecec69aaee6ed47a85ba9d72960301bf762)
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);}
124d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
125d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
126d59c15a7SBarry Smith   if (a->mult_work)   {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
12749b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
128c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
1291a3463dfSHong Zhang 
1301a3463dfSHong Zhang   if (a->inew){
1311a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1321a3463dfSHong Zhang     a->inew = 0;
1331a3463dfSHong Zhang   }
13449b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13549b5e25fSSatish Balay   PetscFunctionReturn(0);
13649b5e25fSSatish Balay }
13749b5e25fSSatish Balay 
1384a2ae208SSatish Balay #undef __FUNCT__
1394a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
14049b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14149b5e25fSSatish Balay {
142045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14349b5e25fSSatish Balay 
14449b5e25fSSatish Balay   PetscFunctionBegin;
1454d9d31abSKris Buschelman   switch (op) {
1464d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1474d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1484d9d31abSKris Buschelman     break;
1494d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1504d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1514d9d31abSKris Buschelman     break;
1524d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1534d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1564d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1594d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1624d9d31abSKris Buschelman     a->nonew          = 1;
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1654d9d31abSKris Buschelman     a->nonew          = -1;
1664d9d31abSKris Buschelman     break;
1674d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1684d9d31abSKris Buschelman     a->nonew          = -2;
1694d9d31abSKris Buschelman     break;
1704d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1714d9d31abSKris Buschelman     a->nonew          = 0;
1724d9d31abSKris Buschelman     break;
1734d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1744d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1754d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1764d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1774d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
178b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1794d9d31abSKris Buschelman     break;
1804d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
18129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1824d9d31abSKris Buschelman   default:
18329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
18449b5e25fSSatish Balay   }
18549b5e25fSSatish Balay   PetscFunctionReturn(0);
18649b5e25fSSatish Balay }
18749b5e25fSSatish Balay 
1884a2ae208SSatish Balay #undef __FUNCT__
1894a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19087828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
19149b5e25fSSatish Balay {
19249b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
19382502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
19449b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
19587828ca2SBarry Smith   PetscScalar  *v_i;
19649b5e25fSSatish Balay 
19749b5e25fSSatish Balay   PetscFunctionBegin;
19849b5e25fSSatish Balay   bs  = a->bs;
19949b5e25fSSatish Balay   ai  = a->i;
20049b5e25fSSatish Balay   aj  = a->j;
20149b5e25fSSatish Balay   aa  = a->a;
20249b5e25fSSatish Balay   bs2 = a->bs2;
20349b5e25fSSatish Balay 
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;
75527f1aa60SHong Zhang #if defined(PETSC_HAVE_SPOOLES) || defined(PETSC_HAVE_MUMPS)
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
80327f1aa60SHong Zhang #if defined(PETSC_HAVE_MUMPS)
80427f1aa60SHong Zhang   ierr = PetscOptionsHasName(A->prefix,"-mat_sbaij_mumps",&flag);CHKERRQ(ierr);
80527f1aa60SHong Zhang   if (flag) { ierr = MatUseMUMPS_MPIAIJ(A);CHKERRQ(ierr); }
80627f1aa60SHong Zhang #endif
807cf4441caSHong Zhang 
80849b5e25fSSatish Balay   PetscFunctionReturn(0);
80949b5e25fSSatish Balay }
81049b5e25fSSatish Balay 
81149b5e25fSSatish Balay /*
81249b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
81349b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
81449b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
81549b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
81649b5e25fSSatish Balay */
8174a2ae208SSatish Balay #undef __FUNCT__
8184a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
819db2f66daSBarry Smith int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
82049b5e25fSSatish Balay {
82149b5e25fSSatish Balay   int        i,j,k,row;
82249b5e25fSSatish Balay   PetscTruth flg;
82349b5e25fSSatish Balay 
82449b5e25fSSatish Balay   PetscFunctionBegin;
82549b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
82649b5e25fSSatish Balay     row = idx[i];
82749b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
82849b5e25fSSatish Balay       sizes[j] = 1;
82949b5e25fSSatish Balay       i++;
83049b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
83149b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
83249b5e25fSSatish Balay       i++;
83349b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
83449b5e25fSSatish Balay       flg = PETSC_TRUE;
83549b5e25fSSatish Balay       for (k=1; k<bs; k++) {
83649b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
83749b5e25fSSatish Balay           flg = PETSC_FALSE;
83849b5e25fSSatish Balay           break;
83949b5e25fSSatish Balay         }
84049b5e25fSSatish Balay       }
84149b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
84249b5e25fSSatish Balay         sizes[j] = bs;
84349b5e25fSSatish Balay         i+= bs;
84449b5e25fSSatish Balay       } else {
84549b5e25fSSatish Balay         sizes[j] = 1;
84649b5e25fSSatish Balay         i++;
84749b5e25fSSatish Balay       }
84849b5e25fSSatish Balay     }
84949b5e25fSSatish Balay   }
85049b5e25fSSatish Balay   *bs_max = j;
85149b5e25fSSatish Balay   PetscFunctionReturn(0);
85249b5e25fSSatish Balay }
85349b5e25fSSatish Balay 
8544a2ae208SSatish Balay #undef __FUNCT__
8554a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
85687828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag)
85749b5e25fSSatish Balay {
85849b5e25fSSatish Balay   PetscFunctionBegin;
859c0f24835SHong Zhang   SETERRQ(PETSC_ERR_SUP,"No support for this function yet");
86049b5e25fSSatish Balay }
86149b5e25fSSatish Balay 
86249b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
86349b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
86449b5e25fSSatish Balay */
86549b5e25fSSatish Balay 
8664a2ae208SSatish Balay #undef __FUNCT__
8674a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
86887828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
86949b5e25fSSatish Balay {
87049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
87149b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
87249b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
87349b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
87449b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
87549b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
87649b5e25fSSatish Balay 
87749b5e25fSSatish Balay   PetscFunctionBegin;
87849b5e25fSSatish Balay 
87949b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
88049b5e25fSSatish Balay     row  = im[k];       /* row number */
88149b5e25fSSatish Balay     brow = row/bs;      /* block row number */
88249b5e25fSSatish Balay     if (row < 0) continue;
88349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
884c464158bSHong Zhang     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
88549b5e25fSSatish Balay #endif
88649b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
88749b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
88849b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
88949b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
89049b5e25fSSatish Balay     low  = 0;
89149b5e25fSSatish Balay 
89249b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
89349b5e25fSSatish Balay       if (in[l] < 0) continue;
89449b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
895c464158bSHong Zhang       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m);
89649b5e25fSSatish Balay #endif
89749b5e25fSSatish Balay       col = in[l];
89849b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
89949b5e25fSSatish Balay 
90049b5e25fSSatish Balay       if (brow <= bcol){
90149b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
9028549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
90349b5e25fSSatish Balay           /* element value a(k,l) */
90449b5e25fSSatish Balay           if (roworiented) {
90549b5e25fSSatish Balay             value = v[l + k*n];
90649b5e25fSSatish Balay           } else {
90749b5e25fSSatish Balay             value = v[k + l*m];
90849b5e25fSSatish Balay           }
90949b5e25fSSatish Balay 
91049b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
91149b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
91249b5e25fSSatish Balay           while (high-low > 7) {
91349b5e25fSSatish Balay             t = (low+high)/2;
91449b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
91549b5e25fSSatish Balay             else              low  = t;
91649b5e25fSSatish Balay           }
91749b5e25fSSatish Balay           for (i=low; i<high; i++) {
91849b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
91949b5e25fSSatish Balay             if (rp[i] > bcol) break;
92049b5e25fSSatish Balay             if (rp[i] == bcol) {
92149b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
92249b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
92349b5e25fSSatish Balay               else                  *bap  = value;
9248549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9258549e402SHong Zhang               if (brow == bcol && ridx < cidx){
9268549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
9278549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
9288549e402SHong Zhang                 else                  *bap  = value;
9298549e402SHong Zhang               }
93049b5e25fSSatish Balay               goto noinsert1;
93149b5e25fSSatish Balay             }
93249b5e25fSSatish Balay           }
93349b5e25fSSatish Balay 
93449b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
93529bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
93649b5e25fSSatish Balay       if (nrow >= rmax) {
93749b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
93849b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
93949b5e25fSSatish Balay         MatScalar *new_a;
94049b5e25fSSatish Balay 
94129bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
94249b5e25fSSatish Balay 
94349b5e25fSSatish Balay         /* Malloc new storage space */
94449b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
94582502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
94649b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
94749b5e25fSSatish Balay         new_i = new_j + new_nz;
94849b5e25fSSatish Balay 
94949b5e25fSSatish Balay         /* copy over old data into new slots */
95049b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
95149b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
95249b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
95349b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
95449b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
95549b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
95649b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
95749b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
95849b5e25fSSatish Balay         /* free up old matrix storage */
95949b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
96049b5e25fSSatish Balay         if (!a->singlemalloc) {
96149b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
96249b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
96349b5e25fSSatish Balay         }
96449b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
96549b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
96649b5e25fSSatish Balay 
96749b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
96849b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
969b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
97049b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
97149b5e25fSSatish Balay         a->reallocs++;
97249b5e25fSSatish Balay         a->s_nz++;
97349b5e25fSSatish Balay       }
97449b5e25fSSatish Balay 
97549b5e25fSSatish Balay       N = nrow++ - 1;
97649b5e25fSSatish Balay       /* shift up all the later entries in this row */
97749b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
97849b5e25fSSatish Balay         rp[ii+1] = rp[ii];
97949b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
98049b5e25fSSatish Balay       }
98149b5e25fSSatish Balay       if (N>=i) {
98249b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
98349b5e25fSSatish Balay       }
98449b5e25fSSatish Balay       rp[i]                      = bcol;
98549b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
98649b5e25fSSatish Balay       noinsert1:;
98749b5e25fSSatish Balay       low = i;
98849b5e25fSSatish Balay       /* } */
9898549e402SHong Zhang         }
99049b5e25fSSatish Balay       } /* end of if .. if..  */
99149b5e25fSSatish Balay     }                     /* end of loop over added columns */
99249b5e25fSSatish Balay     ailen[brow] = nrow;
99349b5e25fSSatish Balay   }                       /* end of loop over added rows */
99449b5e25fSSatish Balay 
99549b5e25fSSatish Balay   PetscFunctionReturn(0);
99649b5e25fSSatish Balay }
99749b5e25fSSatish Balay 
99815e8a5b3SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
99915e8a5b3SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
100049b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
100149b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
100249b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
100349b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
100449b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
100587828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat);
100649b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
100749b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
100849b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
100949b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
101049b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
101149b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
101213a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
101349b5e25fSSatish Balay 
101449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
101549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
101649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
101749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
101849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
101949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
102049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
102149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
102249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
102349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
102449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
102549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
102649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
102749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
102849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
102949b5e25fSSatish Balay 
1030d59c15a7SBarry Smith extern int MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
1031d59c15a7SBarry Smith 
103207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
103307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
103407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
103507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
103607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
103707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
103807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
103907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
104007f98182SSatish Balay 
10415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
10425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
10435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
10445f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
10455f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
10465f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
10475f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
104857d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
10493e0d88b5SBarry Smith extern int MatGetInertia_SeqSBAIJ(Mat,int*,int*,int*);
105049b5e25fSSatish Balay 
105149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
105249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
105349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
105449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
105549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
105649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
105749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
105849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
105949b5e25fSSatish Balay 
106049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
106149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
106249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
106349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
106449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
106549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
106649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
106749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
106849b5e25fSSatish Balay 
10694d101231SSatish Balay #ifdef HAVE_MatICCFactor
10701a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
10714a2ae208SSatish Balay #undef __FUNCT__
10724d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
107315e8a5b3SHong Zhang int MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
107449b5e25fSSatish Balay {
10754ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
107649b5e25fSSatish Balay   Mat         outA;
107749b5e25fSSatish Balay   int         ierr;
107849b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
107949b5e25fSSatish Balay 
108049b5e25fSSatish Balay   PetscFunctionBegin;
10811a3463dfSHong Zhang   /*
108229bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
108349b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
108449b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
108549b5e25fSSatish Balay   if (!row_identity || !col_identity) {
108629bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
108749b5e25fSSatish Balay   }
10881a3463dfSHong Zhang   */
108949b5e25fSSatish Balay 
109049b5e25fSSatish Balay   outA          = inA;
10911260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
109249b5e25fSSatish Balay 
109349b5e25fSSatish Balay   if (!a->diag) {
10941a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
109549b5e25fSSatish Balay   }
109649b5e25fSSatish Balay   /*
109749b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
109849b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
109949b5e25fSSatish Balay   */
110049b5e25fSSatish Balay   switch (a->bs) {
110149b5e25fSSatish Balay   case 1:
11020fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
110307f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1104d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
11054d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
110649b5e25fSSatish Balay   case 2:
11071a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
11081a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11091a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11104d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
111149b5e25fSSatish Balay     break;
111249b5e25fSSatish Balay   case 3:
11131a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
11141a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11151a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11164d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
111749b5e25fSSatish Balay     break;
111849b5e25fSSatish Balay   case 4:
11191a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
11201a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11211a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11224d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
112349b5e25fSSatish Balay     break;
112449b5e25fSSatish Balay   case 5:
11251a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
11261a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11271a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11284d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
112949b5e25fSSatish Balay     break;
113049b5e25fSSatish Balay   case 6:
11311a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
11321a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11331a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11344d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
113549b5e25fSSatish Balay     break;
113649b5e25fSSatish Balay   case 7:
11371a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
11381a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11391a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11404d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
114149b5e25fSSatish Balay     break;
114249b5e25fSSatish Balay   default:
114349b5e25fSSatish Balay     a->row        = row;
11441a3463dfSHong Zhang     a->icol       = col;
114549b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
114649b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
114749b5e25fSSatish Balay 
114849b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
114949b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1150b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
115149b5e25fSSatish Balay 
115249b5e25fSSatish Balay     if (!a->solve_work) {
115387828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
115487828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
115549b5e25fSSatish Balay     }
115649b5e25fSSatish Balay   }
115749b5e25fSSatish Balay 
11581a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
115949b5e25fSSatish Balay 
116049b5e25fSSatish Balay   PetscFunctionReturn(0);
116149b5e25fSSatish Balay }
1162950f1e5bSHong Zhang #endif
1163950f1e5bSHong Zhang 
11644a2ae208SSatish Balay #undef __FUNCT__
11654a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
116649b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
116749b5e25fSSatish Balay {
116849b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
116949b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
117049b5e25fSSatish Balay   int               ierr;
117149b5e25fSSatish Balay 
117249b5e25fSSatish Balay   PetscFunctionBegin;
117349b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
117449b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
117549b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
117649b5e25fSSatish Balay   PetscFunctionReturn(0);
117749b5e25fSSatish Balay }
117849b5e25fSSatish Balay 
117949b5e25fSSatish Balay EXTERN_C_BEGIN
11804a2ae208SSatish Balay #undef __FUNCT__
11814a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
118249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
118349b5e25fSSatish Balay {
1184045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
118549b5e25fSSatish Balay   int         i,nz,n;
118649b5e25fSSatish Balay 
118749b5e25fSSatish Balay   PetscFunctionBegin;
1188045c9aa0SHong Zhang   nz = baij->s_maxnz;
1189c464158bSHong Zhang   n  = mat->n;
119049b5e25fSSatish Balay   for (i=0; i<nz; i++) {
119149b5e25fSSatish Balay     baij->j[i] = indices[i];
119249b5e25fSSatish Balay   }
1193045c9aa0SHong Zhang   baij->s_nz = nz;
119449b5e25fSSatish Balay   for (i=0; i<n; i++) {
119549b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
119649b5e25fSSatish Balay   }
119749b5e25fSSatish Balay 
119849b5e25fSSatish Balay   PetscFunctionReturn(0);
119949b5e25fSSatish Balay }
120049b5e25fSSatish Balay EXTERN_C_END
120149b5e25fSSatish Balay 
12024a2ae208SSatish Balay #undef __FUNCT__
12034a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
120449b5e25fSSatish Balay /*@
120519585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
120649b5e25fSSatish Balay        in the matrix.
120749b5e25fSSatish Balay 
120849b5e25fSSatish Balay   Input Parameters:
120919585528SSatish Balay +  mat     - the SeqSBAIJ matrix
121049b5e25fSSatish Balay -  indices - the column indices
121149b5e25fSSatish Balay 
121249b5e25fSSatish Balay   Level: advanced
121349b5e25fSSatish Balay 
121449b5e25fSSatish Balay   Notes:
121549b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
121649b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
121749b5e25fSSatish Balay   of the MatSetValues() operation.
121849b5e25fSSatish Balay 
121949b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
122019585528SSatish Balay   MatCreateSeqSBAIJ().
122149b5e25fSSatish Balay 
1222ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
122349b5e25fSSatish Balay 
1224ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
122549b5e25fSSatish Balay @*/
122649b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
122749b5e25fSSatish Balay {
122849b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
122949b5e25fSSatish Balay 
123049b5e25fSSatish Balay   PetscFunctionBegin;
123149b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1232c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
123349b5e25fSSatish Balay   if (f) {
123449b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
123549b5e25fSSatish Balay   } else {
123629bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
123749b5e25fSSatish Balay   }
123849b5e25fSSatish Balay   PetscFunctionReturn(0);
123949b5e25fSSatish Balay }
124049b5e25fSSatish Balay 
12414a2ae208SSatish Balay #undef __FUNCT__
12424a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1243273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1244273d9f13SBarry Smith {
1245273d9f13SBarry Smith   int        ierr;
1246273d9f13SBarry Smith 
1247273d9f13SBarry Smith   PetscFunctionBegin;
1248273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1249273d9f13SBarry Smith   PetscFunctionReturn(0);
1250273d9f13SBarry Smith }
1251273d9f13SBarry Smith 
1252a6ece127SHong Zhang #undef __FUNCT__
1253a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1254a6ece127SHong Zhang int MatGetArray_SeqSBAIJ(Mat A,PetscScalar **array)
1255a6ece127SHong Zhang {
1256a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1257a6ece127SHong Zhang   PetscFunctionBegin;
1258a6ece127SHong Zhang   *array = a->a;
1259a6ece127SHong Zhang   PetscFunctionReturn(0);
1260a6ece127SHong Zhang }
1261a6ece127SHong Zhang 
1262a6ece127SHong Zhang #undef __FUNCT__
1263a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1264a6ece127SHong Zhang int MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar **array)
1265a6ece127SHong Zhang {
1266a6ece127SHong Zhang   PetscFunctionBegin;
1267a6ece127SHong Zhang   PetscFunctionReturn(0);
1268a6ece127SHong Zhang }
1269a6ece127SHong Zhang 
127042ee4b1aSHong Zhang #include "petscblaslapack.h"
127142ee4b1aSHong Zhang #undef __FUNCT__
127242ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
127342ee4b1aSHong Zhang int MatAXPY_SeqSBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
127442ee4b1aSHong Zhang {
127542ee4b1aSHong Zhang   Mat_SeqSBAIJ *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1276c4319e64SHong Zhang   int          ierr,one=1,i,bs=y->bs,bs2,j;
127742ee4b1aSHong Zhang 
127842ee4b1aSHong Zhang   PetscFunctionBegin;
127942ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
128042ee4b1aSHong Zhang     BLaxpy_(&x->s_nz,a,x->a,&one,y->a,&one);
1281c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1282c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1283c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1284c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1285c537a176SHong Zhang     }
1286c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1287c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1288c4319e64SHong Zhang       y->XtoY = X;
1289c537a176SHong Zhang     }
1290c4319e64SHong Zhang     bs2 = bs*bs;
1291c537a176SHong Zhang     for (i=0; i<x->s_nz; i++) {
1292c4319e64SHong Zhang       j = 0;
1293c4319e64SHong Zhang       while (j < bs2){
12946550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1295c4319e64SHong Zhang         j++;
1296c537a176SHong Zhang       }
1297c4319e64SHong Zhang     }
1298c4319e64SHong Zhang     PetscLogInfo(0,"MatAXPY_SeqSBAIJ: ratio of nnz_s(X)/nnz_s(Y): %d/%d = %g\n",bs2*x->s_nz,bs2*y->s_nz,(PetscReal)(bs2*x->s_nz)/(bs2*y->s_nz));
129942ee4b1aSHong Zhang   } else {
130042ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
130142ee4b1aSHong Zhang   }
130242ee4b1aSHong Zhang   PetscFunctionReturn(0);
130342ee4b1aSHong Zhang }
130442ee4b1aSHong Zhang 
130549b5e25fSSatish Balay /* -------------------------------------------------------------------*/
130649b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
130749b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
130849b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
130949b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
131049b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
131149b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
131249b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
131349b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
131449b5e25fSSatish Balay        0,
131549b5e25fSSatish Balay        0,
131649b5e25fSSatish Balay        0,
131749b5e25fSSatish Balay        0,
13185f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1319d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
132049b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
132149b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
132249b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
132349b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
132449b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
132549b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
132649b5e25fSSatish Balay        0,
132749b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
132849b5e25fSSatish Balay        0,
132949b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
133049b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
133149b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
133249b5e25fSSatish Balay        0,
133349b5e25fSSatish Balay        0,
13345f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13355f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
1336273d9f13SBarry Smith        MatSetUpPreallocation_SeqSBAIJ,
1337c464158bSHong Zhang        0,
13384d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1339a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1340a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
134149b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
134249b5e25fSSatish Balay        0,
134349b5e25fSSatish Balay        0,
134449b5e25fSSatish Balay        0,
1345950f1e5bSHong Zhang        0,
134642ee4b1aSHong Zhang        MatAXPY_SeqSBAIJ,
134749b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
134849b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
134949b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
135049b5e25fSSatish Balay        0,
135149b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
135249b5e25fSSatish Balay        MatScale_SeqSBAIJ,
135349b5e25fSSatish Balay        0,
135449b5e25fSSatish Balay        0,
135549b5e25fSSatish Balay        0,
135649b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
135749b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
135849b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
135949b5e25fSSatish Balay        0,
136049b5e25fSSatish Balay        0,
136149b5e25fSSatish Balay        0,
136249b5e25fSSatish Balay        0,
136349b5e25fSSatish Balay        0,
136449b5e25fSSatish Balay        0,
136549b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
136649b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
136749b5e25fSSatish Balay        0,
136849b5e25fSSatish Balay        0,
13698a124369SBarry Smith        MatGetPetscMaps_Petsc,
1370d959ec07SHong Zhang        0,
1371d959ec07SHong Zhang        0,
1372d959ec07SHong Zhang        0,
1373d959ec07SHong Zhang        0,
1374d959ec07SHong Zhang        0,
1375d959ec07SHong Zhang        0,
13763e0d88b5SBarry Smith        MatGetRowMax_SeqSBAIJ,
13773e0d88b5SBarry Smith        0,
13783e0d88b5SBarry Smith        0,
13793e0d88b5SBarry Smith        0,
13803e0d88b5SBarry Smith        0,
13813e0d88b5SBarry Smith        0,
13823e0d88b5SBarry Smith        0,
13833e0d88b5SBarry Smith        0,
13843e0d88b5SBarry Smith        0,
13853e0d88b5SBarry Smith        0,
13863e0d88b5SBarry Smith        0,
13873e0d88b5SBarry Smith        0,
13883e0d88b5SBarry Smith        0,
13893e0d88b5SBarry Smith        0,
13903e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
13913e0d88b5SBarry Smith        MatGetInertia_SeqSBAIJ
13923e0d88b5SBarry Smith #else
13933e0d88b5SBarry Smith        0
13943e0d88b5SBarry Smith #endif
13953e0d88b5SBarry Smith };
139649b5e25fSSatish Balay 
139749b5e25fSSatish Balay EXTERN_C_BEGIN
13984a2ae208SSatish Balay #undef __FUNCT__
13994a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
140049b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
140149b5e25fSSatish Balay {
14024afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1403c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
140449b5e25fSSatish Balay   int         ierr;
140549b5e25fSSatish Balay 
140649b5e25fSSatish Balay   PetscFunctionBegin;
140749b5e25fSSatish Balay   if (aij->nonew != 1) {
140829bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
140949b5e25fSSatish Balay   }
141049b5e25fSSatish Balay 
141149b5e25fSSatish Balay   /* allocate space for values if not already there */
141249b5e25fSSatish Balay   if (!aij->saved_values) {
141387828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
141449b5e25fSSatish Balay   }
141549b5e25fSSatish Balay 
141649b5e25fSSatish Balay   /* copy values over */
141787828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
141849b5e25fSSatish Balay   PetscFunctionReturn(0);
141949b5e25fSSatish Balay }
142049b5e25fSSatish Balay EXTERN_C_END
142149b5e25fSSatish Balay 
142249b5e25fSSatish Balay EXTERN_C_BEGIN
14234a2ae208SSatish Balay #undef __FUNCT__
14244a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
142549b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
142649b5e25fSSatish Balay {
14274afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1428c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
142949b5e25fSSatish Balay 
143049b5e25fSSatish Balay   PetscFunctionBegin;
143149b5e25fSSatish Balay   if (aij->nonew != 1) {
143229bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
143349b5e25fSSatish Balay   }
143449b5e25fSSatish Balay   if (!aij->saved_values) {
143529bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
143649b5e25fSSatish Balay   }
143749b5e25fSSatish Balay 
143849b5e25fSSatish Balay   /* copy values over */
143987828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
144049b5e25fSSatish Balay   PetscFunctionReturn(0);
144149b5e25fSSatish Balay }
144249b5e25fSSatish Balay EXTERN_C_END
144349b5e25fSSatish Balay 
14448549e402SHong Zhang EXTERN_C_BEGIN
14454a2ae208SSatish Balay #undef __FUNCT__
1446*a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1447*a23d5eceSKris Buschelman int MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,int bs,int nz,int *nnz)
144849b5e25fSSatish Balay {
1449c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
14500c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
145149b5e25fSSatish Balay   PetscTruth   flg;
14524afc71dfSHong Zhang   int          s_nz;
145349b5e25fSSatish Balay 
145449b5e25fSSatish Balay   PetscFunctionBegin;
1455273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1456e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1457c464158bSHong Zhang   mbs  = B->m/bs;
145849b5e25fSSatish Balay   bs2  = bs*bs;
145949b5e25fSSatish Balay 
1460c464158bSHong Zhang   if (mbs*bs != B->m) {
146129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
146249b5e25fSSatish Balay   }
146349b5e25fSSatish Balay 
1464435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1465435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
146649b5e25fSSatish Balay   if (nnz) {
146749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
146829bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
146980fe4e49SBarry 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);
147049b5e25fSSatish Balay     }
147149b5e25fSSatish Balay   }
147249b5e25fSSatish Balay 
1473e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
147449b5e25fSSatish Balay   if (!flg) {
147549b5e25fSSatish Balay     switch (bs) {
147649b5e25fSSatish Balay     case 1:
14774ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
147849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1479d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
148049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
148149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
148249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
148349b5e25fSSatish Balay       break;
148449b5e25fSSatish Balay     case 2:
14854ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
148649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
148749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
148849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
148949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
149049b5e25fSSatish Balay       break;
149149b5e25fSSatish Balay     case 3:
14925f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
149349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
149449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
149549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
149649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
149749b5e25fSSatish Balay       break;
149849b5e25fSSatish Balay     case 4:
14995f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
150049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
150149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
150249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
150349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
150449b5e25fSSatish Balay       break;
150549b5e25fSSatish Balay     case 5:
15065f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
150749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
150849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
150949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
151049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
151149b5e25fSSatish Balay       break;
151249b5e25fSSatish Balay     case 6:
15135f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
151449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
151549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
151649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
151749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
151849b5e25fSSatish Balay       break;
151949b5e25fSSatish Balay     case 7:
1520de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
152149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
152249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1523de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
152449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
152549b5e25fSSatish Balay       break;
152649b5e25fSSatish Balay     }
152749b5e25fSSatish Balay   }
152849b5e25fSSatish Balay 
152949b5e25fSSatish Balay   b->mbs = mbs;
15304afc71dfSHong Zhang   b->nbs = mbs;
1531b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
153249b5e25fSSatish Balay   if (!nnz) {
1533435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
153449b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
153549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
15368cef66ccSHong Zhang       b->imax[i] = nz;
153749b5e25fSSatish Balay     }
1538153ea458SHong Zhang     nz = nz*mbs; /* total nz */
153949b5e25fSSatish Balay   } else {
154049b5e25fSSatish Balay     nz = 0;
15418cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
154249b5e25fSSatish Balay   }
15438cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
15448cef66ccSHong Zhang   s_nz = nz;
154549b5e25fSSatish Balay 
154649b5e25fSSatish Balay   /* allocate the matrix space */
1547c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
154882502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
154949b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
155049b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
155149b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
155249b5e25fSSatish Balay   b->i = b->j + s_nz;
155349b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
155449b5e25fSSatish Balay 
155549b5e25fSSatish Balay   /* pointer to beginning of each row */
155649b5e25fSSatish Balay   b->i[0] = 0;
155749b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
155849b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
155949b5e25fSSatish Balay   }
156049b5e25fSSatish Balay 
156149b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
15625033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1563b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
156449b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
156549b5e25fSSatish Balay 
156649b5e25fSSatish Balay   b->bs               = bs;
156749b5e25fSSatish Balay   b->bs2              = bs2;
156849b5e25fSSatish Balay   b->s_nz             = 0;
156949b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1570153ea458SHong Zhang 
157116cdd363SHong Zhang   b->inew             = 0;
157216cdd363SHong Zhang   b->jnew             = 0;
157316cdd363SHong Zhang   b->anew             = 0;
157416cdd363SHong Zhang   b->a2anew           = 0;
15751a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1576c464158bSHong Zhang   PetscFunctionReturn(0);
1577c464158bSHong Zhang }
1578*a23d5eceSKris Buschelman EXTERN_C_END
1579153ea458SHong Zhang 
1580*a23d5eceSKris Buschelman EXTERN_C_BEGIN
1581*a23d5eceSKris Buschelman #undef __FUNCT__
1582*a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1583*a23d5eceSKris Buschelman int MatCreate_SeqSBAIJ(Mat B)
1584*a23d5eceSKris Buschelman {
1585*a23d5eceSKris Buschelman   Mat_SeqSBAIJ *b;
1586*a23d5eceSKris Buschelman   int          ierr,size;
1587*a23d5eceSKris Buschelman 
1588*a23d5eceSKris Buschelman   PetscFunctionBegin;
1589*a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1590*a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1591*a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1592*a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1593*a23d5eceSKris Buschelman 
1594*a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1595*a23d5eceSKris Buschelman   B->data = (void*)b;
1596*a23d5eceSKris Buschelman   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1597*a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1598*a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1599*a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1600*a23d5eceSKris Buschelman   B->factor           = 0;
1601*a23d5eceSKris Buschelman   B->lupivotthreshold = 1.0;
1602*a23d5eceSKris Buschelman   B->mapping          = 0;
1603*a23d5eceSKris Buschelman   b->row              = 0;
1604*a23d5eceSKris Buschelman   b->icol             = 0;
1605*a23d5eceSKris Buschelman   b->reallocs         = 0;
1606*a23d5eceSKris Buschelman   b->saved_values     = 0;
1607*a23d5eceSKris Buschelman 
1608*a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1609*a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1610*a23d5eceSKris Buschelman 
1611*a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1612*a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1613*a23d5eceSKris Buschelman   b->nonew            = 0;
1614*a23d5eceSKris Buschelman   b->diag             = 0;
1615*a23d5eceSKris Buschelman   b->solve_work       = 0;
1616*a23d5eceSKris Buschelman   b->mult_work        = 0;
1617*a23d5eceSKris Buschelman   B->spptr            = 0;
1618*a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1619*a23d5eceSKris Buschelman   b->xtoy             = 0;
1620*a23d5eceSKris Buschelman   b->XtoY             = 0;
1621*a23d5eceSKris Buschelman 
1622*a23d5eceSKris Buschelman   b->inew             = 0;
1623*a23d5eceSKris Buschelman   b->jnew             = 0;
1624*a23d5eceSKris Buschelman   b->anew             = 0;
1625*a23d5eceSKris Buschelman   b->a2anew           = 0;
1626*a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1627*a23d5eceSKris Buschelman 
1628*a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1629*a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1630*a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1631*a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1632*a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1633*a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1634*a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1635*a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1636*a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1637*a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1638*a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1639*a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
1640*a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1641*a23d5eceSKris Buschelman }
1642*a23d5eceSKris Buschelman EXTERN_C_END
1643*a23d5eceSKris Buschelman 
1644*a23d5eceSKris Buschelman #undef __FUNCT__
1645*a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1646*a23d5eceSKris Buschelman /*@C
1647*a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1648*a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1649*a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1650*a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1651*a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1652*a23d5eceSKris Buschelman 
1653*a23d5eceSKris Buschelman    Collective on Mat
1654*a23d5eceSKris Buschelman 
1655*a23d5eceSKris Buschelman    Input Parameters:
1656*a23d5eceSKris Buschelman +  A - the symmetric matrix
1657*a23d5eceSKris Buschelman .  bs - size of block
1658*a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1659*a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1660*a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1661*a23d5eceSKris Buschelman 
1662*a23d5eceSKris Buschelman    Options Database Keys:
1663*a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1664*a23d5eceSKris Buschelman                      block calculations (much slower)
1665*a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1666*a23d5eceSKris Buschelman 
1667*a23d5eceSKris Buschelman    Level: intermediate
1668*a23d5eceSKris Buschelman 
1669*a23d5eceSKris Buschelman    Notes:
1670*a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1671*a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1672*a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1673*a23d5eceSKris Buschelman    matrices.
1674*a23d5eceSKris Buschelman 
1675*a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1676*a23d5eceSKris Buschelman @*/
1677*a23d5eceSKris Buschelman int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz) {
1678*a23d5eceSKris Buschelman   int ierr,(*f)(Mat,int,int,int *);
1679*a23d5eceSKris Buschelman 
1680*a23d5eceSKris Buschelman   PetscFunctionBegin;
1681*a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1682*a23d5eceSKris Buschelman   if (f) {
1683*a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1684*a23d5eceSKris Buschelman   }
1685*a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1686*a23d5eceSKris Buschelman }
168749b5e25fSSatish Balay 
16884a2ae208SSatish Balay #undef __FUNCT__
16894a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1690c464158bSHong Zhang /*@C
1691c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1692c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1693c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1694c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1695c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
169649b5e25fSSatish Balay 
1697c464158bSHong Zhang    Collective on MPI_Comm
1698c464158bSHong Zhang 
1699c464158bSHong Zhang    Input Parameters:
1700c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1701c464158bSHong Zhang .  bs - size of block
1702c464158bSHong Zhang .  m - number of rows, or number of columns
1703c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1704744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1705744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1706c464158bSHong Zhang 
1707c464158bSHong Zhang    Output Parameter:
1708c464158bSHong Zhang .  A - the symmetric matrix
1709c464158bSHong Zhang 
1710c464158bSHong Zhang    Options Database Keys:
1711c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1712c464158bSHong Zhang                      block calculations (much slower)
1713c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1714c464158bSHong Zhang 
1715c464158bSHong Zhang    Level: intermediate
1716c464158bSHong Zhang 
1717c464158bSHong Zhang    Notes:
1718c464158bSHong Zhang 
1719c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1720c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1721c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1722c464158bSHong Zhang    matrices.
1723c464158bSHong Zhang 
1724c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1725c464158bSHong Zhang @*/
1726c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1727c464158bSHong Zhang {
1728c464158bSHong Zhang   int ierr;
1729c464158bSHong Zhang 
1730c464158bSHong Zhang   PetscFunctionBegin;
1731c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1732c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1733273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
173449b5e25fSSatish Balay   PetscFunctionReturn(0);
173549b5e25fSSatish Balay }
173649b5e25fSSatish Balay 
17374a2ae208SSatish Balay #undef __FUNCT__
17384a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
173949b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
174049b5e25fSSatish Balay {
174149b5e25fSSatish Balay   Mat         C;
174249b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
174349b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
174449b5e25fSSatish Balay 
174549b5e25fSSatish Balay   PetscFunctionBegin;
174629bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
174749b5e25fSSatish Balay 
174849b5e25fSSatish Balay   *B = 0;
1749692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1750692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1751692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1752692f9cbeSHong Zhang 
175349b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1754273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
175549b5e25fSSatish Balay   C->factor         = A->factor;
175649b5e25fSSatish Balay   c->row            = 0;
175749b5e25fSSatish Balay   c->icol           = 0;
175849b5e25fSSatish Balay   c->saved_values   = 0;
175949b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
176049b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
176149b5e25fSSatish Balay 
176249b5e25fSSatish Balay   c->bs         = a->bs;
176349b5e25fSSatish Balay   c->bs2        = a->bs2;
176449b5e25fSSatish Balay   c->mbs        = a->mbs;
176549b5e25fSSatish Balay   c->nbs        = a->nbs;
176649b5e25fSSatish Balay 
1767b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1768b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
176949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
177049b5e25fSSatish Balay     c->imax[i] = a->imax[i];
177149b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
177249b5e25fSSatish Balay   }
177349b5e25fSSatish Balay 
177449b5e25fSSatish Balay   /* allocate the matrix space */
177549b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
177649b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
177782502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
177849b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
177949b5e25fSSatish Balay   c->i = c->j + nz;
178049b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
178149b5e25fSSatish Balay   if (mbs > 0) {
178249b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
178349b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
178449b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178549b5e25fSSatish Balay     } else {
178649b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178749b5e25fSSatish Balay     }
178849b5e25fSSatish Balay   }
178949b5e25fSSatish Balay 
1790b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
179149b5e25fSSatish Balay   c->sorted      = a->sorted;
179249b5e25fSSatish Balay   c->roworiented = a->roworiented;
179349b5e25fSSatish Balay   c->nonew       = a->nonew;
179449b5e25fSSatish Balay 
179549b5e25fSSatish Balay   if (a->diag) {
1796b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1797b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
179849b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
179949b5e25fSSatish Balay       c->diag[i] = a->diag[i];
180049b5e25fSSatish Balay     }
180149b5e25fSSatish Balay   } else c->diag        = 0;
180249b5e25fSSatish Balay   c->s_nz               = a->s_nz;
180349b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
180449b5e25fSSatish Balay   c->solve_work         = 0;
18052a1b7f2aSHong Zhang   C->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
180649b5e25fSSatish Balay   c->mult_work          = 0;
180749b5e25fSSatish Balay   *B = C;
1808b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
180949b5e25fSSatish Balay   PetscFunctionReturn(0);
181049b5e25fSSatish Balay }
181149b5e25fSSatish Balay 
18128549e402SHong Zhang EXTERN_C_BEGIN
18134a2ae208SSatish Balay #undef __FUNCT__
18144a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1815b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
181649b5e25fSSatish Balay {
181749b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
181849b5e25fSSatish Balay   Mat          B;
181949b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18203f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
182149b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
182249b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
182387828ca2SBarry Smith   PetscScalar  *aa;
182449b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
182527f1aa60SHong Zhang #if defined(PETSC_HAVE_SPOOLES) || defined(PETSC_HAVE_MUMPS)
1826ecc4ab6dSBarry Smith   PetscTruth   flag;
1827ecc4ab6dSBarry Smith #endif
182849b5e25fSSatish Balay 
182949b5e25fSSatish Balay   PetscFunctionBegin;
1830b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
183149b5e25fSSatish Balay   bs2  = bs*bs;
183249b5e25fSSatish Balay 
183349b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
183429bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1835b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
183649b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1837552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
183849b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
183949b5e25fSSatish Balay 
184049b5e25fSSatish Balay   if (header[3] < 0) {
184129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
184249b5e25fSSatish Balay   }
184349b5e25fSSatish Balay 
184429bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
184549b5e25fSSatish Balay 
184649b5e25fSSatish Balay   /*
184749b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
184849b5e25fSSatish Balay     divisible by the blocksize
184949b5e25fSSatish Balay   */
185049b5e25fSSatish Balay   mbs        = M/bs;
185149b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
185249b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
185349b5e25fSSatish Balay   else                  mbs++;
185449b5e25fSSatish Balay   if (extra_rows) {
1855b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
185649b5e25fSSatish Balay   }
185749b5e25fSSatish Balay 
185849b5e25fSSatish Balay   /* read in row lengths */
1859b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
186049b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
186149b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
186249b5e25fSSatish Balay 
186349b5e25fSSatish Balay   /* read in column indices */
1864b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
186549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
186649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
186749b5e25fSSatish Balay 
186849b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
186982502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1870a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
187182502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
187249b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
187349b5e25fSSatish Balay   masked   = mask + mbs;
187449b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
187549b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
187649b5e25fSSatish Balay     nmask = 0;
187749b5e25fSSatish Balay     for (j=0; j<bs; j++) {
187849b5e25fSSatish Balay       kmax = rowlengths[rowcount];
187949b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18802d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
188103630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
188249b5e25fSSatish Balay       }
188349b5e25fSSatish Balay       rowcount++;
188449b5e25fSSatish Balay     }
1885574b2666SHong Zhang     s_browlengths[i] += nmask;
1886574b2666SHong Zhang 
188749b5e25fSSatish Balay     /* zero out the mask elements we set */
188849b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
188949b5e25fSSatish Balay   }
189049b5e25fSSatish Balay 
189149b5e25fSSatish Balay   /* create our matrix */
18927fe2be48SHong Zhang   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr);
189349b5e25fSSatish Balay   B = *A;
189449b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
189549b5e25fSSatish Balay 
189649b5e25fSSatish Balay   /* set matrix "i" values */
189749b5e25fSSatish Balay   a->i[0] = 0;
189849b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1899574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1900574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
190149b5e25fSSatish Balay   }
19027fe2be48SHong Zhang   a->s_nz = a->i[mbs];
190349b5e25fSSatish Balay 
190449b5e25fSSatish Balay   /* read in nonzero values */
190587828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
190649b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
190749b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
190849b5e25fSSatish Balay 
190949b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
191049b5e25fSSatish Balay   nzcount = 0; jcount = 0;
191149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
191249b5e25fSSatish Balay     nzcountb = nzcount;
191349b5e25fSSatish Balay     nmask    = 0;
191449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
191549b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
191649b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19172d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
191803630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19192d703238SHong Zhang       }
19202d703238SHong Zhang     }
19212d703238SHong Zhang     /* sort the masked values */
19222d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19232d703238SHong Zhang 
19242d703238SHong Zhang     /* set "j" values into matrix */
19252d703238SHong Zhang     maskcount = 1;
19262d703238SHong Zhang     for (j=0; j<nmask; j++) {
192749b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
192849b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
192949b5e25fSSatish Balay     }
1930574b2666SHong Zhang 
193149b5e25fSSatish Balay     /* set "a" values into matrix */
193249b5e25fSSatish Balay     ishift = bs2*a->i[i];
193349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
193449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
193549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1936574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1937574b2666SHong Zhang         if (tmp >= i){
193849b5e25fSSatish Balay           block     = mask[tmp] - 1;
193949b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
194049b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1941574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1942574b2666SHong Zhang         }
1943574b2666SHong Zhang         nzcountb++;
194449b5e25fSSatish Balay       }
194549b5e25fSSatish Balay     }
194649b5e25fSSatish Balay     /* zero out the mask elements we set */
194749b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
194849b5e25fSSatish Balay   }
194929bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
195049b5e25fSSatish Balay 
195149b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1952574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
195349b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
195449b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
195549b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
195649b5e25fSSatish Balay 
195749b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
195827fa2452SHong Zhang #if defined(PETSC_HAVE_SPOOLES)
195927fa2452SHong Zhang   ierr = PetscOptionsHasName(B->prefix,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr);
196027fa2452SHong Zhang   if (flag) { ierr = MatUseSpooles_SeqSBAIJ(B);CHKERRQ(ierr); }
1961ecc4ab6dSBarry Smith #endif
196227f1aa60SHong Zhang #if defined(PETSC_HAVE_MUMPS)
196327f1aa60SHong Zhang   ierr = PetscOptionsHasName(B->prefix,"-mat_sbaij_mumps",&flag);CHKERRQ(ierr);
196427f1aa60SHong Zhang   if (flag) { ierr = MatUseMUMPS_MPIAIJ(B);CHKERRQ(ierr); }
196527f1aa60SHong Zhang #endif
196649b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
196749b5e25fSSatish Balay   PetscFunctionReturn(0);
196849b5e25fSSatish Balay }
19698549e402SHong Zhang EXTERN_C_END
1970574b2666SHong Zhang 
1971d06b337dSHong Zhang #undef __FUNCT__
1972d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1973c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1974d06b337dSHong Zhang {
1975d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1976d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1977d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
1978d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
1979d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1980d06b337dSHong Zhang 
1981d06b337dSHong Zhang   PetscFunctionBegin;
1982b965ef7fSBarry Smith   its = its*lits;
198391723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1984d06b337dSHong Zhang 
1985d06b337dSHong Zhang   if (bs > 1)
1986d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1987d06b337dSHong Zhang 
1988d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1989d06b337dSHong Zhang   if (xx != bb) {
1990d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1991d06b337dSHong Zhang   } else {
1992d06b337dSHong Zhang     b = x;
1993d06b337dSHong Zhang   }
1994d06b337dSHong Zhang 
1995d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1996d06b337dSHong Zhang 
1997d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19982798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1999d06b337dSHong Zhang       for (i=0; i<m; i++)
2000d06b337dSHong Zhang         t[i] = b[i];
2001d06b337dSHong Zhang 
2002d06b337dSHong Zhang       for (i=0; i<m; i++){
200344706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2004d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2005d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2006d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2007d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2008d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
200944706875SHong Zhang         PetscLogFlops(2*nz-1);
2010d06b337dSHong Zhang       }
2011d06b337dSHong Zhang     }
2012d06b337dSHong Zhang 
20132798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2014d06b337dSHong Zhang       for (i=0; i<m; i++)
2015d06b337dSHong Zhang         t[i] = b[i];
2016d06b337dSHong Zhang 
2017d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2018d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2019d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2020d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2021d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
202244706875SHong Zhang         PetscLogFlops(2*nz-1);
2023d06b337dSHong Zhang       }
2024d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2025d06b337dSHong Zhang         d  = *(aa + ai[i]);
2026d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2027d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2028d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2029d06b337dSHong Zhang         sum = t[i];
2030d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
203144706875SHong Zhang         PetscLogFlops(2*nz-1);
2032d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2033d06b337dSHong Zhang       }
2034d06b337dSHong Zhang     }
2035d06b337dSHong Zhang     its--;
2036d06b337dSHong Zhang   }
2037d06b337dSHong Zhang 
2038d06b337dSHong Zhang   while (its--) {
203944706875SHong Zhang     /*
204044706875SHong Zhang        forward sweep:
204144706875SHong Zhang        for i=0,...,m-1:
204244706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
204344706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
204444706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2045d06b337dSHong Zhang 
204644706875SHong Zhang     */
20472798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2048d06b337dSHong Zhang       for (i=0; i<m; i++)
2049d06b337dSHong Zhang         t[i] = b[i];
2050d06b337dSHong Zhang 
2051d06b337dSHong Zhang       for (i=0; i<m; i++){
205244706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2053d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2054d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2055d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2056d06b337dSHong Zhang         sum = t[i];
2057d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2058d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2059d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
206044706875SHong Zhang         PetscLogFlops(4*nz-2);
2061d06b337dSHong Zhang       }
2062d06b337dSHong Zhang     }
2063d06b337dSHong Zhang 
20642798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
206544706875SHong Zhang       /*
206644706875SHong Zhang        backward sweep:
206744706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
206844706875SHong Zhang        for i=m-1,...,0:
206944706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
207044706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
207144706875SHong Zhang       */
2072d06b337dSHong Zhang       for (i=0; i<m; i++)
2073d06b337dSHong Zhang         t[i] = b[i];
2074d06b337dSHong Zhang 
2075d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2076d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2077d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2078d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2079d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
208044706875SHong Zhang         PetscLogFlops(2*nz-1);
2081d06b337dSHong Zhang       }
2082d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2083d06b337dSHong Zhang         d  = *(aa + ai[i]);
2084d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2085d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2086d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2087d06b337dSHong Zhang         sum = t[i];
2088d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
208944706875SHong Zhang         PetscLogFlops(2*nz-1);
2090d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2091d06b337dSHong Zhang       }
2092d06b337dSHong Zhang     }
2093d06b337dSHong Zhang   }
2094d06b337dSHong Zhang 
2095d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
2096d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2097d06b337dSHong Zhang   if (bb != xx) {
2098d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2099d06b337dSHong Zhang   }
21002798e883SHong Zhang 
2101d06b337dSHong Zhang   PetscFunctionReturn(0);
2102d06b337dSHong Zhang }
2103d06b337dSHong Zhang 
2104d06b337dSHong Zhang 
2105d06b337dSHong Zhang 
2106d06b337dSHong Zhang 
210749b5e25fSSatish Balay 
210849b5e25fSSatish Balay 
2109