xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 4d9d31ab498596cb7ef6b2cbe5750038ce732ce3)
1*4d9d31abSKris Buschelman /*$Id: sbaij.c,v 1.57 2001/06/16 02:48:46 bsmith Exp buschelm $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
449b5e25fSSatish Balay     Defines the basic matrix operations for the BAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
7ab9f2c04SSatish Balay #include "petscsys.h"                            /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
949b5e25fSSatish Balay #include "src/vec/vecimpl.h"
1049b5e25fSSatish Balay #include "src/inline/spops.h"
11aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1249b5e25fSSatish Balay 
1349b5e25fSSatish Balay #define CHUNKSIZE  10
1449b5e25fSSatish Balay 
1549b5e25fSSatish Balay /*
1649b5e25fSSatish Balay      Checks for missing diagonals
1749b5e25fSSatish Balay */
184a2ae208SSatish Balay #undef __FUNCT__
194a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
2049b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A)
2149b5e25fSSatish Balay {
22045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2349b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2449b5e25fSSatish Balay 
2549b5e25fSSatish Balay   PetscFunctionBegin;
26045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2749b5e25fSSatish Balay   diag = a->diag;
2849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
2949b5e25fSSatish Balay     if (jj[diag[i]] != i) {
3029bbc08cSBarry Smith       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
3149b5e25fSSatish Balay     }
3249b5e25fSSatish Balay   }
3349b5e25fSSatish Balay   PetscFunctionReturn(0);
3449b5e25fSSatish Balay }
3549b5e25fSSatish Balay 
364a2ae208SSatish Balay #undef __FUNCT__
374a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
3849b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A)
3949b5e25fSSatish Balay {
40045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
4182502324SSatish Balay   int          i,j,*diag,m = a->mbs,ierr;
4249b5e25fSSatish Balay 
4349b5e25fSSatish Balay   PetscFunctionBegin;
4449b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4549b5e25fSSatish Balay 
46b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr);
47b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
4849b5e25fSSatish Balay   for (i=0; i<m; i++) {
4949b5e25fSSatish Balay     diag[i] = a->i[i+1];
5049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
5149b5e25fSSatish Balay       if (a->j[j] == i) {
5249b5e25fSSatish Balay         diag[i] = j;
5349b5e25fSSatish Balay         break;
5449b5e25fSSatish Balay       }
5549b5e25fSSatish Balay     }
5649b5e25fSSatish Balay   }
5749b5e25fSSatish Balay   a->diag = diag;
5849b5e25fSSatish Balay   PetscFunctionReturn(0);
5949b5e25fSSatish Balay }
6049b5e25fSSatish Balay 
6149b5e25fSSatish Balay 
6249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
6349b5e25fSSatish Balay 
644a2ae208SSatish Balay #undef __FUNCT__
654a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
6649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
6749b5e25fSSatish Balay {
6849b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6949b5e25fSSatish Balay 
7049b5e25fSSatish Balay   PetscFunctionBegin;
71045c9aa0SHong Zhang   if (ia) {
7229bbc08cSBarry Smith     SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
7349b5e25fSSatish Balay   }
74045c9aa0SHong Zhang   *nn = a->mbs;
7549b5e25fSSatish Balay   PetscFunctionReturn(0);
7649b5e25fSSatish Balay }
7749b5e25fSSatish Balay 
784a2ae208SSatish Balay #undef __FUNCT__
794a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
80045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
8149b5e25fSSatish Balay {
8249b5e25fSSatish Balay   PetscFunctionBegin;
8349b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
8429bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
8516cdd363SHong Zhang   /* PetscFunctionReturn(0); */
8649b5e25fSSatish Balay }
8749b5e25fSSatish Balay 
884a2ae208SSatish Balay #undef __FUNCT__
894a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ"
9049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
9149b5e25fSSatish Balay {
9249b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
9349b5e25fSSatish Balay 
9449b5e25fSSatish Balay   PetscFunctionBegin;
9549b5e25fSSatish Balay   *bs = sbaij->bs;
9649b5e25fSSatish Balay   PetscFunctionReturn(0);
9749b5e25fSSatish Balay }
9849b5e25fSSatish Balay 
994a2ae208SSatish Balay #undef __FUNCT__
1004a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
10149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
10249b5e25fSSatish Balay {
10349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10449b5e25fSSatish Balay   int         ierr;
10549b5e25fSSatish Balay 
10649b5e25fSSatish Balay   PetscFunctionBegin;
10749b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
108b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz);
10949b5e25fSSatish Balay #endif
11049b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
11149b5e25fSSatish Balay   if (!a->singlemalloc) {
11249b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
11363c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
11449b5e25fSSatish Balay   }
11549b5e25fSSatish Balay   if (a->row) {
11649b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
11749b5e25fSSatish Balay   }
11849b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
11949b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
12049b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
12149b5e25fSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
12249b5e25fSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
12349b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
12449b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
1251a3463dfSHong Zhang 
1261a3463dfSHong Zhang   if (a->inew){
1271a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1281a3463dfSHong Zhang     a->inew = 0;
1291a3463dfSHong Zhang   }
13049b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13149b5e25fSSatish Balay   PetscFunctionReturn(0);
13249b5e25fSSatish Balay }
13349b5e25fSSatish Balay 
1344a2ae208SSatish Balay #undef __FUNCT__
1354a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
13649b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
13749b5e25fSSatish Balay {
138045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
13949b5e25fSSatish Balay 
14049b5e25fSSatish Balay   PetscFunctionBegin;
141*4d9d31abSKris Buschelman   switch (op) {
142*4d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
143*4d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
144*4d9d31abSKris Buschelman     break;
145*4d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
146*4d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
147*4d9d31abSKris Buschelman     break;
148*4d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
149*4d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
150*4d9d31abSKris Buschelman     break;
151*4d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
152*4d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
153*4d9d31abSKris Buschelman     break;
154*4d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
155*4d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
156*4d9d31abSKris Buschelman     break;
157*4d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
158*4d9d31abSKris Buschelman     a->nonew          = 1;
159*4d9d31abSKris Buschelman     break;
160*4d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
161*4d9d31abSKris Buschelman     a->nonew          = -1;
162*4d9d31abSKris Buschelman     break;
163*4d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
164*4d9d31abSKris Buschelman     a->nonew          = -2;
165*4d9d31abSKris Buschelman     break;
166*4d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
167*4d9d31abSKris Buschelman     a->nonew          = 0;
168*4d9d31abSKris Buschelman     break;
169*4d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
170*4d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
171*4d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
172*4d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
173*4d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
174b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
175*4d9d31abSKris Buschelman     break;
176*4d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
17729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
178*4d9d31abSKris Buschelman     break;
179*4d9d31abSKris Buschelman   default:
18029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
181*4d9d31abSKris Buschelman     break;
18249b5e25fSSatish Balay   }
18349b5e25fSSatish Balay   PetscFunctionReturn(0);
18449b5e25fSSatish Balay }
18549b5e25fSSatish Balay 
1864a2ae208SSatish Balay #undef __FUNCT__
1874a2ae208SSatish Balay #define __FUNCT__ "MatGetOwnershipRange_SeqSBAIJ"
18849b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n)
18949b5e25fSSatish Balay {
19049b5e25fSSatish Balay   PetscFunctionBegin;
191c464158bSHong Zhang   *m = 0;
1927e8a2ce2SHong Zhang   *n = A->m;
19349b5e25fSSatish Balay   PetscFunctionReturn(0);
19449b5e25fSSatish Balay }
19549b5e25fSSatish Balay 
1964a2ae208SSatish Balay #undef __FUNCT__
1974a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19849b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v)
19949b5e25fSSatish Balay {
20049b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
20182502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
20249b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
20349b5e25fSSatish Balay   Scalar       *v_i;
20449b5e25fSSatish Balay 
20549b5e25fSSatish Balay   PetscFunctionBegin;
20649b5e25fSSatish Balay   bs  = a->bs;
20749b5e25fSSatish Balay   ai  = a->i;
20849b5e25fSSatish Balay   aj  = a->j;
20949b5e25fSSatish Balay   aa  = a->a;
21049b5e25fSSatish Balay   bs2 = a->bs2;
21149b5e25fSSatish Balay 
212c464158bSHong Zhang   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
21349b5e25fSSatish Balay 
21449b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21549b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21649b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
21749b5e25fSSatish Balay   *ncols = bs*M;
21849b5e25fSSatish Balay 
21949b5e25fSSatish Balay   if (v) {
22049b5e25fSSatish Balay     *v = 0;
22149b5e25fSSatish Balay     if (*ncols) {
22282502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(Scalar),v);CHKERRQ(ierr);
22349b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22449b5e25fSSatish Balay         v_i  = *v + i*bs;
22549b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22649b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
22749b5e25fSSatish Balay       }
22849b5e25fSSatish Balay     }
22949b5e25fSSatish Balay   }
23049b5e25fSSatish Balay 
23149b5e25fSSatish Balay   if (cols) {
23249b5e25fSSatish Balay     *cols = 0;
23349b5e25fSSatish Balay     if (*ncols) {
23482502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
23549b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23649b5e25fSSatish Balay         cols_i = *cols + i*bs;
23749b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
23849b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
23949b5e25fSSatish Balay       }
24049b5e25fSSatish Balay     }
24149b5e25fSSatish Balay   }
24249b5e25fSSatish Balay 
24349b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2445ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2455ddb2528SHong Zhang #ifdef column_search
24649b5e25fSSatish Balay   v_i    = *v    + M*bs;
24749b5e25fSSatish Balay   cols_i = *cols + M*bs;
24849b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
24949b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25049b5e25fSSatish Balay     for (j=0; j<M; j++){
25149b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25249b5e25fSSatish Balay       if (itmp == bn){
25349b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25449b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25549b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25649b5e25fSSatish Balay           *v_i++    = aa_i[k];
25749b5e25fSSatish Balay         }
25849b5e25fSSatish Balay         *ncols += bs;
25949b5e25fSSatish Balay         break;
26049b5e25fSSatish Balay       }
26149b5e25fSSatish Balay     }
26249b5e25fSSatish Balay   }
2635ddb2528SHong Zhang #endif
26449b5e25fSSatish Balay 
26549b5e25fSSatish Balay   PetscFunctionReturn(0);
26649b5e25fSSatish Balay }
26749b5e25fSSatish Balay 
2684a2ae208SSatish Balay #undef __FUNCT__
2694a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
27049b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
27149b5e25fSSatish Balay {
27249b5e25fSSatish Balay   int ierr;
27349b5e25fSSatish Balay 
27449b5e25fSSatish Balay   PetscFunctionBegin;
27549b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27649b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
27749b5e25fSSatish Balay   PetscFunctionReturn(0);
27849b5e25fSSatish Balay }
27949b5e25fSSatish Balay 
2804a2ae208SSatish Balay #undef __FUNCT__
2814a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
28249b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
28349b5e25fSSatish Balay {
28449b5e25fSSatish Balay   PetscFunctionBegin;
28529bbc08cSBarry Smith   SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called");
28616cdd363SHong Zhang   /* PetscFunctionReturn(0); */
28749b5e25fSSatish Balay }
28849b5e25fSSatish Balay 
2894a2ae208SSatish Balay #undef __FUNCT__
2904a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary"
291b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
29249b5e25fSSatish Balay {
29349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
29449b5e25fSSatish Balay   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
29549b5e25fSSatish Balay   Scalar       *aa;
29649b5e25fSSatish Balay   FILE         *file;
29749b5e25fSSatish Balay 
29849b5e25fSSatish Balay   PetscFunctionBegin;
299b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
30082502324SSatish Balay   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
30149b5e25fSSatish Balay   col_lens[0] = MAT_COOKIE;
30249b5e25fSSatish Balay 
303c464158bSHong Zhang   col_lens[1] = A->m;
304c464158bSHong Zhang   col_lens[2] = A->m;
30549b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
30649b5e25fSSatish Balay 
30749b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
30849b5e25fSSatish Balay   count = 0;
30949b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31149b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
31249b5e25fSSatish Balay     }
31349b5e25fSSatish Balay   }
314c464158bSHong Zhang   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
31549b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
31649b5e25fSSatish Balay 
31749b5e25fSSatish Balay   /* store column indices (zero start index) */
31882502324SSatish Balay   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
31949b5e25fSSatish Balay   count = 0;
32049b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
32149b5e25fSSatish Balay     for (j=0; j<bs; j++) {
32249b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
32349b5e25fSSatish Balay         for (l=0; l<bs; l++) {
32449b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
32549b5e25fSSatish Balay         }
32649b5e25fSSatish Balay       }
32749b5e25fSSatish Balay     }
32849b5e25fSSatish Balay   }
32949b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
33049b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
33149b5e25fSSatish Balay 
33249b5e25fSSatish Balay   /* store nonzero values */
33382502324SSatish Balay   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar),&aa);CHKERRQ(ierr);
33449b5e25fSSatish Balay   count = 0;
33549b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
33649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
33749b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
33849b5e25fSSatish Balay         for (l=0; l<bs; l++) {
33949b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
34049b5e25fSSatish Balay         }
34149b5e25fSSatish Balay       }
34249b5e25fSSatish Balay     }
34349b5e25fSSatish Balay   }
34449b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
34549b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
34649b5e25fSSatish Balay 
347b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
34849b5e25fSSatish Balay   if (file) {
34949b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
35049b5e25fSSatish Balay   }
35149b5e25fSSatish Balay   PetscFunctionReturn(0);
35249b5e25fSSatish Balay }
35349b5e25fSSatish Balay 
3544a2ae208SSatish Balay #undef __FUNCT__
3554a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
356b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
35749b5e25fSSatish Balay {
35849b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
359fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
360fb9695e5SSatish Balay   char              *name;
361f3ef73ceSBarry Smith   PetscViewerFormat format;
36249b5e25fSSatish Balay 
36349b5e25fSSatish Balay   PetscFunctionBegin;
36480fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
365b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
366fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
367b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
368fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
36929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
370fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
371b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
37249b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
37349b5e25fSSatish Balay       for (j=0; j<bs; j++) {
374b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
37549b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
37649b5e25fSSatish Balay           for (l=0; l<bs; l++) {
37749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
37849b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
379b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
38049b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38149b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
382b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
38349b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38449b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
385b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38649b5e25fSSatish Balay             }
38749b5e25fSSatish Balay #else
38849b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
389b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
39049b5e25fSSatish Balay             }
39149b5e25fSSatish Balay #endif
39249b5e25fSSatish Balay           }
39349b5e25fSSatish Balay         }
394b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39549b5e25fSSatish Balay       }
39649b5e25fSSatish Balay     }
397b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
39849b5e25fSSatish Balay   } else {
399b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
40049b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
40149b5e25fSSatish Balay       for (j=0; j<bs; j++) {
402b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
40349b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
40449b5e25fSSatish Balay           for (l=0; l<bs; l++) {
40549b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
40649b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
407b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
40849b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40949b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
410b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
41149b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41249b5e25fSSatish Balay             } else {
413b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41449b5e25fSSatish Balay             }
41549b5e25fSSatish Balay #else
416b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
41749b5e25fSSatish Balay #endif
41849b5e25fSSatish Balay           }
41949b5e25fSSatish Balay         }
420b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
42149b5e25fSSatish Balay       }
42249b5e25fSSatish Balay     }
423b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
42449b5e25fSSatish Balay   }
425b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
42649b5e25fSSatish Balay   PetscFunctionReturn(0);
42749b5e25fSSatish Balay }
42849b5e25fSSatish Balay 
4294a2ae208SSatish Balay #undef __FUNCT__
4304a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
431b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
43249b5e25fSSatish Balay {
43349b5e25fSSatish Balay   Mat          A = (Mat) Aa;
43449b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
43549b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
43649b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
43749b5e25fSSatish Balay   MatScalar    *aa;
43849b5e25fSSatish Balay   MPI_Comm     comm;
439b0a32e0cSBarry Smith   PetscViewer  viewer;
44049b5e25fSSatish Balay 
44149b5e25fSSatish Balay   PetscFunctionBegin;
44249b5e25fSSatish Balay   /*
44349b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
44449b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
44549b5e25fSSatish Balay    rest should return immediately.
44649b5e25fSSatish Balay   */
44749b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
44849b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
44949b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
45049b5e25fSSatish Balay 
45149b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
45249b5e25fSSatish Balay 
453b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
454b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
45549b5e25fSSatish Balay 
45649b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
457b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
45849b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45949b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
460c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
46149b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46249b5e25fSSatish Balay       aa = a->a + j*bs2;
46349b5e25fSSatish Balay       for (k=0; k<bs; k++) {
46449b5e25fSSatish Balay         for (l=0; l<bs; l++) {
46549b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
466b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
46749b5e25fSSatish Balay         }
46849b5e25fSSatish Balay       }
46949b5e25fSSatish Balay     }
47049b5e25fSSatish Balay   }
471b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
47249b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
47349b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
474c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
47549b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
47649b5e25fSSatish Balay       aa = a->a + j*bs2;
47749b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47849b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47949b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
480b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
48149b5e25fSSatish Balay         }
48249b5e25fSSatish Balay       }
48349b5e25fSSatish Balay     }
48449b5e25fSSatish Balay   }
48549b5e25fSSatish Balay 
486b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
48749b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
48849b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
489c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
49049b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
49149b5e25fSSatish Balay       aa = a->a + j*bs2;
49249b5e25fSSatish Balay       for (k=0; k<bs; k++) {
49349b5e25fSSatish Balay         for (l=0; l<bs; l++) {
49449b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
495b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
49649b5e25fSSatish Balay         }
49749b5e25fSSatish Balay       }
49849b5e25fSSatish Balay     }
49949b5e25fSSatish Balay   }
50049b5e25fSSatish Balay   PetscFunctionReturn(0);
50149b5e25fSSatish Balay }
50249b5e25fSSatish Balay 
5034a2ae208SSatish Balay #undef __FUNCT__
5044a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
505b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
50649b5e25fSSatish Balay {
50749b5e25fSSatish Balay   int          ierr;
50849b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
509b0a32e0cSBarry Smith   PetscDraw    draw;
51049b5e25fSSatish Balay   PetscTruth   isnull;
51149b5e25fSSatish Balay 
51249b5e25fSSatish Balay   PetscFunctionBegin;
51349b5e25fSSatish Balay 
514b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
515b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
51649b5e25fSSatish Balay 
51749b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
518c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
51949b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
520b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
521b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
52249b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
52349b5e25fSSatish Balay   PetscFunctionReturn(0);
52449b5e25fSSatish Balay }
52549b5e25fSSatish Balay 
5264a2ae208SSatish Balay #undef __FUNCT__
5274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
528b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
52949b5e25fSSatish Balay {
53049b5e25fSSatish Balay   int        ierr;
53149b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
53249b5e25fSSatish Balay 
53349b5e25fSSatish Balay   PetscFunctionBegin;
534b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
535b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
536fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
537fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
53849b5e25fSSatish Balay   if (issocket) {
53929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
54049b5e25fSSatish Balay   } else if (isascii){
54149b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
54249b5e25fSSatish Balay   } else if (isbinary) {
54349b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
54449b5e25fSSatish Balay   } else if (isdraw) {
54549b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
54649b5e25fSSatish Balay   } else {
54729bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
54849b5e25fSSatish Balay   }
54949b5e25fSSatish Balay   PetscFunctionReturn(0);
55049b5e25fSSatish Balay }
55149b5e25fSSatish Balay 
55249b5e25fSSatish Balay 
5534a2ae208SSatish Balay #undef __FUNCT__
5544a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
55549b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
55649b5e25fSSatish Balay {
557045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
55849b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
55949b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
56049b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
56149b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
56249b5e25fSSatish Balay 
56349b5e25fSSatish Balay   PetscFunctionBegin;
56449b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
56549b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
56629bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
567c464158bSHong Zhang     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
56849b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
56949b5e25fSSatish Balay     nrow = ailen[brow];
57049b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
57129bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
572c464158bSHong Zhang       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
57349b5e25fSSatish Balay       col  = in[l] ;
57449b5e25fSSatish Balay       bcol = col/bs;
57549b5e25fSSatish Balay       cidx = col%bs;
57649b5e25fSSatish Balay       ridx = row%bs;
57749b5e25fSSatish Balay       high = nrow;
57849b5e25fSSatish Balay       low  = 0; /* assume unsorted */
57949b5e25fSSatish Balay       while (high-low > 5) {
58049b5e25fSSatish Balay         t = (low+high)/2;
58149b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
58249b5e25fSSatish Balay         else             low  = t;
58349b5e25fSSatish Balay       }
58449b5e25fSSatish Balay       for (i=low; i<high; i++) {
58549b5e25fSSatish Balay         if (rp[i] > bcol) break;
58649b5e25fSSatish Balay         if (rp[i] == bcol) {
58749b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
58849b5e25fSSatish Balay           goto finished;
58949b5e25fSSatish Balay         }
59049b5e25fSSatish Balay       }
59149b5e25fSSatish Balay       *v++ = zero;
59249b5e25fSSatish Balay       finished:;
59349b5e25fSSatish Balay     }
59449b5e25fSSatish Balay   }
59549b5e25fSSatish Balay   PetscFunctionReturn(0);
59649b5e25fSSatish Balay }
59749b5e25fSSatish Balay 
59849b5e25fSSatish Balay 
5994a2ae208SSatish Balay #undef __FUNCT__
6004a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
60149b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
60249b5e25fSSatish Balay {
60349b5e25fSSatish Balay   PetscFunctionBegin;
60429bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format");
60516cdd363SHong Zhang   /* PetscFunctionReturn(0); */
60649b5e25fSSatish Balay }
60749b5e25fSSatish Balay 
6084a2ae208SSatish Balay #undef __FUNCT__
6094a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
61049b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
61149b5e25fSSatish Balay {
61249b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
61349b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
614c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
61549b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
61649b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
61749b5e25fSSatish Balay 
61849b5e25fSSatish Balay   PetscFunctionBegin;
61949b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
62049b5e25fSSatish Balay 
62149b5e25fSSatish Balay   if (m) rmax = ailen[0];
62249b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
62349b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
62449b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
62549b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
62649b5e25fSSatish Balay     if (fshift) {
62749b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
62849b5e25fSSatish Balay       N = ailen[i];
62949b5e25fSSatish Balay       for (j=0; j<N; j++) {
63049b5e25fSSatish Balay         ip[j-fshift] = ip[j];
63149b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
63249b5e25fSSatish Balay       }
63349b5e25fSSatish Balay     }
63449b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
63549b5e25fSSatish Balay   }
63649b5e25fSSatish Balay   if (mbs) {
63749b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
63849b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
63949b5e25fSSatish Balay   }
64049b5e25fSSatish Balay   /* reset ilen and imax for each row */
64149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
64249b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
64349b5e25fSSatish Balay   }
64449b5e25fSSatish Balay   a->s_nz = ai[mbs];
64549b5e25fSSatish Balay 
64649b5e25fSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
64749b5e25fSSatish Balay   if (fshift && a->diag) {
64849b5e25fSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
649b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-(m+1)*sizeof(int));
65049b5e25fSSatish Balay     a->diag = 0;
65149b5e25fSSatish Balay   }
652b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
653c464158bSHong Zhang            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
654b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
65549b5e25fSSatish Balay            a->reallocs);
656b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
65749b5e25fSSatish Balay   a->reallocs          = 0;
65849b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
65949b5e25fSSatish Balay 
66049b5e25fSSatish Balay   PetscFunctionReturn(0);
66149b5e25fSSatish Balay }
66249b5e25fSSatish Balay 
66349b5e25fSSatish Balay /*
66449b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
66549b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
66649b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
66749b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
66849b5e25fSSatish Balay */
6694a2ae208SSatish Balay #undef __FUNCT__
6704a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
67149b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
67249b5e25fSSatish Balay {
67349b5e25fSSatish Balay   int        i,j,k,row;
67449b5e25fSSatish Balay   PetscTruth flg;
67549b5e25fSSatish Balay 
67649b5e25fSSatish Balay   PetscFunctionBegin;
67749b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
67849b5e25fSSatish Balay     row = idx[i];
67949b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
68049b5e25fSSatish Balay       sizes[j] = 1;
68149b5e25fSSatish Balay       i++;
68249b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
68349b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
68449b5e25fSSatish Balay       i++;
68549b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
68649b5e25fSSatish Balay       flg = PETSC_TRUE;
68749b5e25fSSatish Balay       for (k=1; k<bs; k++) {
68849b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
68949b5e25fSSatish Balay           flg = PETSC_FALSE;
69049b5e25fSSatish Balay           break;
69149b5e25fSSatish Balay         }
69249b5e25fSSatish Balay       }
69349b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
69449b5e25fSSatish Balay         sizes[j] = bs;
69549b5e25fSSatish Balay         i+= bs;
69649b5e25fSSatish Balay       } else {
69749b5e25fSSatish Balay         sizes[j] = 1;
69849b5e25fSSatish Balay         i++;
69949b5e25fSSatish Balay       }
70049b5e25fSSatish Balay     }
70149b5e25fSSatish Balay   }
70249b5e25fSSatish Balay   *bs_max = j;
70349b5e25fSSatish Balay   PetscFunctionReturn(0);
70449b5e25fSSatish Balay }
70549b5e25fSSatish Balay 
7064a2ae208SSatish Balay #undef __FUNCT__
7074a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
70849b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag)
70949b5e25fSSatish Balay {
71049b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data;
71149b5e25fSSatish Balay   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
71249b5e25fSSatish Balay   int         bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
71349b5e25fSSatish Balay   Scalar      zero = 0.0;
71449b5e25fSSatish Balay   MatScalar   *aa;
71549b5e25fSSatish Balay 
71649b5e25fSSatish Balay   PetscFunctionBegin;
71749b5e25fSSatish Balay   /* Make a copy of the IS and  sort it */
71849b5e25fSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
71949b5e25fSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
72049b5e25fSSatish Balay 
72149b5e25fSSatish Balay   /* allocate memory for rows,sizes */
722b0a32e0cSBarry Smith   ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
72349b5e25fSSatish Balay   sizes = rows + is_n;
72449b5e25fSSatish Balay 
72549b5e25fSSatish Balay   /* initialize copy IS values to rows, and sort them */
72649b5e25fSSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
72749b5e25fSSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
72849b5e25fSSatish Balay   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
72949b5e25fSSatish Balay     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
73049b5e25fSSatish Balay     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
73149b5e25fSSatish Balay   } else {
73249b5e25fSSatish Balay     ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
73349b5e25fSSatish Balay   }
73449b5e25fSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
73549b5e25fSSatish Balay 
73649b5e25fSSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
73749b5e25fSSatish Balay     row   = rows[j];                  /* row to be zeroed */
738c464158bSHong Zhang     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
73949b5e25fSSatish Balay     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
74049b5e25fSSatish Balay     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
74149b5e25fSSatish Balay     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
74249b5e25fSSatish Balay       if (diag) {
74349b5e25fSSatish Balay         if (sbaij->ilen[row/bs] > 0) {
74449b5e25fSSatish Balay           sbaij->ilen[row/bs] = 1;
74549b5e25fSSatish Balay           sbaij->j[sbaij->i[row/bs]] = row/bs;
74649b5e25fSSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
74749b5e25fSSatish Balay         }
74849b5e25fSSatish Balay         /* Now insert all the diagoanl values for this bs */
74949b5e25fSSatish Balay         for (k=0; k<bs; k++) {
75049b5e25fSSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
75149b5e25fSSatish Balay         }
75249b5e25fSSatish Balay       } else { /* (!diag) */
75349b5e25fSSatish Balay         sbaij->ilen[row/bs] = 0;
75449b5e25fSSatish Balay       } /* end (!diag) */
75549b5e25fSSatish Balay     } else { /* (sizes[i] != bs), broken block */
75649b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG)
75729bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
75849b5e25fSSatish Balay #endif
75949b5e25fSSatish Balay       for (k=0; k<count; k++) {
76049b5e25fSSatish Balay         aa[0] = zero;
76149b5e25fSSatish Balay         aa+=bs;
76249b5e25fSSatish Balay       }
76349b5e25fSSatish Balay       if (diag) {
76449b5e25fSSatish Balay         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
76549b5e25fSSatish Balay       }
76649b5e25fSSatish Balay     }
76749b5e25fSSatish Balay   }
76849b5e25fSSatish Balay 
76949b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
77049b5e25fSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77149b5e25fSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77249b5e25fSSatish Balay   PetscFunctionReturn(0);
77349b5e25fSSatish Balay }
77449b5e25fSSatish Balay 
77549b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
77649b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
77749b5e25fSSatish Balay */
77849b5e25fSSatish Balay 
7794a2ae208SSatish Balay #undef __FUNCT__
7804a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
78149b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
78249b5e25fSSatish Balay {
78349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
78449b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
78549b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
78649b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
78749b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
78849b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
78949b5e25fSSatish Balay 
79049b5e25fSSatish Balay   PetscFunctionBegin;
79149b5e25fSSatish Balay 
79249b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
79349b5e25fSSatish Balay     row  = im[k];       /* row number */
79449b5e25fSSatish Balay     brow = row/bs;      /* block row number */
79549b5e25fSSatish Balay     if (row < 0) continue;
79649b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
797c464158bSHong Zhang     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
79849b5e25fSSatish Balay #endif
79949b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
80049b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
80149b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
80249b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
80349b5e25fSSatish Balay     low  = 0;
80449b5e25fSSatish Balay 
80549b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
80649b5e25fSSatish Balay       if (in[l] < 0) continue;
80749b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
808c464158bSHong Zhang       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m);
80949b5e25fSSatish Balay #endif
81049b5e25fSSatish Balay       col = in[l];
81149b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
81249b5e25fSSatish Balay 
81349b5e25fSSatish Balay       if (brow <= bcol){
81449b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8158549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
81649b5e25fSSatish Balay           /* element value a(k,l) */
81749b5e25fSSatish Balay           if (roworiented) {
81849b5e25fSSatish Balay             value = v[l + k*n];
81949b5e25fSSatish Balay           } else {
82049b5e25fSSatish Balay             value = v[k + l*m];
82149b5e25fSSatish Balay           }
82249b5e25fSSatish Balay 
82349b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
82449b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
82549b5e25fSSatish Balay           while (high-low > 7) {
82649b5e25fSSatish Balay             t = (low+high)/2;
82749b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
82849b5e25fSSatish Balay             else              low  = t;
82949b5e25fSSatish Balay           }
83049b5e25fSSatish Balay           for (i=low; i<high; i++) {
83149b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
83249b5e25fSSatish Balay             if (rp[i] > bcol) break;
83349b5e25fSSatish Balay             if (rp[i] == bcol) {
83449b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
83549b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
83649b5e25fSSatish Balay               else                  *bap  = value;
8378549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8388549e402SHong Zhang               if (brow == bcol && ridx < cidx){
8398549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
8408549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
8418549e402SHong Zhang                 else                  *bap  = value;
8428549e402SHong Zhang               }
84349b5e25fSSatish Balay               goto noinsert1;
84449b5e25fSSatish Balay             }
84549b5e25fSSatish Balay           }
84649b5e25fSSatish Balay 
84749b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
84829bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
84949b5e25fSSatish Balay       if (nrow >= rmax) {
85049b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
85149b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
85249b5e25fSSatish Balay         MatScalar *new_a;
85349b5e25fSSatish Balay 
85429bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
85549b5e25fSSatish Balay 
85649b5e25fSSatish Balay         /* Malloc new storage space */
85749b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
85882502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
85949b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
86049b5e25fSSatish Balay         new_i = new_j + new_nz;
86149b5e25fSSatish Balay 
86249b5e25fSSatish Balay         /* copy over old data into new slots */
86349b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
86449b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
86549b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
86649b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
86749b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
86849b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
86949b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
87049b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
87149b5e25fSSatish Balay         /* free up old matrix storage */
87249b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
87349b5e25fSSatish Balay         if (!a->singlemalloc) {
87449b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
87549b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
87649b5e25fSSatish Balay         }
87749b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
87849b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
87949b5e25fSSatish Balay 
88049b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
88149b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
882b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
88349b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
88449b5e25fSSatish Balay         a->reallocs++;
88549b5e25fSSatish Balay         a->s_nz++;
88649b5e25fSSatish Balay       }
88749b5e25fSSatish Balay 
88849b5e25fSSatish Balay       N = nrow++ - 1;
88949b5e25fSSatish Balay       /* shift up all the later entries in this row */
89049b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
89149b5e25fSSatish Balay         rp[ii+1] = rp[ii];
89249b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
89349b5e25fSSatish Balay       }
89449b5e25fSSatish Balay       if (N>=i) {
89549b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
89649b5e25fSSatish Balay       }
89749b5e25fSSatish Balay       rp[i]                      = bcol;
89849b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
89949b5e25fSSatish Balay       noinsert1:;
90049b5e25fSSatish Balay       low = i;
90149b5e25fSSatish Balay       /* } */
9028549e402SHong Zhang         }
90349b5e25fSSatish Balay       } /* end of if .. if..  */
90449b5e25fSSatish Balay     }                     /* end of loop over added columns */
90549b5e25fSSatish Balay     ailen[brow] = nrow;
90649b5e25fSSatish Balay   }                       /* end of loop over added rows */
90749b5e25fSSatish Balay 
90849b5e25fSSatish Balay   PetscFunctionReturn(0);
90949b5e25fSSatish Balay }
91049b5e25fSSatish Balay 
911915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
912915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
91349b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
91449b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
91549b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
91649b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
91749b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
91849b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat);
91949b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
92049b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
92149b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
92249b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
92349b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
92449b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
92513a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
92649b5e25fSSatish Balay 
92749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
92849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
92949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
93049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
93149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
93249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
93349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
93449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
93549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
93649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
93749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
93849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
93949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
94049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
94149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
94249b5e25fSSatish Balay 
94307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
94407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
94507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
94607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
94707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
94807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
94907f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
95007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
95107f98182SSatish Balay 
9525f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
9535f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
9545f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
9555f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
9565f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
9575f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
9585f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
95957d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
96049b5e25fSSatish Balay 
96149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
96249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
96349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
96449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
96549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
96649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
96749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
96849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
96949b5e25fSSatish Balay 
97049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
97149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
97249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
97349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
97449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
97549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
97649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
97749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
97849b5e25fSSatish Balay 
9794d101231SSatish Balay #ifdef HAVE_MatICCFactor
9801a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
9814a2ae208SSatish Balay #undef __FUNCT__
9824d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
9834d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level)
98449b5e25fSSatish Balay {
9854ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
98649b5e25fSSatish Balay   Mat         outA;
98749b5e25fSSatish Balay   int         ierr;
98849b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
98949b5e25fSSatish Balay 
99049b5e25fSSatish Balay   PetscFunctionBegin;
9911a3463dfSHong Zhang   /*
99229bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
99349b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
99449b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
99549b5e25fSSatish Balay   if (!row_identity || !col_identity) {
99629bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
99749b5e25fSSatish Balay   }
9981a3463dfSHong Zhang   */
99949b5e25fSSatish Balay 
100049b5e25fSSatish Balay   outA          = inA;
10011260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
100249b5e25fSSatish Balay 
100349b5e25fSSatish Balay   if (!a->diag) {
10041a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
100549b5e25fSSatish Balay   }
100649b5e25fSSatish Balay   /*
100749b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
100849b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
100949b5e25fSSatish Balay   */
101049b5e25fSSatish Balay   switch (a->bs) {
101149b5e25fSSatish Balay   case 1:
10120fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
101307f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
10144d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
101549b5e25fSSatish Balay   case 2:
10161a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
10171a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10181a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10194d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
102049b5e25fSSatish Balay     break;
102149b5e25fSSatish Balay   case 3:
10221a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
10231a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10241a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10254d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
102649b5e25fSSatish Balay     break;
102749b5e25fSSatish Balay   case 4:
10281a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
10291a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10301a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10314d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
103249b5e25fSSatish Balay     break;
103349b5e25fSSatish Balay   case 5:
10341a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
10351a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10361a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10374d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
103849b5e25fSSatish Balay     break;
103949b5e25fSSatish Balay   case 6:
10401a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
10411a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10421a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10434d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
104449b5e25fSSatish Balay     break;
104549b5e25fSSatish Balay   case 7:
10461a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
10471a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10481a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10494d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
105049b5e25fSSatish Balay     break;
105149b5e25fSSatish Balay   default:
105249b5e25fSSatish Balay     a->row        = row;
10531a3463dfSHong Zhang     a->icol       = col;
105449b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
105549b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
105649b5e25fSSatish Balay 
105749b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
105849b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1059b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
106049b5e25fSSatish Balay 
106149b5e25fSSatish Balay     if (!a->solve_work) {
106282502324SSatish Balay       ierr = PetscMalloc((A->m+a->bs)*sizeof(Scalar),&a->solve_work);CHKERRQ(ierr);
1063b0a32e0cSBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(Scalar));
106449b5e25fSSatish Balay     }
106549b5e25fSSatish Balay   }
106649b5e25fSSatish Balay 
10671a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
106849b5e25fSSatish Balay 
106949b5e25fSSatish Balay   PetscFunctionReturn(0);
107049b5e25fSSatish Balay }
1071950f1e5bSHong Zhang #endif
1072950f1e5bSHong Zhang 
10734a2ae208SSatish Balay #undef __FUNCT__
10744a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
107549b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
107649b5e25fSSatish Balay {
107749b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
107849b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
107949b5e25fSSatish Balay   int               ierr;
108049b5e25fSSatish Balay 
108149b5e25fSSatish Balay   PetscFunctionBegin;
108249b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
108349b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
108449b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
108549b5e25fSSatish Balay   PetscFunctionReturn(0);
108649b5e25fSSatish Balay }
108749b5e25fSSatish Balay 
108849b5e25fSSatish Balay EXTERN_C_BEGIN
10894a2ae208SSatish Balay #undef __FUNCT__
10904a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
109149b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
109249b5e25fSSatish Balay {
1093045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
109449b5e25fSSatish Balay   int         i,nz,n;
109549b5e25fSSatish Balay 
109649b5e25fSSatish Balay   PetscFunctionBegin;
1097045c9aa0SHong Zhang   nz = baij->s_maxnz;
1098c464158bSHong Zhang   n  = mat->n;
109949b5e25fSSatish Balay   for (i=0; i<nz; i++) {
110049b5e25fSSatish Balay     baij->j[i] = indices[i];
110149b5e25fSSatish Balay   }
1102045c9aa0SHong Zhang   baij->s_nz = nz;
110349b5e25fSSatish Balay   for (i=0; i<n; i++) {
110449b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
110549b5e25fSSatish Balay   }
110649b5e25fSSatish Balay 
110749b5e25fSSatish Balay   PetscFunctionReturn(0);
110849b5e25fSSatish Balay }
110949b5e25fSSatish Balay EXTERN_C_END
111049b5e25fSSatish Balay 
11114a2ae208SSatish Balay #undef __FUNCT__
11124a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
111349b5e25fSSatish Balay /*@
111419585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
111549b5e25fSSatish Balay        in the matrix.
111649b5e25fSSatish Balay 
111749b5e25fSSatish Balay   Input Parameters:
111819585528SSatish Balay +  mat     - the SeqSBAIJ matrix
111949b5e25fSSatish Balay -  indices - the column indices
112049b5e25fSSatish Balay 
112149b5e25fSSatish Balay   Level: advanced
112249b5e25fSSatish Balay 
112349b5e25fSSatish Balay   Notes:
112449b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
112549b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
112649b5e25fSSatish Balay   of the MatSetValues() operation.
112749b5e25fSSatish Balay 
112849b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
112919585528SSatish Balay   MatCreateSeqSBAIJ().
113049b5e25fSSatish Balay 
1131ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
113249b5e25fSSatish Balay 
1133ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
113449b5e25fSSatish Balay @*/
113549b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
113649b5e25fSSatish Balay {
113749b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
113849b5e25fSSatish Balay 
113949b5e25fSSatish Balay   PetscFunctionBegin;
114049b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1141b9617806SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)())&f);CHKERRQ(ierr);
114249b5e25fSSatish Balay   if (f) {
114349b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
114449b5e25fSSatish Balay   } else {
114529bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
114649b5e25fSSatish Balay   }
114749b5e25fSSatish Balay   PetscFunctionReturn(0);
114849b5e25fSSatish Balay }
114949b5e25fSSatish Balay 
11504a2ae208SSatish Balay #undef __FUNCT__
11514a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1152273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1153273d9f13SBarry Smith {
1154273d9f13SBarry Smith   int        ierr;
1155273d9f13SBarry Smith 
1156273d9f13SBarry Smith   PetscFunctionBegin;
1157273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1158273d9f13SBarry Smith   PetscFunctionReturn(0);
1159273d9f13SBarry Smith }
1160273d9f13SBarry Smith 
116149b5e25fSSatish Balay /* -------------------------------------------------------------------*/
116249b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
116349b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
116449b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
116549b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
116649b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
116749b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
116849b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
116949b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
117049b5e25fSSatish Balay        0,
117149b5e25fSSatish Balay        0,
117249b5e25fSSatish Balay        0,
117349b5e25fSSatish Balay        0,
11745f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
117549b5e25fSSatish Balay        0,
117649b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
117749b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
117849b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
117949b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
118049b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
118149b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
118249b5e25fSSatish Balay        0,
118349b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
118449b5e25fSSatish Balay        0,
118549b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
118649b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
118749b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
118849b5e25fSSatish Balay        0,
118949b5e25fSSatish Balay        0,
11905f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
11915f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
1192273d9f13SBarry Smith        MatSetUpPreallocation_SeqSBAIJ,
1193c464158bSHong Zhang        0,
119449b5e25fSSatish Balay        MatGetOwnershipRange_SeqSBAIJ,
119549b5e25fSSatish Balay        0,
11964d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
119749b5e25fSSatish Balay        0,
119849b5e25fSSatish Balay        0,
119949b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
120049b5e25fSSatish Balay        0,
120149b5e25fSSatish Balay        0,
120249b5e25fSSatish Balay        0,
1203950f1e5bSHong Zhang        0,
120449b5e25fSSatish Balay        0,
120549b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
120649b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
120749b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
120849b5e25fSSatish Balay        0,
120949b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
121049b5e25fSSatish Balay        MatScale_SeqSBAIJ,
121149b5e25fSSatish Balay        0,
121249b5e25fSSatish Balay        0,
121349b5e25fSSatish Balay        0,
121449b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
121549b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
121649b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
121749b5e25fSSatish Balay        0,
121849b5e25fSSatish Balay        0,
121949b5e25fSSatish Balay        0,
122049b5e25fSSatish Balay        0,
122149b5e25fSSatish Balay        0,
122249b5e25fSSatish Balay        0,
122349b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
122449b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
122549b5e25fSSatish Balay        0,
122649b5e25fSSatish Balay        0,
1227d959ec07SHong Zhang        MatGetMaps_Petsc,
1228d959ec07SHong Zhang        0,
1229d959ec07SHong Zhang        0,
1230d959ec07SHong Zhang        0,
1231d959ec07SHong Zhang        0,
1232d959ec07SHong Zhang        0,
1233d959ec07SHong Zhang        0,
123424d5174aSHong Zhang        MatGetRowMax_SeqSBAIJ};
123549b5e25fSSatish Balay 
123649b5e25fSSatish Balay EXTERN_C_BEGIN
12374a2ae208SSatish Balay #undef __FUNCT__
12384a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
123949b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
124049b5e25fSSatish Balay {
12414afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1242c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
124349b5e25fSSatish Balay   int         ierr;
124449b5e25fSSatish Balay 
124549b5e25fSSatish Balay   PetscFunctionBegin;
124649b5e25fSSatish Balay   if (aij->nonew != 1) {
124729bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
124849b5e25fSSatish Balay   }
124949b5e25fSSatish Balay 
125049b5e25fSSatish Balay   /* allocate space for values if not already there */
125149b5e25fSSatish Balay   if (!aij->saved_values) {
1252eb7adc28SSatish Balay     ierr = PetscMalloc((nz+1)*sizeof(Scalar),&aij->saved_values);CHKERRQ(ierr);
125349b5e25fSSatish Balay   }
125449b5e25fSSatish Balay 
125549b5e25fSSatish Balay   /* copy values over */
125649b5e25fSSatish Balay   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
125749b5e25fSSatish Balay   PetscFunctionReturn(0);
125849b5e25fSSatish Balay }
125949b5e25fSSatish Balay EXTERN_C_END
126049b5e25fSSatish Balay 
126149b5e25fSSatish Balay EXTERN_C_BEGIN
12624a2ae208SSatish Balay #undef __FUNCT__
12634a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
126449b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
126549b5e25fSSatish Balay {
12664afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1267c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
126849b5e25fSSatish Balay 
126949b5e25fSSatish Balay   PetscFunctionBegin;
127049b5e25fSSatish Balay   if (aij->nonew != 1) {
127129bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
127249b5e25fSSatish Balay   }
127349b5e25fSSatish Balay   if (!aij->saved_values) {
127429bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
127549b5e25fSSatish Balay   }
127649b5e25fSSatish Balay 
127749b5e25fSSatish Balay   /* copy values over */
127849b5e25fSSatish Balay   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
127949b5e25fSSatish Balay   PetscFunctionReturn(0);
128049b5e25fSSatish Balay }
128149b5e25fSSatish Balay EXTERN_C_END
128249b5e25fSSatish Balay 
12838549e402SHong Zhang EXTERN_C_BEGIN
12844a2ae208SSatish Balay #undef __FUNCT__
12854a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ"
1286c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B)
1287c464158bSHong Zhang {
1288c464158bSHong Zhang   Mat_SeqSBAIJ *b;
12890c955e93SHong Zhang   int          ierr,size;
1290c464158bSHong Zhang 
1291c464158bSHong Zhang   PetscFunctionBegin;
1292c464158bSHong Zhang   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1293c464158bSHong Zhang   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1294273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1295273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1296c464158bSHong Zhang 
1297b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1298b0a32e0cSBarry Smith   B->data = (void*)b;
1299c464158bSHong Zhang   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1300c464158bSHong Zhang   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1301c464158bSHong Zhang   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1302c464158bSHong Zhang   B->ops->view        = MatView_SeqSBAIJ;
1303c464158bSHong Zhang   B->factor           = 0;
1304c464158bSHong Zhang   B->lupivotthreshold = 1.0;
1305c464158bSHong Zhang   B->mapping          = 0;
1306c464158bSHong Zhang   b->row              = 0;
1307c464158bSHong Zhang   b->icol             = 0;
1308c464158bSHong Zhang   b->reallocs         = 0;
1309c464158bSHong Zhang   b->saved_values     = 0;
1310c464158bSHong Zhang 
1311c464158bSHong Zhang   ierr = MapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1312c464158bSHong Zhang   ierr = MapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1313c464158bSHong Zhang 
1314c464158bSHong Zhang   b->sorted           = PETSC_FALSE;
1315c464158bSHong Zhang   b->roworiented      = PETSC_TRUE;
1316c464158bSHong Zhang   b->nonew            = 0;
1317c464158bSHong Zhang   b->diag             = 0;
1318c464158bSHong Zhang   b->solve_work       = 0;
1319c464158bSHong Zhang   b->mult_work        = 0;
1320c464158bSHong Zhang   b->spptr            = 0;
1321c464158bSHong Zhang   b->keepzeroedrows   = PETSC_FALSE;
1322c464158bSHong Zhang 
1323c464158bSHong Zhang   b->inew             = 0;
1324c464158bSHong Zhang   b->jnew             = 0;
1325c464158bSHong Zhang   b->anew             = 0;
1326c464158bSHong Zhang   b->a2anew           = 0;
1327c464158bSHong Zhang   b->permute          = PETSC_FALSE;
1328c464158bSHong Zhang 
1329c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1330c464158bSHong Zhang                                      "MatStoreValues_SeqSBAIJ",
1331c464158bSHong Zhang                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1332c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1333c464158bSHong Zhang                                      "MatRetrieveValues_SeqSBAIJ",
1334c464158bSHong Zhang                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1335c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1336c464158bSHong Zhang                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1337c464158bSHong Zhang                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1338c464158bSHong Zhang   PetscFunctionReturn(0);
1339c464158bSHong Zhang }
13408549e402SHong Zhang EXTERN_C_END
1341c464158bSHong Zhang 
13424a2ae208SSatish Balay #undef __FUNCT__
13434a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
134449b5e25fSSatish Balay /*@C
1345273d9f13SBarry Smith    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
134649b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
134749b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
134849b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
134949b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
135049b5e25fSSatish Balay 
1351c464158bSHong Zhang    Collective on Mat
135249b5e25fSSatish Balay 
135349b5e25fSSatish Balay    Input Parameters:
1354c464158bSHong Zhang +  A - the symmetric matrix
135549b5e25fSSatish Balay .  bs - size of block
135649b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
135749b5e25fSSatish Balay -  nnz - array containing the number of block nonzeros in the various block rows
135849b5e25fSSatish Balay          (possibly different for each block row) or PETSC_NULL
135949b5e25fSSatish Balay 
136049b5e25fSSatish Balay    Options Database Keys:
136149b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
136249b5e25fSSatish Balay                      block calculations (much slower)
136349b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
136449b5e25fSSatish Balay 
136549b5e25fSSatish Balay    Level: intermediate
136649b5e25fSSatish Balay 
136749b5e25fSSatish Balay    Notes:
1368c464158bSHong Zhang    The block AIJ format is compatible with standard Fortran 77
136949b5e25fSSatish Balay    storage.  That is, the stored row and column indices can begin at
137049b5e25fSSatish Balay    either one (as in Fortran) or zero.  See the users' manual for details.
137149b5e25fSSatish Balay 
137249b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
137349b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
137449b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
137549b5e25fSSatish Balay    matrices.
137649b5e25fSSatish Balay 
1377a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
137849b5e25fSSatish Balay @*/
1379273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
138049b5e25fSSatish Balay {
1381c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
13820c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
138349b5e25fSSatish Balay   PetscTruth   flg;
13844afc71dfSHong Zhang   int          s_nz;
138549b5e25fSSatish Balay 
138649b5e25fSSatish Balay   PetscFunctionBegin;
1387273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1388b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1389c464158bSHong Zhang   mbs  = B->m/bs;
139049b5e25fSSatish Balay   bs2  = bs*bs;
139149b5e25fSSatish Balay 
1392c464158bSHong Zhang   if (mbs*bs != B->m) {
139329bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
139449b5e25fSSatish Balay   }
139549b5e25fSSatish Balay 
1396435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1397435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
139849b5e25fSSatish Balay   if (nnz) {
139949b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
140029bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
140180fe4e49SBarry 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);
140249b5e25fSSatish Balay     }
140349b5e25fSSatish Balay   }
140449b5e25fSSatish Balay 
1405b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
140649b5e25fSSatish Balay   if (!flg) {
140749b5e25fSSatish Balay     switch (bs) {
140849b5e25fSSatish Balay     case 1:
14094ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
141049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
141149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
141249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
141349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
141449b5e25fSSatish Balay       break;
141549b5e25fSSatish Balay     case 2:
14164ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
141749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
141849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
141949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
142049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
142149b5e25fSSatish Balay       break;
142249b5e25fSSatish Balay     case 3:
14235f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
142449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
142549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
142649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
142749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
142849b5e25fSSatish Balay       break;
142949b5e25fSSatish Balay     case 4:
14305f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
143149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
143249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
143349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
143449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
143549b5e25fSSatish Balay       break;
143649b5e25fSSatish Balay     case 5:
14375f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
143849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
143949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
144049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
144149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
144249b5e25fSSatish Balay       break;
144349b5e25fSSatish Balay     case 6:
14445f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
144549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
144649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
144749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
144849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
144949b5e25fSSatish Balay       break;
145049b5e25fSSatish Balay     case 7:
1451de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
145249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
145349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1454de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
145549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
145649b5e25fSSatish Balay       break;
145749b5e25fSSatish Balay     }
145849b5e25fSSatish Balay   }
145949b5e25fSSatish Balay 
146049b5e25fSSatish Balay   b->mbs = mbs;
14614afc71dfSHong Zhang   b->nbs = mbs;
1462b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
146349b5e25fSSatish Balay   if (!nnz) {
1464435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
146549b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
146649b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
14678cef66ccSHong Zhang       b->imax[i] = nz;
146849b5e25fSSatish Balay     }
1469153ea458SHong Zhang     nz = nz*mbs; /* total nz */
147049b5e25fSSatish Balay   } else {
147149b5e25fSSatish Balay     nz = 0;
14728cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
147349b5e25fSSatish Balay   }
14748cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
14758cef66ccSHong Zhang   s_nz = nz;
147649b5e25fSSatish Balay 
147749b5e25fSSatish Balay   /* allocate the matrix space */
1478c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
147982502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
148049b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
148149b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
148249b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
148349b5e25fSSatish Balay   b->i = b->j + s_nz;
148449b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
148549b5e25fSSatish Balay 
148649b5e25fSSatish Balay   /* pointer to beginning of each row */
148749b5e25fSSatish Balay   b->i[0] = 0;
148849b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
148949b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
149049b5e25fSSatish Balay   }
149149b5e25fSSatish Balay 
149249b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
14935033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1494b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
149549b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
149649b5e25fSSatish Balay 
149749b5e25fSSatish Balay   b->bs               = bs;
149849b5e25fSSatish Balay   b->bs2              = bs2;
149949b5e25fSSatish Balay   b->s_nz             = 0;
150049b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1501153ea458SHong Zhang 
150216cdd363SHong Zhang   b->inew             = 0;
150316cdd363SHong Zhang   b->jnew             = 0;
150416cdd363SHong Zhang   b->anew             = 0;
150516cdd363SHong Zhang   b->a2anew           = 0;
15061a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1507c464158bSHong Zhang   PetscFunctionReturn(0);
1508c464158bSHong Zhang }
1509153ea458SHong Zhang 
151049b5e25fSSatish Balay 
15114a2ae208SSatish Balay #undef __FUNCT__
15124a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1513c464158bSHong Zhang /*@C
1514c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1515c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1516c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1517c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1518c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
151949b5e25fSSatish Balay 
1520c464158bSHong Zhang    Collective on MPI_Comm
1521c464158bSHong Zhang 
1522c464158bSHong Zhang    Input Parameters:
1523c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1524c464158bSHong Zhang .  bs - size of block
1525c464158bSHong Zhang .  m - number of rows, or number of columns
1526c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1527c464158bSHong Zhang -  nnz - array containing the number of block nonzeros in the various block rows
1528c464158bSHong Zhang          (possibly different for each block row) or PETSC_NULL
1529c464158bSHong Zhang 
1530c464158bSHong Zhang    Output Parameter:
1531c464158bSHong Zhang .  A - the symmetric matrix
1532c464158bSHong Zhang 
1533c464158bSHong Zhang    Options Database Keys:
1534c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1535c464158bSHong Zhang                      block calculations (much slower)
1536c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1537c464158bSHong Zhang 
1538c464158bSHong Zhang    Level: intermediate
1539c464158bSHong Zhang 
1540c464158bSHong Zhang    Notes:
1541c464158bSHong Zhang    The block AIJ format is compatible with standard Fortran 77
1542c464158bSHong Zhang    storage.  That is, the stored row and column indices can begin at
1543c464158bSHong Zhang    either one (as in Fortran) or zero.  See the users' manual for details.
1544c464158bSHong Zhang 
1545c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1546c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1547c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1548c464158bSHong Zhang    matrices.
1549c464158bSHong Zhang 
1550c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1551c464158bSHong Zhang @*/
1552c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1553c464158bSHong Zhang {
1554c464158bSHong Zhang   int ierr;
1555c464158bSHong Zhang 
1556c464158bSHong Zhang   PetscFunctionBegin;
1557c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1558c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1559273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
156049b5e25fSSatish Balay   PetscFunctionReturn(0);
156149b5e25fSSatish Balay }
156249b5e25fSSatish Balay 
15634a2ae208SSatish Balay #undef __FUNCT__
15644a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
156549b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
156649b5e25fSSatish Balay {
156749b5e25fSSatish Balay   Mat         C;
156849b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
156949b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
157049b5e25fSSatish Balay 
157149b5e25fSSatish Balay   PetscFunctionBegin;
157229bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
157349b5e25fSSatish Balay 
157449b5e25fSSatish Balay   *B = 0;
1575692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1576692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1577692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1578692f9cbeSHong Zhang 
157949b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1580273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
158149b5e25fSSatish Balay   C->factor         = A->factor;
158249b5e25fSSatish Balay   c->row            = 0;
158349b5e25fSSatish Balay   c->icol           = 0;
158449b5e25fSSatish Balay   c->saved_values   = 0;
158549b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
158649b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
158749b5e25fSSatish Balay 
158849b5e25fSSatish Balay   c->bs         = a->bs;
158949b5e25fSSatish Balay   c->bs2        = a->bs2;
159049b5e25fSSatish Balay   c->mbs        = a->mbs;
159149b5e25fSSatish Balay   c->nbs        = a->nbs;
159249b5e25fSSatish Balay 
1593b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1594b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
159549b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
159649b5e25fSSatish Balay     c->imax[i] = a->imax[i];
159749b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
159849b5e25fSSatish Balay   }
159949b5e25fSSatish Balay 
160049b5e25fSSatish Balay   /* allocate the matrix space */
160149b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
160249b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
160382502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
160449b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
160549b5e25fSSatish Balay   c->i = c->j + nz;
160649b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
160749b5e25fSSatish Balay   if (mbs > 0) {
160849b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
160949b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
161049b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
161149b5e25fSSatish Balay     } else {
161249b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
161349b5e25fSSatish Balay     }
161449b5e25fSSatish Balay   }
161549b5e25fSSatish Balay 
1616b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
161749b5e25fSSatish Balay   c->sorted      = a->sorted;
161849b5e25fSSatish Balay   c->roworiented = a->roworiented;
161949b5e25fSSatish Balay   c->nonew       = a->nonew;
162049b5e25fSSatish Balay 
162149b5e25fSSatish Balay   if (a->diag) {
1622b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1623b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
162449b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
162549b5e25fSSatish Balay       c->diag[i] = a->diag[i];
162649b5e25fSSatish Balay     }
162749b5e25fSSatish Balay   } else c->diag        = 0;
162849b5e25fSSatish Balay   c->s_nz               = a->s_nz;
162949b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
163049b5e25fSSatish Balay   c->solve_work         = 0;
163149b5e25fSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
163249b5e25fSSatish Balay   c->mult_work          = 0;
163349b5e25fSSatish Balay   *B = C;
1634b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
163549b5e25fSSatish Balay   PetscFunctionReturn(0);
163649b5e25fSSatish Balay }
163749b5e25fSSatish Balay 
16388549e402SHong Zhang EXTERN_C_BEGIN
16394a2ae208SSatish Balay #undef __FUNCT__
16404a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1641b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
164249b5e25fSSatish Balay {
164349b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
164449b5e25fSSatish Balay   Mat          B;
164549b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1646574b2666SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount;
164749b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
164849b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
164949b5e25fSSatish Balay   Scalar       *aa;
165049b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
165149b5e25fSSatish Balay 
165249b5e25fSSatish Balay   PetscFunctionBegin;
1653b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
165449b5e25fSSatish Balay   bs2  = bs*bs;
165549b5e25fSSatish Balay 
165649b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
165729bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1658b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
165949b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
166029bbc08cSBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
166149b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
166249b5e25fSSatish Balay 
166349b5e25fSSatish Balay   if (header[3] < 0) {
166429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
166549b5e25fSSatish Balay   }
166649b5e25fSSatish Balay 
166729bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
166849b5e25fSSatish Balay 
166949b5e25fSSatish Balay   /*
167049b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
167149b5e25fSSatish Balay     divisible by the blocksize
167249b5e25fSSatish Balay   */
167349b5e25fSSatish Balay   mbs        = M/bs;
167449b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
167549b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
167649b5e25fSSatish Balay   else                  mbs++;
167749b5e25fSSatish Balay   if (extra_rows) {
1678b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
167949b5e25fSSatish Balay   }
168049b5e25fSSatish Balay 
168149b5e25fSSatish Balay   /* read in row lengths */
1682b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
168349b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
168449b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
168549b5e25fSSatish Balay 
168649b5e25fSSatish Balay   /* read in column indices */
1687b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
168849b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
168949b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
169049b5e25fSSatish Balay 
169149b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
169282502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
169382502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1694574b2666SHong Zhang   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
169582502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
169649b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
169749b5e25fSSatish Balay   masked   = mask + mbs;
169849b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
169949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
170049b5e25fSSatish Balay     nmask = 0;
170149b5e25fSSatish Balay     for (j=0; j<bs; j++) {
170249b5e25fSSatish Balay       kmax = rowlengths[rowcount];
170349b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
17042d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
170503630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
170649b5e25fSSatish Balay       }
170749b5e25fSSatish Balay       rowcount++;
170849b5e25fSSatish Balay     }
1709574b2666SHong Zhang     s_browlengths[i] += nmask;
1710574b2666SHong Zhang     browlengths[i]    = 2*s_browlengths[i];
1711574b2666SHong Zhang 
171249b5e25fSSatish Balay     /* zero out the mask elements we set */
171349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
171449b5e25fSSatish Balay   }
171549b5e25fSSatish Balay 
171649b5e25fSSatish Balay   /* create our matrix */
171749b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
171849b5e25fSSatish Balay   B = *A;
171949b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
172049b5e25fSSatish Balay 
172149b5e25fSSatish Balay   /* set matrix "i" values */
172249b5e25fSSatish Balay   a->i[0] = 0;
172349b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1724574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1725574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
172649b5e25fSSatish Balay   }
172749b5e25fSSatish Balay   a->s_nz         = 0;
1728574b2666SHong Zhang   for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i];
172949b5e25fSSatish Balay 
173049b5e25fSSatish Balay   /* read in nonzero values */
1731b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(Scalar),&aa);CHKERRQ(ierr);
173249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
173349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
173449b5e25fSSatish Balay 
173549b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
173649b5e25fSSatish Balay   nzcount = 0; jcount = 0;
173749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
173849b5e25fSSatish Balay     nzcountb = nzcount;
173949b5e25fSSatish Balay     nmask    = 0;
174049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
174149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
174249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
17432d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
174403630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
17452d703238SHong Zhang       }
17462d703238SHong Zhang     }
17472d703238SHong Zhang     /* sort the masked values */
17482d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
17492d703238SHong Zhang 
17502d703238SHong Zhang     /* set "j" values into matrix */
17512d703238SHong Zhang     maskcount = 1;
17522d703238SHong Zhang     for (j=0; j<nmask; j++) {
175349b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
175449b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
175549b5e25fSSatish Balay     }
1756574b2666SHong Zhang 
175749b5e25fSSatish Balay     /* set "a" values into matrix */
175849b5e25fSSatish Balay     ishift = bs2*a->i[i];
175949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
176049b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
176149b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1762574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1763574b2666SHong Zhang         if (tmp >= i){
176449b5e25fSSatish Balay           block     = mask[tmp] - 1;
176549b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
176649b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1767574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1768574b2666SHong Zhang         }
1769574b2666SHong Zhang         nzcountb++;
177049b5e25fSSatish Balay       }
177149b5e25fSSatish Balay     }
177249b5e25fSSatish Balay     /* zero out the mask elements we set */
177349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
177449b5e25fSSatish Balay   }
177529bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
177649b5e25fSSatish Balay 
177749b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
177849b5e25fSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1779574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
178049b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
178149b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
178249b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
178349b5e25fSSatish Balay 
178449b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
178549b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
178649b5e25fSSatish Balay   PetscFunctionReturn(0);
178749b5e25fSSatish Balay }
17888549e402SHong Zhang EXTERN_C_END
1789574b2666SHong Zhang 
179049b5e25fSSatish Balay 
179149b5e25fSSatish Balay 
1792