xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 4e0d8c253c5baa1c48327ddb2b78bf7a43f58497)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
4a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h"         /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/inline/spops.h"
9aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1049b5e25fSSatish Balay 
1149b5e25fSSatish Balay #define CHUNKSIZE  10
1249b5e25fSSatish Balay 
1349b5e25fSSatish Balay /*
1449b5e25fSSatish Balay      Checks for missing diagonals
1549b5e25fSSatish Balay */
164a2ae208SSatish Balay #undef __FUNCT__
174a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
18dfbe8321SBarry Smith PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A)
1949b5e25fSSatish Balay {
20045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
216849ba73SBarry Smith   PetscErrorCode ierr;
2213f74950SBarry Smith   PetscInt       *diag,*jj = a->j,i;
2349b5e25fSSatish Balay 
2449b5e25fSSatish Balay   PetscFunctionBegin;
25045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2649b5e25fSSatish Balay   diag = a->diag;
2749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
2877431f27SBarry Smith     if (jj[diag[i]] != i) SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Matrix is missing diagonal number %D",i);
2949b5e25fSSatish Balay   }
3049b5e25fSSatish Balay   PetscFunctionReturn(0);
3149b5e25fSSatish Balay }
3249b5e25fSSatish Balay 
334a2ae208SSatish Balay #undef __FUNCT__
344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
35dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
3649b5e25fSSatish Balay {
37045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
386849ba73SBarry Smith   PetscErrorCode ierr;
3909f38230SBarry Smith   PetscInt       i;
4049b5e25fSSatish Balay 
4149b5e25fSSatish Balay   PetscFunctionBegin;
4209f38230SBarry Smith   if (!a->diag) {
4309f38230SBarry Smith     ierr = PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
4409f38230SBarry Smith   }
4509f38230SBarry Smith   for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i];
4649b5e25fSSatish Balay   PetscFunctionReturn(0);
4749b5e25fSSatish Balay }
4849b5e25fSSatish Balay 
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
518f7157efSSatish Balay static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
5249b5e25fSSatish Balay {
53a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
548f7157efSSatish Balay   PetscInt     i,j,n = a->mbs,nz = a->i[n],bs = A->rmap.bs;
558f7157efSSatish Balay   PetscErrorCode ierr;
5649b5e25fSSatish Balay 
5749b5e25fSSatish Balay   PetscFunctionBegin;
58d3e5a4abSHong Zhang   *nn = n;
59a1373b80SHong Zhang   if (!ia) PetscFunctionReturn(0);
608f7157efSSatish Balay   if (!blockcompressed) {
618f7157efSSatish Balay     /* malloc & create the natural set of indices */
62f1d0d59dSSatish Balay     ierr = PetscMalloc2((n+1)*bs,PetscInt,ia,nz*bs,PetscInt,ja);CHKERRQ(ierr);
638f7157efSSatish Balay     for (i=0; i<n+1; i++) {
648f7157efSSatish Balay       for (j=0; j<bs; j++) {
658f7157efSSatish Balay         *ia[i*bs+j] = a->i[i]*bs+j+oshift;
668f7157efSSatish Balay       }
678f7157efSSatish Balay     }
688f7157efSSatish Balay     for (i=0; i<nz; i++) {
698f7157efSSatish Balay       for (j=0; j<bs; j++) {
708f7157efSSatish Balay         *ja[i*bs+j] = a->j[i]*bs+j+oshift;
718f7157efSSatish Balay       }
728f7157efSSatish Balay     }
738f7157efSSatish Balay   } else { /* blockcompressed */
74a6ece127SHong Zhang     if (oshift == 1) {
75a6ece127SHong Zhang       /* temporarily add 1 to i and j indices */
766c6c5352SBarry Smith       for (i=0; i<nz; i++) a->j[i]++;
77a1373b80SHong Zhang       for (i=0; i<n+1; i++) a->i[i]++;
788f7157efSSatish Balay     }
79a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
80a6ece127SHong Zhang   }
818f7157efSSatish Balay 
8249b5e25fSSatish Balay   PetscFunctionReturn(0);
8349b5e25fSSatish Balay }
8449b5e25fSSatish Balay 
854a2ae208SSatish Balay #undef __FUNCT__
864a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
878f7157efSSatish Balay static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscTruth blockcompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
8849b5e25fSSatish Balay {
89b7aaefc3SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
908f7157efSSatish Balay   PetscInt     i,n = a->mbs,nz = a->i[n];
918f7157efSSatish Balay   PetscErrorCode ierr;
92a6ece127SHong Zhang 
9349b5e25fSSatish Balay   PetscFunctionBegin;
9449b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
95a6ece127SHong Zhang 
968f7157efSSatish Balay   if (!blockcompressed) {
978f7157efSSatish Balay     ierr = PetscFree2(*ia,*ja);CHKERRQ(ierr);
988f7157efSSatish Balay   } else if (oshift == 1) { /* blockcompressed */
996c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
100a6ece127SHong Zhang     for (i=0; i<n+1; i++) a->i[i]--;
101a6ece127SHong Zhang   }
1028f7157efSSatish Balay 
103a6ece127SHong Zhang   PetscFunctionReturn(0);
10449b5e25fSSatish Balay }
10549b5e25fSSatish Balay 
1064a2ae208SSatish Balay #undef __FUNCT__
1074a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
108dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
10949b5e25fSSatish Balay {
11049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
111dfbe8321SBarry Smith   PetscErrorCode ierr;
11249b5e25fSSatish Balay 
11349b5e25fSSatish Balay   PetscFunctionBegin;
114a9f03627SSatish Balay #if defined(PETSC_USE_LOG)
115899cda47SBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap.N,a->nz);
116a9f03627SSatish Balay #endif
117e6b907acSBarry Smith   ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr);
1189bfd6278SHong Zhang   if (a->row) {ierr = ISDestroy(a->row);CHKERRQ(ierr);}
1199bfd6278SHong Zhang   if (a->col){ierr = ISDestroy(a->col);CHKERRQ(ierr);}
1209bfd6278SHong Zhang   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
12105b42c5fSBarry Smith   ierr = PetscFree(a->diag);CHKERRQ(ierr);
12205b42c5fSBarry Smith   ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);
12305b42c5fSBarry Smith   ierr = PetscFree(a->solve_work);CHKERRQ(ierr);
124e2ee2a47SBarry Smith   ierr = PetscFree(a->relax_work);CHKERRQ(ierr);
12505b42c5fSBarry Smith   ierr = PetscFree(a->solves_work);CHKERRQ(ierr);
12605b42c5fSBarry Smith   ierr = PetscFree(a->mult_work);CHKERRQ(ierr);
12705b42c5fSBarry Smith   ierr = PetscFree(a->saved_values);CHKERRQ(ierr);
12805b42c5fSBarry Smith   ierr = PetscFree(a->xtoy);CHKERRQ(ierr);
1291a3463dfSHong Zhang 
1301a3463dfSHong Zhang   ierr = PetscFree(a->inew);CHKERRQ(ierr);
13149b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
132901853e0SKris Buschelman 
133dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
134901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr);
135901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr);
136901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr);
137901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
138901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr);
139901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
14049b5e25fSSatish Balay   PetscFunctionReturn(0);
14149b5e25fSSatish Balay }
14249b5e25fSSatish Balay 
1434a2ae208SSatish Balay #undef __FUNCT__
1444a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
145*4e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscTruth flg)
14649b5e25fSSatish Balay {
147045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
14863ba0a88SBarry Smith   PetscErrorCode ierr;
14949b5e25fSSatish Balay 
15049b5e25fSSatish Balay   PetscFunctionBegin;
1514d9d31abSKris Buschelman   switch (op) {
1524d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
153*4e0d8c25SBarry Smith     a->roworiented = flg;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
156*4e0d8c25SBarry Smith     a->sorted = flg;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
159*4e0d8c25SBarry Smith     a->keepzeroedrows = flg;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
162*4e0d8c25SBarry Smith     a->nonew = (flg ? 1 : 0);
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
165*4e0d8c25SBarry Smith     a->nonew = (flg ? -1 : 0);
1664d9d31abSKris Buschelman     break;
1674d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
168*4e0d8c25SBarry Smith     a->nonew = (flg ? -2 : 0);
1694d9d31abSKris Buschelman     break;
1704d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
171*4e0d8c25SBarry Smith   case MAT_NEW_DIAGONALS:
1724d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1734d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
174290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
1754d9d31abSKris Buschelman     break;
1769a4540c5SBarry Smith   case MAT_HERMITIAN:
177*4e0d8c25SBarry Smith     if (flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
17877e54ba9SKris Buschelman   case MAT_SYMMETRIC:
17977e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
1809a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
181*4e0d8c25SBarry Smith     if (!flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
182290bbb0aSBarry Smith     ierr = PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);CHKERRQ(ierr);
183290bbb0aSBarry Smith     break;
184941593c8SHong Zhang   case MAT_IGNORE_LOWER_TRIANGULAR:
185*4e0d8c25SBarry Smith     a->ignore_ltriangular = flg;
186941593c8SHong Zhang     break;
187941593c8SHong Zhang   case MAT_ERROR_LOWER_TRIANGULAR:
188*4e0d8c25SBarry Smith     a->ignore_ltriangular = flg;
18977e54ba9SKris Buschelman     break;
190f5edf698SHong Zhang   case MAT_GETROW_UPPERTRIANGULAR:
191*4e0d8c25SBarry Smith     a->getrow_utriangular = flg;
192f5edf698SHong Zhang     break;
1934d9d31abSKris Buschelman   default:
194ad86a440SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
19549b5e25fSSatish Balay   }
19649b5e25fSSatish Balay   PetscFunctionReturn(0);
19749b5e25fSSatish Balay }
19849b5e25fSSatish Balay 
1994a2ae208SSatish Balay #undef __FUNCT__
2004a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
20113f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v)
20249b5e25fSSatish Balay {
20349b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
2046849ba73SBarry Smith   PetscErrorCode ierr;
20513f74950SBarry Smith   PetscInt       itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
20649b5e25fSSatish Balay   MatScalar      *aa,*aa_i;
20787828ca2SBarry Smith   PetscScalar    *v_i;
20849b5e25fSSatish Balay 
20949b5e25fSSatish Balay   PetscFunctionBegin;
210*4e0d8c25SBarry Smith   if (A && !a->getrow_utriangular) SETERRQ(PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");
211f5edf698SHong Zhang   /* Get the upper triangular part of the row */
212899cda47SBarry Smith   bs  = A->rmap.bs;
21349b5e25fSSatish Balay   ai  = a->i;
21449b5e25fSSatish Balay   aj  = a->j;
21549b5e25fSSatish Balay   aa  = a->a;
21649b5e25fSSatish Balay   bs2 = a->bs2;
21749b5e25fSSatish Balay 
218899cda47SBarry Smith   if (row < 0 || row >= A->rmap.N) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row);
21949b5e25fSSatish Balay 
22049b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
22149b5e25fSSatish Balay   bp  = row % bs; /* Block position */
22249b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
22349b5e25fSSatish Balay   *ncols = bs*M;
22449b5e25fSSatish Balay 
22549b5e25fSSatish Balay   if (v) {
22649b5e25fSSatish Balay     *v = 0;
22749b5e25fSSatish Balay     if (*ncols) {
22887828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
22949b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23049b5e25fSSatish Balay         v_i  = *v + i*bs;
23149b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
23249b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
23349b5e25fSSatish Balay       }
23449b5e25fSSatish Balay     }
23549b5e25fSSatish Balay   }
23649b5e25fSSatish Balay 
23749b5e25fSSatish Balay   if (cols) {
23849b5e25fSSatish Balay     *cols = 0;
23949b5e25fSSatish Balay     if (*ncols) {
24013f74950SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr);
24149b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
24249b5e25fSSatish Balay         cols_i = *cols + i*bs;
24349b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24449b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24549b5e25fSSatish Balay       }
24649b5e25fSSatish Balay     }
24749b5e25fSSatish Balay   }
24849b5e25fSSatish Balay 
24949b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2505ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2515ddb2528SHong Zhang #ifdef column_search
25249b5e25fSSatish Balay   v_i    = *v    + M*bs;
25349b5e25fSSatish Balay   cols_i = *cols + M*bs;
25449b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25549b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25649b5e25fSSatish Balay     for (j=0; j<M; j++){
25749b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25849b5e25fSSatish Balay       if (itmp == bn){
25949b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
26049b5e25fSSatish Balay         for (k=0; k<bs; k++) {
26149b5e25fSSatish Balay           *cols_i++ = i*bs+k;
26249b5e25fSSatish Balay           *v_i++    = aa_i[k];
26349b5e25fSSatish Balay         }
26449b5e25fSSatish Balay         *ncols += bs;
26549b5e25fSSatish Balay         break;
26649b5e25fSSatish Balay       }
26749b5e25fSSatish Balay     }
26849b5e25fSSatish Balay   }
2695ddb2528SHong Zhang #endif
27049b5e25fSSatish Balay   PetscFunctionReturn(0);
27149b5e25fSSatish Balay }
27249b5e25fSSatish Balay 
2734a2ae208SSatish Balay #undef __FUNCT__
2744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
27513f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
27649b5e25fSSatish Balay {
277dfbe8321SBarry Smith   PetscErrorCode ierr;
27849b5e25fSSatish Balay 
27949b5e25fSSatish Balay   PetscFunctionBegin;
28005b42c5fSBarry Smith   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
28105b42c5fSBarry Smith   if (v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}
28249b5e25fSSatish Balay   PetscFunctionReturn(0);
28349b5e25fSSatish Balay }
28449b5e25fSSatish Balay 
2854a2ae208SSatish Balay #undef __FUNCT__
286f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ"
287f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
288f5edf698SHong Zhang {
289f5edf698SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
290f5edf698SHong Zhang 
291f5edf698SHong Zhang   PetscFunctionBegin;
292f5edf698SHong Zhang   a->getrow_utriangular = PETSC_TRUE;
293f5edf698SHong Zhang   PetscFunctionReturn(0);
294f5edf698SHong Zhang }
295f5edf698SHong Zhang #undef __FUNCT__
296f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ"
297f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
298f5edf698SHong Zhang {
299f5edf698SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
300f5edf698SHong Zhang 
301f5edf698SHong Zhang   PetscFunctionBegin;
302f5edf698SHong Zhang   a->getrow_utriangular = PETSC_FALSE;
303f5edf698SHong Zhang   PetscFunctionReturn(0);
304f5edf698SHong Zhang }
305f5edf698SHong Zhang 
306f5edf698SHong Zhang #undef __FUNCT__
3074a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
308dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B)
30949b5e25fSSatish Balay {
310dfbe8321SBarry Smith   PetscErrorCode ierr;
31149b5e25fSSatish Balay   PetscFunctionBegin;
312999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
3138115998fSBarry Smith   PetscFunctionReturn(0);
31449b5e25fSSatish Balay }
31549b5e25fSSatish Balay 
3164a2ae208SSatish Balay #undef __FUNCT__
3174a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
3186849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
31949b5e25fSSatish Balay {
32049b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
321dfbe8321SBarry Smith   PetscErrorCode    ierr;
322899cda47SBarry Smith   PetscInt          i,j,bs = A->rmap.bs,k,l,bs2=a->bs2;
3232dcb1b2aSMatthew Knepley   const char        *name;
324f3ef73ceSBarry Smith   PetscViewerFormat format;
32549b5e25fSSatish Balay 
32649b5e25fSSatish Balay   PetscFunctionBegin;
32780fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
328b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
329456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
33077431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);CHKERRQ(ierr);
331fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
33270d5e725SHong Zhang     if (A->factor && bs>1){
33370d5e725SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr);
33470d5e725SHong Zhang       PetscFunctionReturn(0);
33570d5e725SHong Zhang     }
336c9f458caSMatthew Knepley     Mat aij;
337c9f458caSMatthew Knepley     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr);
338c9f458caSMatthew Knepley     ierr = MatView(aij,viewer);CHKERRQ(ierr);
339c9f458caSMatthew Knepley     ierr = MatDestroy(aij);CHKERRQ(ierr);
340fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
341b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
34249b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
34349b5e25fSSatish Balay       for (j=0; j<bs; j++) {
34477431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
34549b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
34649b5e25fSSatish Balay           for (l=0; l<bs; l++) {
34749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
34849b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
349a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l,
35049b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35149b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
352a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l,
35349b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35449b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
355a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35649b5e25fSSatish Balay             }
35749b5e25fSSatish Balay #else
35849b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
359a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
36049b5e25fSSatish Balay             }
36149b5e25fSSatish Balay #endif
36249b5e25fSSatish Balay           }
36349b5e25fSSatish Balay         }
364b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
36549b5e25fSSatish Balay       }
36649b5e25fSSatish Balay     }
367b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
368c1490034SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
369c1490034SHong Zhang      PetscFunctionReturn(0);
37049b5e25fSSatish Balay   } else {
3718608aa04SHong Zhang     if (A->factor && bs>1){
3728608aa04SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored. MatView_SeqSBAIJ_ASCII() may not display complete or logically correct entries!\n");CHKERRQ(ierr);
3738608aa04SHong Zhang     }
374b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
37549b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
37649b5e25fSSatish Balay       for (j=0; j<bs; j++) {
37777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
37849b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
37949b5e25fSSatish Balay           for (l=0; l<bs; l++) {
38049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
38149b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
382a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l,
38349b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38449b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
385a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l,
38649b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38749b5e25fSSatish Balay             } else {
388a83599f4SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38949b5e25fSSatish Balay             }
39049b5e25fSSatish Balay #else
391e9f7bc9eSHong Zhang             ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
39249b5e25fSSatish Balay #endif
39349b5e25fSSatish Balay           }
39449b5e25fSSatish Balay         }
395b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39649b5e25fSSatish Balay       }
39749b5e25fSSatish Balay     }
398b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
39949b5e25fSSatish Balay   }
400b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
40149b5e25fSSatish Balay   PetscFunctionReturn(0);
40249b5e25fSSatish Balay }
40349b5e25fSSatish Balay 
4044a2ae208SSatish Balay #undef __FUNCT__
4054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
4066849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
40749b5e25fSSatish Balay {
40849b5e25fSSatish Balay   Mat            A = (Mat) Aa;
40949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
4106849ba73SBarry Smith   PetscErrorCode ierr;
411899cda47SBarry Smith   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap.bs,bs2=a->bs2;
41213f74950SBarry Smith   PetscMPIInt    rank;
41349b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
41449b5e25fSSatish Balay   MatScalar      *aa;
41549b5e25fSSatish Balay   MPI_Comm       comm;
416b0a32e0cSBarry Smith   PetscViewer    viewer;
41749b5e25fSSatish Balay 
41849b5e25fSSatish Balay   PetscFunctionBegin;
41949b5e25fSSatish Balay   /*
42049b5e25fSSatish Balay     This is nasty. If this is called from an originally parallel matrix
42149b5e25fSSatish Balay     then all processes call this,but only the first has the matrix so the
42249b5e25fSSatish Balay     rest should return immediately.
42349b5e25fSSatish Balay   */
42449b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
42549b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
42649b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
42749b5e25fSSatish Balay 
42849b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
42949b5e25fSSatish Balay 
430b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
431b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
43249b5e25fSSatish Balay 
43349b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
434b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
43549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
43649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
437899cda47SBarry Smith       y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0;
43849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
43949b5e25fSSatish Balay       aa = a->a + j*bs2;
44049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
44149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
44249b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
443b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
44449b5e25fSSatish Balay         }
44549b5e25fSSatish Balay       }
44649b5e25fSSatish Balay     }
44749b5e25fSSatish Balay   }
448b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
44949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
451899cda47SBarry Smith       y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0;
45249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
45349b5e25fSSatish Balay       aa = a->a + j*bs2;
45449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
45549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
45649b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
457b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
45849b5e25fSSatish Balay         }
45949b5e25fSSatish Balay       }
46049b5e25fSSatish Balay     }
46149b5e25fSSatish Balay   }
46249b5e25fSSatish Balay 
463b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
46449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
466899cda47SBarry Smith       y_l = A->rmap.N - row - 1.0; y_r = y_l + 1.0;
46749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46849b5e25fSSatish Balay       aa = a->a + j*bs2;
46949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47149b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
472b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47349b5e25fSSatish Balay         }
47449b5e25fSSatish Balay       }
47549b5e25fSSatish Balay     }
47649b5e25fSSatish Balay   }
47749b5e25fSSatish Balay   PetscFunctionReturn(0);
47849b5e25fSSatish Balay }
47949b5e25fSSatish Balay 
4804a2ae208SSatish Balay #undef __FUNCT__
4814a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
4826849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
48349b5e25fSSatish Balay {
484dfbe8321SBarry Smith   PetscErrorCode ierr;
48549b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,w,h;
486b0a32e0cSBarry Smith   PetscDraw      draw;
48749b5e25fSSatish Balay   PetscTruth     isnull;
48849b5e25fSSatish Balay 
48949b5e25fSSatish Balay   PetscFunctionBegin;
490b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
491b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
49249b5e25fSSatish Balay 
49349b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
494899cda47SBarry Smith   xr  = A->rmap.N; yr = A->rmap.N; h = yr/10.0; w = xr/10.0;
49549b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
496b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
497b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
49849b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
49949b5e25fSSatish Balay   PetscFunctionReturn(0);
50049b5e25fSSatish Balay }
50149b5e25fSSatish Balay 
5024a2ae208SSatish Balay #undef __FUNCT__
5034a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
504dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
50549b5e25fSSatish Balay {
506dfbe8321SBarry Smith   PetscErrorCode ierr;
50732077d6dSBarry Smith   PetscTruth     iascii,isdraw;
50849b5e25fSSatish Balay 
50949b5e25fSSatish Balay   PetscFunctionBegin;
51032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
511fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
51232077d6dSBarry Smith   if (iascii){
51349b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
51449b5e25fSSatish Balay   } else if (isdraw) {
51549b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
51649b5e25fSSatish Balay   } else {
517a5e6ed63SBarry Smith     Mat B;
518ceb03754SKris Buschelman     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
519a5e6ed63SBarry Smith     ierr = MatView(B,viewer);CHKERRQ(ierr);
520a5e6ed63SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
52149b5e25fSSatish Balay   }
52249b5e25fSSatish Balay   PetscFunctionReturn(0);
52349b5e25fSSatish Balay }
52449b5e25fSSatish Balay 
52549b5e25fSSatish Balay 
5264a2ae208SSatish Balay #undef __FUNCT__
5274a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
52813f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
52949b5e25fSSatish Balay {
530045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
53113f74950SBarry Smith   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
53213f74950SBarry Smith   PetscInt     *ai = a->i,*ailen = a->ilen;
533899cda47SBarry Smith   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap.bs,bs2=a->bs2;
53497e567efSBarry Smith   MatScalar    *ap,*aa = a->a;
53549b5e25fSSatish Balay 
53649b5e25fSSatish Balay   PetscFunctionBegin;
53749b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
53849b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
53997e567efSBarry Smith     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
540899cda47SBarry Smith     if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1);
54149b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
54249b5e25fSSatish Balay     nrow = ailen[brow];
54349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
54497e567efSBarry Smith       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
545899cda47SBarry Smith       if (in[l] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1);
54649b5e25fSSatish Balay       col  = in[l] ;
54749b5e25fSSatish Balay       bcol = col/bs;
54849b5e25fSSatish Balay       cidx = col%bs;
54949b5e25fSSatish Balay       ridx = row%bs;
55049b5e25fSSatish Balay       high = nrow;
55149b5e25fSSatish Balay       low  = 0; /* assume unsorted */
55249b5e25fSSatish Balay       while (high-low > 5) {
55349b5e25fSSatish Balay         t = (low+high)/2;
55449b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
55549b5e25fSSatish Balay         else             low  = t;
55649b5e25fSSatish Balay       }
55749b5e25fSSatish Balay       for (i=low; i<high; i++) {
55849b5e25fSSatish Balay         if (rp[i] > bcol) break;
55949b5e25fSSatish Balay         if (rp[i] == bcol) {
56049b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
56149b5e25fSSatish Balay           goto finished;
56249b5e25fSSatish Balay         }
56349b5e25fSSatish Balay       }
56497e567efSBarry Smith       *v++ = 0.0;
56549b5e25fSSatish Balay       finished:;
56649b5e25fSSatish Balay     }
56749b5e25fSSatish Balay   }
56849b5e25fSSatish Balay   PetscFunctionReturn(0);
56949b5e25fSSatish Balay }
57049b5e25fSSatish Balay 
57149b5e25fSSatish Balay 
5724a2ae208SSatish Balay #undef __FUNCT__
5734a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
57413f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
57549b5e25fSSatish Balay {
5760880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
5776849ba73SBarry Smith   PetscErrorCode  ierr;
578e2ee6c50SBarry Smith   PetscInt        *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
57913f74950SBarry Smith   PetscInt        *imax=a->imax,*ai=a->i,*ailen=a->ilen;
580899cda47SBarry Smith   PetscInt        *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap.bs,stepval;
5810880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
582f15d580aSBarry Smith   const MatScalar *value = v;
583f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
5840880e062SHong Zhang 
58549b5e25fSSatish Balay   PetscFunctionBegin;
5860880e062SHong Zhang   if (roworiented) {
5870880e062SHong Zhang     stepval = (n-1)*bs;
5880880e062SHong Zhang   } else {
5890880e062SHong Zhang     stepval = (m-1)*bs;
5900880e062SHong Zhang   }
5910880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
5920880e062SHong Zhang     row  = im[k];
5930880e062SHong Zhang     if (row < 0) continue;
5942515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
59577431f27SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
5960880e062SHong Zhang #endif
5970880e062SHong Zhang     rp   = aj + ai[row];
5980880e062SHong Zhang     ap   = aa + bs2*ai[row];
5990880e062SHong Zhang     rmax = imax[row];
6000880e062SHong Zhang     nrow = ailen[row];
6010880e062SHong Zhang     low  = 0;
602818f2c47SBarry Smith     high = nrow;
6030880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
6040880e062SHong Zhang       if (in[l] < 0) continue;
6050880e062SHong Zhang       col = in[l];
6062515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
60777431f27SBarry Smith       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
608b1823623SSatish Balay #endif
6090880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
6100880e062SHong Zhang       if (roworiented) {
6110880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
6120880e062SHong Zhang       } else {
6130880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
6140880e062SHong Zhang       }
6157cd84e04SBarry Smith       if (col <= lastcol) low = 0; else high = nrow;
616e2ee6c50SBarry Smith       lastcol = col;
6170880e062SHong Zhang       while (high-low > 7) {
6180880e062SHong Zhang         t = (low+high)/2;
6190880e062SHong Zhang         if (rp[t] > col) high = t;
6200880e062SHong Zhang         else             low  = t;
6210880e062SHong Zhang       }
6220880e062SHong Zhang       for (i=low; i<high; i++) {
6230880e062SHong Zhang         if (rp[i] > col) break;
6240880e062SHong Zhang         if (rp[i] == col) {
6250880e062SHong Zhang           bap  = ap +  bs2*i;
6260880e062SHong Zhang           if (roworiented) {
6270880e062SHong Zhang             if (is == ADD_VALUES) {
6280880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6290880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6300880e062SHong Zhang                   bap[jj] += *value++;
6310880e062SHong Zhang                 }
6320880e062SHong Zhang               }
6330880e062SHong Zhang             } else {
6340880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6350880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6360880e062SHong Zhang                   bap[jj] = *value++;
6370880e062SHong Zhang                 }
6380880e062SHong Zhang                }
6390880e062SHong Zhang             }
6400880e062SHong Zhang           } else {
6410880e062SHong Zhang             if (is == ADD_VALUES) {
6420880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6430880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6440880e062SHong Zhang                   *bap++ += *value++;
6450880e062SHong Zhang                 }
6460880e062SHong Zhang               }
6470880e062SHong Zhang             } else {
6480880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6490880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6500880e062SHong Zhang                   *bap++  = *value++;
6510880e062SHong Zhang                 }
6520880e062SHong Zhang               }
6530880e062SHong Zhang             }
6540880e062SHong Zhang           }
6550880e062SHong Zhang           goto noinsert2;
6560880e062SHong Zhang         }
6570880e062SHong Zhang       }
6580880e062SHong Zhang       if (nonew == 1) goto noinsert2;
659085a36d4SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
660421e10b8SBarry Smith       MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
661c03d1d03SSatish Balay       N = nrow++ - 1; high++;
6620880e062SHong Zhang       /* shift up all the later entries in this row */
6630880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
6640880e062SHong Zhang         rp[ii+1] = rp[ii];
6650880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
6660880e062SHong Zhang       }
6670880e062SHong Zhang       if (N >= i) {
6680880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6690880e062SHong Zhang       }
6700880e062SHong Zhang       rp[i] = col;
6710880e062SHong Zhang       bap   = ap +  bs2*i;
6720880e062SHong Zhang       if (roworiented) {
6730880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6740880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
6750880e062SHong Zhang             bap[jj] = *value++;
6760880e062SHong Zhang           }
6770880e062SHong Zhang         }
6780880e062SHong Zhang       } else {
6790880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6800880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
6810880e062SHong Zhang             *bap++  = *value++;
6820880e062SHong Zhang           }
6830880e062SHong Zhang         }
6840880e062SHong Zhang        }
6850880e062SHong Zhang     noinsert2:;
6860880e062SHong Zhang       low = i;
6870880e062SHong Zhang     }
6880880e062SHong Zhang     ailen[row] = nrow;
6890880e062SHong Zhang   }
6900880e062SHong Zhang    PetscFunctionReturn(0);
69149b5e25fSSatish Balay }
69249b5e25fSSatish Balay 
6934a2ae208SSatish Balay #undef __FUNCT__
6944a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
695dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
69649b5e25fSSatish Balay {
69749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
6986849ba73SBarry Smith   PetscErrorCode ierr;
69913f74950SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
700899cda47SBarry Smith   PetscInt       m = A->rmap.N,*ip,N,*ailen = a->ilen;
70113f74950SBarry Smith   PetscInt       mbs = a->mbs,bs2 = a->bs2,rmax = 0;
70249b5e25fSSatish Balay   MatScalar      *aa = a->a,*ap;
70349b5e25fSSatish Balay 
70449b5e25fSSatish Balay   PetscFunctionBegin;
70549b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
70649b5e25fSSatish Balay 
70749b5e25fSSatish Balay   if (m) rmax = ailen[0];
70849b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
70949b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
71049b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
71149b5e25fSSatish Balay      rmax   = PetscMax(rmax,ailen[i]);
71249b5e25fSSatish Balay      if (fshift) {
71349b5e25fSSatish Balay        ip = aj + ai[i]; ap = aa + bs2*ai[i];
71449b5e25fSSatish Balay        N = ailen[i];
71549b5e25fSSatish Balay        for (j=0; j<N; j++) {
71649b5e25fSSatish Balay          ip[j-fshift] = ip[j];
71749b5e25fSSatish Balay          ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
71849b5e25fSSatish Balay        }
71949b5e25fSSatish Balay      }
72049b5e25fSSatish Balay      ai[i] = ai[i-1] + ailen[i-1];
72149b5e25fSSatish Balay   }
72249b5e25fSSatish Balay   if (mbs) {
72349b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
72449b5e25fSSatish Balay      ai[mbs] = ai[mbs-1] + ailen[mbs-1];
72549b5e25fSSatish Balay   }
72649b5e25fSSatish Balay   /* reset ilen and imax for each row */
72749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
72849b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
72949b5e25fSSatish Balay   }
7306c6c5352SBarry Smith   a->nz = ai[mbs];
73149b5e25fSSatish Balay 
732b424e231SHong Zhang   /* diagonals may have moved, reset it */
733b424e231SHong Zhang   if (a->diag) {
73413f74950SBarry Smith     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
73549b5e25fSSatish Balay   }
736899cda47SBarry Smith   ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap.N,A->rmap.bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr);
737ae15b995SBarry Smith   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr);
738ae15b995SBarry Smith   ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr);
73949b5e25fSSatish Balay   a->reallocs          = 0;
74049b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
74149b5e25fSSatish Balay   PetscFunctionReturn(0);
74249b5e25fSSatish Balay }
74349b5e25fSSatish Balay 
74449b5e25fSSatish Balay /*
74549b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
74649b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
74749b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
74849b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
74949b5e25fSSatish Balay */
7504a2ae208SSatish Balay #undef __FUNCT__
7514a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
75213f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
75349b5e25fSSatish Balay {
75413f74950SBarry Smith   PetscInt   i,j,k,row;
75549b5e25fSSatish Balay   PetscTruth flg;
75649b5e25fSSatish Balay 
75749b5e25fSSatish Balay   PetscFunctionBegin;
75849b5e25fSSatish Balay    for (i=0,j=0; i<n; j++) {
75949b5e25fSSatish Balay      row = idx[i];
76049b5e25fSSatish Balay      if (row%bs!=0) { /* Not the begining of a block */
76149b5e25fSSatish Balay        sizes[j] = 1;
76249b5e25fSSatish Balay        i++;
76349b5e25fSSatish Balay      } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
76449b5e25fSSatish Balay        sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
76549b5e25fSSatish Balay        i++;
76649b5e25fSSatish Balay      } else { /* Begining of the block, so check if the complete block exists */
76749b5e25fSSatish Balay        flg = PETSC_TRUE;
76849b5e25fSSatish Balay        for (k=1; k<bs; k++) {
76949b5e25fSSatish Balay          if (row+k != idx[i+k]) { /* break in the block */
77049b5e25fSSatish Balay            flg = PETSC_FALSE;
77149b5e25fSSatish Balay            break;
77249b5e25fSSatish Balay          }
77349b5e25fSSatish Balay        }
774abc0a331SBarry Smith        if (flg) { /* No break in the bs */
77549b5e25fSSatish Balay          sizes[j] = bs;
77649b5e25fSSatish Balay          i+= bs;
77749b5e25fSSatish Balay        } else {
77849b5e25fSSatish Balay          sizes[j] = 1;
77949b5e25fSSatish Balay          i++;
78049b5e25fSSatish Balay        }
78149b5e25fSSatish Balay      }
78249b5e25fSSatish Balay    }
78349b5e25fSSatish Balay    *bs_max = j;
78449b5e25fSSatish Balay    PetscFunctionReturn(0);
78549b5e25fSSatish Balay }
78649b5e25fSSatish Balay 
78749b5e25fSSatish Balay 
78849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
78949b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
79049b5e25fSSatish Balay */
79149b5e25fSSatish Balay 
7924a2ae208SSatish Balay #undef __FUNCT__
7934a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
79413f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
79549b5e25fSSatish Balay {
79649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
7976849ba73SBarry Smith   PetscErrorCode ierr;
798e2ee6c50SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
79913f74950SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
800899cda47SBarry Smith   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->rmap.bs,brow,bcol;
80113f74950SBarry Smith   PetscInt       ridx,cidx,bs2=a->bs2;
80249b5e25fSSatish Balay   MatScalar      *ap,value,*aa=a->a,*bap;
80349b5e25fSSatish Balay 
80449b5e25fSSatish Balay   PetscFunctionBegin;
80549b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
80649b5e25fSSatish Balay     row  = im[k];       /* row number */
80749b5e25fSSatish Balay     brow = row/bs;      /* block row number */
80849b5e25fSSatish Balay     if (row < 0) continue;
8092515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
810899cda47SBarry Smith     if (row >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.N-1);
81149b5e25fSSatish Balay #endif
81249b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
81349b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
81449b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
81549b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
81649b5e25fSSatish Balay     low  = 0;
81749b5e25fSSatish Balay 
81849b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
81949b5e25fSSatish Balay       if (in[l] < 0) continue;
8202515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
821899cda47SBarry Smith       if (in[l] >= A->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap.N-1);
82249b5e25fSSatish Balay #endif
82349b5e25fSSatish Balay       col = in[l];
82449b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
82549b5e25fSSatish Balay 
826941593c8SHong Zhang       if (brow > bcol) {
827941593c8SHong Zhang         if (a->ignore_ltriangular){
828941593c8SHong Zhang           continue; /* ignore lower triangular values */
829941593c8SHong Zhang         } else {
830*4e0d8c25SBarry Smith           SETERRQ(PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
831941593c8SHong Zhang         }
832941593c8SHong Zhang       }
833f4989cb3SHong Zhang 
83449b5e25fSSatish Balay       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8358549e402SHong Zhang       if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
83649b5e25fSSatish Balay         /* element value a(k,l) */
83749b5e25fSSatish Balay         if (roworiented) {
83849b5e25fSSatish Balay           value = v[l + k*n];
83949b5e25fSSatish Balay         } else {
84049b5e25fSSatish Balay           value = v[k + l*m];
84149b5e25fSSatish Balay         }
84249b5e25fSSatish Balay 
84349b5e25fSSatish Balay         /* move pointer bap to a(k,l) quickly and add/insert value */
8447cd84e04SBarry Smith         if (col <= lastcol) low = 0; high = nrow;
845e2ee6c50SBarry Smith         lastcol = col;
84649b5e25fSSatish Balay         while (high-low > 7) {
84749b5e25fSSatish Balay           t = (low+high)/2;
84849b5e25fSSatish Balay           if (rp[t] > bcol) high = t;
84949b5e25fSSatish Balay           else              low  = t;
85049b5e25fSSatish Balay         }
85149b5e25fSSatish Balay         for (i=low; i<high; i++) {
85249b5e25fSSatish Balay           if (rp[i] > bcol) break;
85349b5e25fSSatish Balay           if (rp[i] == bcol) {
85449b5e25fSSatish Balay             bap  = ap +  bs2*i + bs*cidx + ridx;
85549b5e25fSSatish Balay             if (is == ADD_VALUES) *bap += value;
85649b5e25fSSatish Balay             else                  *bap  = value;
8578549e402SHong Zhang             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8588549e402SHong Zhang             if (brow == bcol && ridx < cidx){
8598549e402SHong Zhang               bap  = ap +  bs2*i + bs*ridx + cidx;
8608549e402SHong Zhang               if (is == ADD_VALUES) *bap += value;
8618549e402SHong Zhang               else                  *bap  = value;
8628549e402SHong Zhang             }
86349b5e25fSSatish Balay             goto noinsert1;
86449b5e25fSSatish Balay           }
86549b5e25fSSatish Balay         }
86649b5e25fSSatish Balay 
86749b5e25fSSatish Balay         if (nonew == 1) goto noinsert1;
868085a36d4SBarry Smith         if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
869421e10b8SBarry Smith         MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
87049b5e25fSSatish Balay 
871c03d1d03SSatish Balay         N = nrow++ - 1; high++;
87249b5e25fSSatish Balay         /* shift up all the later entries in this row */
87349b5e25fSSatish Balay         for (ii=N; ii>=i; ii--) {
87449b5e25fSSatish Balay           rp[ii+1] = rp[ii];
87549b5e25fSSatish Balay           ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
87649b5e25fSSatish Balay         }
87749b5e25fSSatish Balay         if (N>=i) {
87849b5e25fSSatish Balay           ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
87949b5e25fSSatish Balay         }
88049b5e25fSSatish Balay         rp[i]                      = bcol;
88149b5e25fSSatish Balay         ap[bs2*i + bs*cidx + ridx] = value;
88249b5e25fSSatish Balay       noinsert1:;
88349b5e25fSSatish Balay         low = i;
8848549e402SHong Zhang       }
88549b5e25fSSatish Balay     }   /* end of loop over added columns */
88649b5e25fSSatish Balay     ailen[brow] = nrow;
88749b5e25fSSatish Balay   }   /* end of loop over added rows */
88849b5e25fSSatish Balay   PetscFunctionReturn(0);
88949b5e25fSSatish Balay }
89049b5e25fSSatish Balay 
8914a2ae208SSatish Balay #undef __FUNCT__
8924d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
893dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
89449b5e25fSSatish Balay {
8954ccecd49SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
89649b5e25fSSatish Balay   Mat            outA;
897dfbe8321SBarry Smith   PetscErrorCode ierr;
898c84f5b01SHong Zhang   PetscTruth     row_identity;
89949b5e25fSSatish Balay 
90049b5e25fSSatish Balay   PetscFunctionBegin;
901c84f5b01SHong Zhang   if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
902c84f5b01SHong Zhang   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
903c84f5b01SHong Zhang   if (!row_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
904c84f5b01SHong Zhang   if (inA->rmap.bs != 1) SETERRQ1(PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap.bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */
905c84f5b01SHong Zhang 
90649b5e25fSSatish Balay   outA        = inA;
9071260e1f8SHong Zhang   inA->factor = FACTOR_CHOLESKY;
90849b5e25fSSatish Balay 
9091a3463dfSHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
91049b5e25fSSatish Balay   /*
911c84f5b01SHong Zhang     Blocksize < 8 have a special faster factorization/solver
912c84f5b01SHong Zhang     for ICC(0) factorization with natural ordering
91349b5e25fSSatish Balay   */
914c84f5b01SHong Zhang   switch (inA->rmap.bs){ /* Note: row_identity = PETSC_TRUE! */
91549b5e25fSSatish Balay   case 1:
916c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
9170fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
91807f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
919d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
920759322afSHong Zhang     inA->ops->forwardsolve     = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
921759322afSHong Zhang     inA->ops->backwardsolve    = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
922c84f5b01SHong Zhang     ierr = PetscInfo(inA,"Using special in-place natural ordering solvetrans BS=1\n");CHKERRQ(ierr);
923c84f5b01SHong Zhang     break;
92449b5e25fSSatish Balay   case 2:
925c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
9261a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
9271a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
928759322afSHong Zhang     inA->ops->forwardsolve     = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering;
929759322afSHong Zhang     inA->ops->backwardsolve    = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering;
930ae15b995SBarry Smith     ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr);
93149b5e25fSSatish Balay     break;
93249b5e25fSSatish Balay   case 3:
933c84f5b01SHong Zhang      inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
9341a3463dfSHong Zhang      inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
9351a3463dfSHong Zhang      inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
936759322afSHong Zhang      inA->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering;
937759322afSHong Zhang      inA->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering;
938ae15b995SBarry Smith      ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr);
93949b5e25fSSatish Balay      break;
94049b5e25fSSatish Balay   case 4:
941c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
9421a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
9431a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
944759322afSHong Zhang     inA->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering;
945759322afSHong Zhang     inA->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering;
946ae15b995SBarry Smith     ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr);
94749b5e25fSSatish Balay     break;
94849b5e25fSSatish Balay   case 5:
949c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
9501a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
9511a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
952759322afSHong Zhang     inA->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering;
953759322afSHong Zhang     inA->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering;
954ae15b995SBarry Smith     ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr);
95549b5e25fSSatish Balay     break;
95649b5e25fSSatish Balay   case 6:
957c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
9581a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
9591a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
960759322afSHong Zhang     inA->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering;
961759322afSHong Zhang     inA->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering;
962ae15b995SBarry Smith     ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr);
96349b5e25fSSatish Balay     break;
96449b5e25fSSatish Balay   case 7:
965c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
9661a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
9671a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
968759322afSHong Zhang     inA->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering;
969759322afSHong Zhang     inA->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering;
970ae15b995SBarry Smith     ierr = PetscInfo(inA,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr);
97149b5e25fSSatish Balay     break;
97249b5e25fSSatish Balay   default:
973c84f5b01SHong Zhang     inA->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
974c84f5b01SHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_N_NaturalOrdering;
975c84f5b01SHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
976759322afSHong Zhang     inA->ops->forwardsolve    = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering;
977759322afSHong Zhang     inA->ops->backwardsolve   = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering;
978c84f5b01SHong Zhang     break;
979c84f5b01SHong Zhang   }
98049b5e25fSSatish Balay 
981c3122656SLisandro Dalcin   ierr   = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
982c3122656SLisandro Dalcin   if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); }
983c84f5b01SHong Zhang   a->row = row;
984c3122656SLisandro Dalcin   ierr   = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
985c3122656SLisandro Dalcin   if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); }
986c84f5b01SHong Zhang   a->col = row;
987c84f5b01SHong Zhang 
988c84f5b01SHong Zhang   /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
989c84f5b01SHong Zhang   if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);}
99052e6d16bSBarry Smith   ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr);
99149b5e25fSSatish Balay 
99249b5e25fSSatish Balay   if (!a->solve_work) {
993c84f5b01SHong Zhang     ierr = PetscMalloc((inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
994c84f5b01SHong Zhang     ierr = PetscLogObjectMemory(inA,(inA->rmap.N+inA->rmap.bs)*sizeof(PetscScalar));CHKERRQ(ierr);
99549b5e25fSSatish Balay   }
99649b5e25fSSatish Balay 
997af281ebdSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr);
99849b5e25fSSatish Balay   PetscFunctionReturn(0);
99949b5e25fSSatish Balay }
1000950f1e5bSHong Zhang 
100149b5e25fSSatish Balay EXTERN_C_BEGIN
10024a2ae208SSatish Balay #undef __FUNCT__
10034a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
1004be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
100549b5e25fSSatish Balay {
1006045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
100713f74950SBarry Smith   PetscInt     i,nz,n;
100849b5e25fSSatish Balay 
100949b5e25fSSatish Balay   PetscFunctionBegin;
10106c6c5352SBarry Smith   nz = baij->maxnz;
1011899cda47SBarry Smith   n  = mat->cmap.n;
101249b5e25fSSatish Balay   for (i=0; i<nz; i++) {
101349b5e25fSSatish Balay     baij->j[i] = indices[i];
101449b5e25fSSatish Balay   }
10156c6c5352SBarry Smith    baij->nz = nz;
101649b5e25fSSatish Balay    for (i=0; i<n; i++) {
101749b5e25fSSatish Balay      baij->ilen[i] = baij->imax[i];
101849b5e25fSSatish Balay    }
101949b5e25fSSatish Balay    PetscFunctionReturn(0);
102049b5e25fSSatish Balay }
102149b5e25fSSatish Balay EXTERN_C_END
102249b5e25fSSatish Balay 
10234a2ae208SSatish Balay #undef __FUNCT__
10244a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
102549b5e25fSSatish Balay /*@
102619585528SSatish Balay   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
102749b5e25fSSatish Balay   in the matrix.
102849b5e25fSSatish Balay 
102949b5e25fSSatish Balay   Input Parameters:
103019585528SSatish Balay   +  mat     - the SeqSBAIJ matrix
103149b5e25fSSatish Balay   -  indices - the column indices
103249b5e25fSSatish Balay 
103349b5e25fSSatish Balay   Level: advanced
103449b5e25fSSatish Balay 
103549b5e25fSSatish Balay   Notes:
103649b5e25fSSatish Balay   This can be called if you have precomputed the nonzero structure of the
103749b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
103849b5e25fSSatish Balay   of the MatSetValues() operation.
103949b5e25fSSatish Balay 
104049b5e25fSSatish Balay   You MUST have set the correct numbers of nonzeros per row in the call to
1041d1be2dadSMatthew Knepley   MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
104249b5e25fSSatish Balay 
1043ab9f2c04SSatish Balay   MUST be called before any calls to MatSetValues()
104449b5e25fSSatish Balay 
1045ab9f2c04SSatish Balay   .seealso: MatCreateSeqSBAIJ
104649b5e25fSSatish Balay @*/
1047be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
104849b5e25fSSatish Balay {
104913f74950SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt *);
105049b5e25fSSatish Balay 
105149b5e25fSSatish Balay   PetscFunctionBegin;
10524482741eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
10534482741eSBarry Smith   PetscValidPointer(indices,2);
1054c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
105549b5e25fSSatish Balay   if (f) {
105649b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
105749b5e25fSSatish Balay   } else {
1058e005ede5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
105949b5e25fSSatish Balay   }
106049b5e25fSSatish Balay   PetscFunctionReturn(0);
106149b5e25fSSatish Balay }
106249b5e25fSSatish Balay 
10634a2ae208SSatish Balay #undef __FUNCT__
10643c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ"
10653c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
10663c896bc6SHong Zhang {
10673c896bc6SHong Zhang   PetscErrorCode ierr;
10683c896bc6SHong Zhang 
10693c896bc6SHong Zhang   PetscFunctionBegin;
10703c896bc6SHong Zhang   /* If the two matrices have the same copy implementation, use fast copy. */
10713c896bc6SHong Zhang   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
10723c896bc6SHong Zhang     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10733c896bc6SHong Zhang     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
10743c896bc6SHong Zhang 
1075899cda47SBarry Smith     if (a->i[A->rmap.N] != b->i[B->rmap.N]) {
10763c896bc6SHong Zhang       SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
10773c896bc6SHong Zhang     }
1078899cda47SBarry Smith     ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.N])*sizeof(PetscScalar));CHKERRQ(ierr);
10793c896bc6SHong Zhang   } else {
1080f5edf698SHong Zhang     ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
10813c896bc6SHong Zhang     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1082f5edf698SHong Zhang     ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
10833c896bc6SHong Zhang   }
10843c896bc6SHong Zhang   PetscFunctionReturn(0);
10853c896bc6SHong Zhang }
10863c896bc6SHong Zhang 
10873c896bc6SHong Zhang #undef __FUNCT__
10884a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1089dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A)
1090273d9f13SBarry Smith {
1091dfbe8321SBarry Smith   PetscErrorCode ierr;
1092273d9f13SBarry Smith 
1093273d9f13SBarry Smith   PetscFunctionBegin;
10947edd0491SSatish Balay   ierr =  MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,PetscMax(A->rmap.bs,1),PETSC_DEFAULT,0);CHKERRQ(ierr);
1095273d9f13SBarry Smith   PetscFunctionReturn(0);
1096273d9f13SBarry Smith }
1097273d9f13SBarry Smith 
1098a6ece127SHong Zhang #undef __FUNCT__
1099a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1100dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1101a6ece127SHong Zhang {
1102a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1103a6ece127SHong Zhang   PetscFunctionBegin;
1104a6ece127SHong Zhang   *array = a->a;
1105a6ece127SHong Zhang   PetscFunctionReturn(0);
1106a6ece127SHong Zhang }
1107a6ece127SHong Zhang 
1108a6ece127SHong Zhang #undef __FUNCT__
1109a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1110dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1111a6ece127SHong Zhang {
1112a6ece127SHong Zhang   PetscFunctionBegin;
1113a6ece127SHong Zhang   PetscFunctionReturn(0);
1114a6ece127SHong Zhang  }
1115a6ece127SHong Zhang 
111642ee4b1aSHong Zhang #include "petscblaslapack.h"
111742ee4b1aSHong Zhang #undef __FUNCT__
111842ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1119f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
112042ee4b1aSHong Zhang {
112142ee4b1aSHong Zhang   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1122dfbe8321SBarry Smith   PetscErrorCode ierr;
1123899cda47SBarry Smith   PetscInt       i,bs=Y->rmap.bs,bs2,j;
11244ce68768SBarry Smith   PetscBLASInt   bnz = (PetscBLASInt)x->nz,one = 1;
112542ee4b1aSHong Zhang 
112642ee4b1aSHong Zhang   PetscFunctionBegin;
112742ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1128f4df32b1SMatthew Knepley     PetscScalar alpha = a;
1129f4df32b1SMatthew Knepley     BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
1130c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1131c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1132c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1133c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1134c537a176SHong Zhang     }
1135c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1136c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1137c4319e64SHong Zhang       y->XtoY = X;
1138c537a176SHong Zhang     }
1139c4319e64SHong Zhang     bs2 = bs*bs;
11406c6c5352SBarry Smith     for (i=0; i<x->nz; i++) {
1141c4319e64SHong Zhang       j = 0;
1142c4319e64SHong Zhang       while (j < bs2){
1143f4df32b1SMatthew Knepley         y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]);
1144c4319e64SHong Zhang         j++;
1145c537a176SHong Zhang       }
1146c4319e64SHong Zhang     }
1147ae15b995SBarry Smith     ierr = PetscInfo3(0,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %G\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));CHKERRQ(ierr);
114842ee4b1aSHong Zhang   } else {
1149f5edf698SHong Zhang     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
1150f4df32b1SMatthew Knepley     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
1151f5edf698SHong Zhang     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
115242ee4b1aSHong Zhang   }
115342ee4b1aSHong Zhang   PetscFunctionReturn(0);
115442ee4b1aSHong Zhang }
115542ee4b1aSHong Zhang 
1156efcf0fc3SBarry Smith #undef __FUNCT__
1157efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1158dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg)
1159efcf0fc3SBarry Smith {
1160efcf0fc3SBarry Smith   PetscFunctionBegin;
1161efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1162efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1163efcf0fc3SBarry Smith }
1164efcf0fc3SBarry Smith 
1165efcf0fc3SBarry Smith #undef __FUNCT__
1166efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1167dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1168efcf0fc3SBarry Smith {
1169efcf0fc3SBarry Smith    PetscFunctionBegin;
1170efcf0fc3SBarry Smith    *flg = PETSC_TRUE;
1171efcf0fc3SBarry Smith    PetscFunctionReturn(0);
1172efcf0fc3SBarry Smith }
1173efcf0fc3SBarry Smith 
1174efcf0fc3SBarry Smith #undef __FUNCT__
1175efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1176ab5e4463SMatthew Knepley PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg)
1177efcf0fc3SBarry Smith  {
1178efcf0fc3SBarry Smith    PetscFunctionBegin;
1179efcf0fc3SBarry Smith    *flg = PETSC_FALSE;
1180efcf0fc3SBarry Smith    PetscFunctionReturn(0);
1181efcf0fc3SBarry Smith  }
1182efcf0fc3SBarry Smith 
118399cafbc1SBarry Smith #undef __FUNCT__
118499cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ"
118599cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
118699cafbc1SBarry Smith {
118799cafbc1SBarry Smith   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
118899cafbc1SBarry Smith   PetscInt       i,nz = a->bs2*a->i[a->mbs];
118999cafbc1SBarry Smith   PetscScalar    *aa = a->a;
119099cafbc1SBarry Smith 
119199cafbc1SBarry Smith   PetscFunctionBegin;
119299cafbc1SBarry Smith   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
119399cafbc1SBarry Smith   PetscFunctionReturn(0);
119499cafbc1SBarry Smith }
119599cafbc1SBarry Smith 
119699cafbc1SBarry Smith #undef __FUNCT__
119799cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ"
119899cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
119999cafbc1SBarry Smith {
120099cafbc1SBarry Smith   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
120199cafbc1SBarry Smith   PetscInt       i,nz = a->bs2*a->i[a->mbs];
120299cafbc1SBarry Smith   PetscScalar    *aa = a->a;
120399cafbc1SBarry Smith 
120499cafbc1SBarry Smith   PetscFunctionBegin;
120599cafbc1SBarry Smith   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
120699cafbc1SBarry Smith   PetscFunctionReturn(0);
120799cafbc1SBarry Smith }
120899cafbc1SBarry Smith 
120949b5e25fSSatish Balay /* -------------------------------------------------------------------*/
121049b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
121149b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
121249b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
121349b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
121497304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
1215431c96f7SBarry Smith        MatMult_SeqSBAIJ_N,       /* transpose versions are same as non-transpose versions */
1216e005ede5SBarry Smith        MatMultAdd_SeqSBAIJ_N,
121749b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
121849b5e25fSSatish Balay        0,
121949b5e25fSSatish Balay        0,
122097304618SKris Buschelman /*10*/ 0,
122149b5e25fSSatish Balay        0,
12225f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1223d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
122449b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
122597304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
122649b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
122749b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
122849b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
122949b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
123097304618SKris Buschelman /*20*/ 0,
123149b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
123249b5e25fSSatish Balay        0,
123349b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
123449b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
1235dcf5cc72SBarry Smith /*25*/ 0,
123649b5e25fSSatish Balay        0,
123749b5e25fSSatish Balay        0,
12385f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
12395f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
124097304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1241c464158bSHong Zhang        0,
12424d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1243a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1244a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
124597304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
12469495ac64SHong Zhang        MatForwardSolve_SeqSBAIJ_N,
12479495ac64SHong Zhang        MatBackwardSolve_SeqSBAIJ_N,
124849b5e25fSSatish Balay        0,
1249c84f5b01SHong Zhang        MatICCFactor_SeqSBAIJ,
125097304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
125149b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
125249b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
125349b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
12543c896bc6SHong Zhang        MatCopy_SeqSBAIJ,
12558c07d4e3SBarry Smith /*45*/ 0,
125649b5e25fSSatish Balay        MatScale_SeqSBAIJ,
125749b5e25fSSatish Balay        0,
125849b5e25fSSatish Balay        0,
125949b5e25fSSatish Balay        0,
1260521d7252SBarry Smith /*50*/ 0,
126149b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
126249b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
126349b5e25fSSatish Balay        0,
126449b5e25fSSatish Balay        0,
126597304618SKris Buschelman /*55*/ 0,
126649b5e25fSSatish Balay        0,
126749b5e25fSSatish Balay        0,
126849b5e25fSSatish Balay        0,
126949b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
127097304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
127149b5e25fSSatish Balay        0,
127249b5e25fSSatish Balay        0,
1273357abbc8SBarry Smith        0,
1274d959ec07SHong Zhang        0,
127597304618SKris Buschelman /*65*/ 0,
1276d959ec07SHong Zhang        0,
1277d959ec07SHong Zhang        0,
1278d959ec07SHong Zhang        0,
1279d959ec07SHong Zhang        0,
1280985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqSBAIJ,
12813e0d88b5SBarry Smith        0,
12823e0d88b5SBarry Smith        0,
12833e0d88b5SBarry Smith        0,
12843e0d88b5SBarry Smith        0,
128597304618SKris Buschelman /*75*/ 0,
12863e0d88b5SBarry Smith        0,
12873e0d88b5SBarry Smith        0,
12883e0d88b5SBarry Smith        0,
12893e0d88b5SBarry Smith        0,
129097304618SKris Buschelman /*80*/ 0,
12913e0d88b5SBarry Smith        0,
12923e0d88b5SBarry Smith        0,
12933e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
129497304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
12953e0d88b5SBarry Smith #else
129697304618SKris Buschelman        0,
12973e0d88b5SBarry Smith #endif
1298865e5f61SKris Buschelman        MatLoad_SeqSBAIJ,
1299865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ,
1300865e5f61SKris Buschelman        MatIsHermitian_SeqSBAIJ,
1301efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1302865e5f61SKris Buschelman        0,
1303865e5f61SKris Buschelman        0,
1304865e5f61SKris Buschelman /*90*/ 0,
1305865e5f61SKris Buschelman        0,
1306865e5f61SKris Buschelman        0,
1307865e5f61SKris Buschelman        0,
1308865e5f61SKris Buschelman        0,
1309865e5f61SKris Buschelman /*95*/ 0,
1310865e5f61SKris Buschelman        0,
1311865e5f61SKris Buschelman        0,
131299cafbc1SBarry Smith        0,
131399cafbc1SBarry Smith        0,
131499cafbc1SBarry Smith /*100*/0,
131599cafbc1SBarry Smith        0,
131699cafbc1SBarry Smith        0,
131799cafbc1SBarry Smith        0,
131899cafbc1SBarry Smith        0,
131999cafbc1SBarry Smith /*105*/0,
132099cafbc1SBarry Smith        MatRealPart_SeqSBAIJ,
1321f5edf698SHong Zhang        MatImaginaryPart_SeqSBAIJ,
1322f5edf698SHong Zhang        MatGetRowUpperTriangular_SeqSBAIJ,
1323f5edf698SHong Zhang        MatRestoreRowUpperTriangular_SeqSBAIJ
132499cafbc1SBarry Smith };
1325be1d678aSKris Buschelman 
132649b5e25fSSatish Balay EXTERN_C_BEGIN
13274a2ae208SSatish Balay #undef __FUNCT__
13284a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
1329be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat)
133049b5e25fSSatish Balay {
13314afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ *)mat->data;
1332899cda47SBarry Smith   PetscInt       nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2;
1333dfbe8321SBarry Smith   PetscErrorCode ierr;
133449b5e25fSSatish Balay 
133549b5e25fSSatish Balay   PetscFunctionBegin;
133649b5e25fSSatish Balay   if (aij->nonew != 1) {
1337*4e0d8c25SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS,PETSC_TRUE);first");
133849b5e25fSSatish Balay   }
133949b5e25fSSatish Balay 
134049b5e25fSSatish Balay   /* allocate space for values if not already there */
134149b5e25fSSatish Balay   if (!aij->saved_values) {
134287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
134349b5e25fSSatish Balay   }
134449b5e25fSSatish Balay 
134549b5e25fSSatish Balay   /* copy values over */
134687828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
134749b5e25fSSatish Balay   PetscFunctionReturn(0);
134849b5e25fSSatish Balay }
134949b5e25fSSatish Balay EXTERN_C_END
135049b5e25fSSatish Balay 
135149b5e25fSSatish Balay EXTERN_C_BEGIN
13524a2ae208SSatish Balay #undef __FUNCT__
13534a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
1354be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat)
135549b5e25fSSatish Balay {
13564afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ *)mat->data;
13576849ba73SBarry Smith   PetscErrorCode ierr;
1358899cda47SBarry Smith   PetscInt       nz = aij->i[mat->rmap.N]*mat->rmap.bs*aij->bs2;
135949b5e25fSSatish Balay 
136049b5e25fSSatish Balay   PetscFunctionBegin;
136149b5e25fSSatish Balay   if (aij->nonew != 1) {
1362*4e0d8c25SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS,PETSC_TRUE);first");
136349b5e25fSSatish Balay   }
136449b5e25fSSatish Balay   if (!aij->saved_values) {
1365e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
136649b5e25fSSatish Balay   }
136749b5e25fSSatish Balay 
136849b5e25fSSatish Balay   /* copy values over */
136987828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
137049b5e25fSSatish Balay   PetscFunctionReturn(0);
137149b5e25fSSatish Balay }
137249b5e25fSSatish Balay EXTERN_C_END
137349b5e25fSSatish Balay 
13748549e402SHong Zhang EXTERN_C_BEGIN
13754a2ae208SSatish Balay #undef __FUNCT__
1376a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1377be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
137849b5e25fSSatish Balay {
1379c464158bSHong Zhang   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
13806849ba73SBarry Smith   PetscErrorCode ierr;
1381ab93d7beSBarry Smith   PetscInt       i,mbs,bs2;
1382ab93d7beSBarry Smith   PetscTruth     skipallocation = PETSC_FALSE,flg;
138349b5e25fSSatish Balay 
138449b5e25fSSatish Balay   PetscFunctionBegin;
1385273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1386e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1387899cda47SBarry Smith   B->rmap.bs = B->cmap.bs = bs;
13886148ca0dSBarry Smith   ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr);
13896148ca0dSBarry Smith   ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr);
1390899cda47SBarry Smith 
1391899cda47SBarry Smith   mbs  = B->rmap.N/bs;
139249b5e25fSSatish Balay   bs2  = bs*bs;
139349b5e25fSSatish Balay 
1394899cda47SBarry Smith   if (mbs*bs != B->rmap.N) {
139529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
139649b5e25fSSatish Balay   }
139749b5e25fSSatish Balay 
1398ab93d7beSBarry Smith   if (nz == MAT_SKIP_ALLOCATION) {
1399ab93d7beSBarry Smith     skipallocation = PETSC_TRUE;
1400ab93d7beSBarry Smith     nz             = 0;
1401ab93d7beSBarry Smith   }
1402ab93d7beSBarry Smith 
1403435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
140477431f27SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
140549b5e25fSSatish Balay   if (nnz) {
140649b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
140777431f27SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
140877431f27SBarry 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);
140949b5e25fSSatish Balay     }
141049b5e25fSSatish Balay   }
141149b5e25fSSatish Balay 
1412e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
141349b5e25fSSatish Balay   if (!flg) {
141449b5e25fSSatish Balay     switch (bs) {
141549b5e25fSSatish Balay     case 1:
14164ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
141749b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_1;
1418d59c15a7SBarry Smith       B->ops->solves           = MatSolves_SeqSBAIJ_1;
1419e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_1;
142049b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_1;
142149b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1422431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1423431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
142449b5e25fSSatish Balay       break;
142549b5e25fSSatish Balay     case 2:
14264ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
142749b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_2;
1428e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_2;
142949b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_2;
143049b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1431431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1432431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
143349b5e25fSSatish Balay       break;
143449b5e25fSSatish Balay     case 3:
14355f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
143649b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_3;
1437e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_3;
143849b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_3;
143949b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1440431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1441431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
144249b5e25fSSatish Balay       break;
144349b5e25fSSatish Balay     case 4:
14445f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
144549b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_4;
1446e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_4;
144749b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_4;
144849b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1449431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1450431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
145149b5e25fSSatish Balay       break;
145249b5e25fSSatish Balay     case 5:
14535f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
145449b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_5;
1455e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_5;
145649b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_5;
145749b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1458431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1459431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
146049b5e25fSSatish Balay       break;
146149b5e25fSSatish Balay     case 6:
14625f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
146349b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_6;
1464e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_6;
146549b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_6;
146649b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1467431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1468431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
146949b5e25fSSatish Balay       break;
147049b5e25fSSatish Balay     case 7:
1471de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
147249b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_7;
1473e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_7;
1474de53e5efSHong Zhang       B->ops->mult             = MatMult_SeqSBAIJ_7;
147549b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1476431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1477431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
147849b5e25fSSatish Balay       break;
147949b5e25fSSatish Balay     }
148049b5e25fSSatish Balay   }
148149b5e25fSSatish Balay 
148249b5e25fSSatish Balay   b->mbs = mbs;
14834afc71dfSHong Zhang   b->nbs = mbs;
1484ab93d7beSBarry Smith   if (!skipallocation) {
1485ab93d7beSBarry Smith     /* b->ilen will count nonzeros in each block row so far. */
1486ab93d7beSBarry Smith     ierr   = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr);
1487ab93d7beSBarry Smith     for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
148849b5e25fSSatish Balay     if (!nnz) {
1489435da068SBarry Smith       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
149049b5e25fSSatish Balay       else if (nz <= 0)        nz = 1;
149149b5e25fSSatish Balay       for (i=0; i<mbs; i++) {
14928cef66ccSHong Zhang         b->imax[i] = nz;
149349b5e25fSSatish Balay       }
1494153ea458SHong Zhang       nz = nz*mbs; /* total nz */
149549b5e25fSSatish Balay     } else {
149649b5e25fSSatish Balay       nz = 0;
14978cef66ccSHong Zhang       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
149849b5e25fSSatish Balay     }
14996c6c5352SBarry Smith     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
150049b5e25fSSatish Balay 
150149b5e25fSSatish Balay     /* allocate the matrix space */
1502899cda47SBarry Smith     ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.N+1,PetscInt,&b->i);CHKERRQ(ierr);
15036c6c5352SBarry Smith     ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
150413f74950SBarry Smith     ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
150549b5e25fSSatish Balay     b->singlemalloc = PETSC_TRUE;
150649b5e25fSSatish Balay 
150749b5e25fSSatish Balay     /* pointer to beginning of each row */
150849b5e25fSSatish Balay     b->i[0] = 0;
150949b5e25fSSatish Balay     for (i=1; i<mbs+1; i++) {
151049b5e25fSSatish Balay       b->i[i] = b->i[i-1] + b->imax[i-1];
151149b5e25fSSatish Balay     }
1512e6b907acSBarry Smith     b->free_a     = PETSC_TRUE;
1513e6b907acSBarry Smith     b->free_ij    = PETSC_TRUE;
1514e811da20SHong Zhang   } else {
1515e6b907acSBarry Smith     b->free_a     = PETSC_FALSE;
1516e6b907acSBarry Smith     b->free_ij    = PETSC_FALSE;
1517ab93d7beSBarry Smith   }
151849b5e25fSSatish Balay 
1519899cda47SBarry Smith   B->rmap.bs               = bs;
152049b5e25fSSatish Balay   b->bs2              = bs2;
15216c6c5352SBarry Smith   b->nz             = 0;
15226c6c5352SBarry Smith   b->maxnz          = nz*bs2;
1523153ea458SHong Zhang 
152416cdd363SHong Zhang   b->inew             = 0;
152516cdd363SHong Zhang   b->jnew             = 0;
152616cdd363SHong Zhang   b->anew             = 0;
152716cdd363SHong Zhang   b->a2anew           = 0;
15281a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1529c464158bSHong Zhang   PetscFunctionReturn(0);
1530c464158bSHong Zhang }
1531a23d5eceSKris Buschelman EXTERN_C_END
1532153ea458SHong Zhang 
1533d769727bSBarry Smith EXTERN_C_BEGIN
1534f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1535f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1536d769727bSBarry Smith EXTERN_C_END
1537d769727bSBarry Smith 
15380bad9183SKris Buschelman /*MC
1539fafad747SKris Buschelman   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
15400bad9183SKris Buschelman   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
15410bad9183SKris Buschelman 
15420bad9183SKris Buschelman   Options Database Keys:
15430bad9183SKris Buschelman   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
15440bad9183SKris Buschelman 
15450bad9183SKris Buschelman   Level: beginner
15460bad9183SKris Buschelman 
15470bad9183SKris Buschelman   .seealso: MatCreateSeqSBAIJ
15480bad9183SKris Buschelman M*/
15490bad9183SKris Buschelman 
1550a23d5eceSKris Buschelman EXTERN_C_BEGIN
1551a23d5eceSKris Buschelman #undef __FUNCT__
1552a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1553be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B)
1554a23d5eceSKris Buschelman {
1555a23d5eceSKris Buschelman   Mat_SeqSBAIJ   *b;
1556dfbe8321SBarry Smith   PetscErrorCode ierr;
155713f74950SBarry Smith   PetscMPIInt    size;
1558941593c8SHong Zhang   PetscTruth     flg;
1559a23d5eceSKris Buschelman 
1560a23d5eceSKris Buschelman   PetscFunctionBegin;
1561a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1562a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1563a23d5eceSKris Buschelman 
156438f2d2fdSLisandro Dalcin   ierr    = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1565a23d5eceSKris Buschelman   B->data = (void*)b;
1566a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1567a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1568a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1569a23d5eceSKris Buschelman   B->factor           = 0;
1570a23d5eceSKris Buschelman   B->mapping          = 0;
1571a23d5eceSKris Buschelman   b->row              = 0;
1572a23d5eceSKris Buschelman   b->icol             = 0;
1573a23d5eceSKris Buschelman   b->reallocs         = 0;
1574a23d5eceSKris Buschelman   b->saved_values     = 0;
1575a23d5eceSKris Buschelman 
1576a23d5eceSKris Buschelman 
1577a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1578a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1579a23d5eceSKris Buschelman   b->nonew            = 0;
1580a23d5eceSKris Buschelman   b->diag             = 0;
1581a23d5eceSKris Buschelman   b->solve_work       = 0;
1582a23d5eceSKris Buschelman   b->mult_work        = 0;
1583a23d5eceSKris Buschelman   B->spptr            = 0;
1584a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1585a23d5eceSKris Buschelman   b->xtoy             = 0;
1586a23d5eceSKris Buschelman   b->XtoY             = 0;
1587a23d5eceSKris Buschelman 
1588a23d5eceSKris Buschelman   b->inew             = 0;
1589a23d5eceSKris Buschelman   b->jnew             = 0;
1590a23d5eceSKris Buschelman   b->anew             = 0;
1591a23d5eceSKris Buschelman   b->a2anew           = 0;
1592a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1593a23d5eceSKris Buschelman 
1594941593c8SHong Zhang   b->ignore_ltriangular = PETSC_FALSE;
1595941593c8SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr);
1596941593c8SHong Zhang   if (flg) b->ignore_ltriangular = PETSC_TRUE;
1597941593c8SHong Zhang 
1598f5edf698SHong Zhang   b->getrow_utriangular = PETSC_FALSE;
1599f5edf698SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr);
1600f5edf698SHong Zhang   if (flg) b->getrow_utriangular = PETSC_TRUE;
1601f5edf698SHong Zhang 
1602a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1603a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1604a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1605a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1606a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1607a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1608a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1609a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1610a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
16114e5e7fe4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",
16124e5e7fe4SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqAIJ",
16134e5e7fe4SHong Zhang                                       MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1614a0e1a404SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",
1615a0e1a404SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqBAIJ",
1616a0e1a404SHong Zhang                                       MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1617a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1618a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1619a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
162023ce1328SBarry Smith 
162123ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
162223ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
162323ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
162423ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
162517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr);
1626a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1627a23d5eceSKris Buschelman }
1628a23d5eceSKris Buschelman EXTERN_C_END
1629a23d5eceSKris Buschelman 
1630a23d5eceSKris Buschelman #undef __FUNCT__
1631a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1632a23d5eceSKris Buschelman /*@C
1633a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1634a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1635a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1636a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1637a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1638a23d5eceSKris Buschelman 
1639a23d5eceSKris Buschelman    Collective on Mat
1640a23d5eceSKris Buschelman 
1641a23d5eceSKris Buschelman    Input Parameters:
1642a23d5eceSKris Buschelman +  A - the symmetric matrix
1643a23d5eceSKris Buschelman .  bs - size of block
1644a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1645a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1646a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1647a23d5eceSKris Buschelman 
1648a23d5eceSKris Buschelman    Options Database Keys:
1649a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1650a23d5eceSKris Buschelman                      block calculations (much slower)
1651a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1652a23d5eceSKris Buschelman 
1653a23d5eceSKris Buschelman    Level: intermediate
1654a23d5eceSKris Buschelman 
1655a23d5eceSKris Buschelman    Notes:
1656a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1657a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1658a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1659a23d5eceSKris Buschelman    matrices.
1660a23d5eceSKris Buschelman 
166149a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
166249a6f317SBarry Smith 
166349a6f317SBarry Smith 
1664a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1665a23d5eceSKris Buschelman @*/
1666be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
166713f74950SBarry Smith {
166813f74950SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]);
1669a23d5eceSKris Buschelman 
1670a23d5eceSKris Buschelman   PetscFunctionBegin;
1671a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1672a23d5eceSKris Buschelman   if (f) {
1673a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1674a23d5eceSKris Buschelman   }
1675a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1676a23d5eceSKris Buschelman }
167749b5e25fSSatish Balay 
16784a2ae208SSatish Balay #undef __FUNCT__
16794a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1680c464158bSHong Zhang /*@C
1681c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1682c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1683c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1684c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1685c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
168649b5e25fSSatish Balay 
1687c464158bSHong Zhang    Collective on MPI_Comm
1688c464158bSHong Zhang 
1689c464158bSHong Zhang    Input Parameters:
1690c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1691c464158bSHong Zhang .  bs - size of block
1692c464158bSHong Zhang .  m - number of rows, or number of columns
1693c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1694744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1695744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1696c464158bSHong Zhang 
1697c464158bSHong Zhang    Output Parameter:
1698c464158bSHong Zhang .  A - the symmetric matrix
1699c464158bSHong Zhang 
1700c464158bSHong Zhang    Options Database Keys:
1701c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1702c464158bSHong Zhang                      block calculations (much slower)
1703c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1704c464158bSHong Zhang 
1705c464158bSHong Zhang    Level: intermediate
1706c464158bSHong Zhang 
1707c464158bSHong Zhang    Notes:
17086d6d819aSHong Zhang    The number of rows and columns must be divisible by blocksize.
17096d6d819aSHong Zhang    This matrix type does not support complex Hermitian operation.
1710c464158bSHong Zhang 
1711c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1712c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1713c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1714c464158bSHong Zhang    matrices.
1715c464158bSHong Zhang 
171649a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
171749a6f317SBarry Smith 
1718c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1719c464158bSHong Zhang @*/
1720be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
1721c464158bSHong Zhang {
1722dfbe8321SBarry Smith   PetscErrorCode ierr;
1723c464158bSHong Zhang 
1724c464158bSHong Zhang   PetscFunctionBegin;
1725f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1726f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1727c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1728ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr);
172949b5e25fSSatish Balay   PetscFunctionReturn(0);
173049b5e25fSSatish Balay }
173149b5e25fSSatish Balay 
17324a2ae208SSatish Balay #undef __FUNCT__
17334a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
1734dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
173549b5e25fSSatish Balay {
173649b5e25fSSatish Balay   Mat            C;
173749b5e25fSSatish Balay   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
17386849ba73SBarry Smith   PetscErrorCode ierr;
1739b40805acSSatish Balay   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
174049b5e25fSSatish Balay 
174149b5e25fSSatish Balay   PetscFunctionBegin;
174229bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
174349b5e25fSSatish Balay 
174449b5e25fSSatish Balay   *B = 0;
1745f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
1746899cda47SBarry Smith   ierr = MatSetSizes(C,A->rmap.N,A->cmap.n,A->rmap.N,A->cmap.n);CHKERRQ(ierr);
1747be5d1d56SKris Buschelman   ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
17481d5dac46SHong Zhang   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1749692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1750692f9cbeSHong Zhang 
1751273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
175249b5e25fSSatish Balay   C->factor         = A->factor;
175349b5e25fSSatish Balay   c->row            = 0;
175449b5e25fSSatish Balay   c->icol           = 0;
175549b5e25fSSatish Balay   c->saved_values   = 0;
175649b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
175749b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
175849b5e25fSSatish Balay 
1759899cda47SBarry Smith   ierr = PetscMapCopy(A->comm,&A->rmap,&C->rmap);CHKERRQ(ierr);
1760899cda47SBarry Smith   ierr = PetscMapCopy(A->comm,&A->cmap,&C->cmap);CHKERRQ(ierr);
176149b5e25fSSatish Balay   c->bs2  = a->bs2;
176249b5e25fSSatish Balay   c->mbs  = a->mbs;
176349b5e25fSSatish Balay   c->nbs  = a->nbs;
176449b5e25fSSatish Balay 
17658777fc3fSSatish Balay   ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr);
176649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
176749b5e25fSSatish Balay     c->imax[i] = a->imax[i];
176849b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
176949b5e25fSSatish Balay   }
177049b5e25fSSatish Balay 
177149b5e25fSSatish Balay   /* allocate the matrix space */
1772b40805acSSatish Balay   ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr);
177349b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
177413f74950SBarry Smith   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
1775b40805acSSatish Balay   ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr);
177649b5e25fSSatish Balay   if (mbs > 0) {
177713f74950SBarry Smith     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
177849b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
177949b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178049b5e25fSSatish Balay     } else {
178149b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
178249b5e25fSSatish Balay     }
178349b5e25fSSatish Balay   }
178449b5e25fSSatish Balay 
178549b5e25fSSatish Balay   c->sorted      = a->sorted;
178649b5e25fSSatish Balay   c->roworiented = a->roworiented;
178749b5e25fSSatish Balay   c->nonew       = a->nonew;
178849b5e25fSSatish Balay 
178949b5e25fSSatish Balay   if (a->diag) {
179013f74950SBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr);
179152e6d16bSBarry Smith     ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
179249b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
179349b5e25fSSatish Balay       c->diag[i] = a->diag[i];
179449b5e25fSSatish Balay     }
179549b5e25fSSatish Balay   } else c->diag  = 0;
17966c6c5352SBarry Smith   c->nz           = a->nz;
17976c6c5352SBarry Smith   c->maxnz        = a->maxnz;
179849b5e25fSSatish Balay   c->solve_work   = 0;
179949b5e25fSSatish Balay   c->mult_work    = 0;
1800e6b907acSBarry Smith   c->free_a       = PETSC_TRUE;
1801e6b907acSBarry Smith   c->free_ij      = PETSC_TRUE;
180249b5e25fSSatish Balay   *B = C;
1803b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
180449b5e25fSSatish Balay   PetscFunctionReturn(0);
180549b5e25fSSatish Balay }
180649b5e25fSSatish Balay 
18074a2ae208SSatish Balay #undef __FUNCT__
18084a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1809f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A)
181049b5e25fSSatish Balay {
181149b5e25fSSatish Balay   Mat_SeqSBAIJ   *a;
181249b5e25fSSatish Balay   Mat            B;
18136849ba73SBarry Smith   PetscErrorCode ierr;
181413f74950SBarry Smith   int            fd;
181513f74950SBarry Smith   PetscMPIInt    size;
181613f74950SBarry Smith   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs=1;
181713f74950SBarry Smith   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
181813f74950SBarry Smith   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows;
181913f74950SBarry Smith   PetscInt       *masked,nmask,tmp,bs2,ishift;
182087828ca2SBarry Smith   PetscScalar    *aa;
182149b5e25fSSatish Balay   MPI_Comm       comm = ((PetscObject)viewer)->comm;
182249b5e25fSSatish Balay 
182349b5e25fSSatish Balay   PetscFunctionBegin;
1824b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
182549b5e25fSSatish Balay   bs2  = bs*bs;
182649b5e25fSSatish Balay 
182749b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
182829bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1829b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
183049b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1831552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
183249b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
183349b5e25fSSatish Balay 
183449b5e25fSSatish Balay   if (header[3] < 0) {
183529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
183649b5e25fSSatish Balay   }
183749b5e25fSSatish Balay 
183829bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
183949b5e25fSSatish Balay 
184049b5e25fSSatish Balay   /*
184149b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
184249b5e25fSSatish Balay     divisible by the blocksize
184349b5e25fSSatish Balay   */
184449b5e25fSSatish Balay   mbs        = M/bs;
184549b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
184649b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
184749b5e25fSSatish Balay   else                  mbs++;
184849b5e25fSSatish Balay   if (extra_rows) {
1849ae15b995SBarry Smith     ierr = PetscInfo(0,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr);
185049b5e25fSSatish Balay   }
185149b5e25fSSatish Balay 
185249b5e25fSSatish Balay   /* read in row lengths */
185313f74950SBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
185449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
185549b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
185649b5e25fSSatish Balay 
185749b5e25fSSatish Balay   /* read in column indices */
185813f74950SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr);
185949b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
186049b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
186149b5e25fSSatish Balay 
186249b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
186313f74950SBarry Smith   ierr     = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr);
186413f74950SBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr);
186513f74950SBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr);
186613f74950SBarry Smith   ierr     = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr);
186749b5e25fSSatish Balay   masked   = mask + mbs;
186849b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
186949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
187049b5e25fSSatish Balay     nmask = 0;
187149b5e25fSSatish Balay     for (j=0; j<bs; j++) {
187249b5e25fSSatish Balay       kmax = rowlengths[rowcount];
187349b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18742d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
187503630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
187649b5e25fSSatish Balay       }
187749b5e25fSSatish Balay       rowcount++;
187849b5e25fSSatish Balay     }
1879574b2666SHong Zhang     s_browlengths[i] += nmask;
1880574b2666SHong Zhang 
188149b5e25fSSatish Balay     /* zero out the mask elements we set */
188249b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
188349b5e25fSSatish Balay   }
188449b5e25fSSatish Balay 
188549b5e25fSSatish Balay   /* create our matrix */
1886f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
1887f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr);
18889abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
1889ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr);
189049b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
189149b5e25fSSatish Balay 
189249b5e25fSSatish Balay   /* set matrix "i" values */
189349b5e25fSSatish Balay   a->i[0] = 0;
189449b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1895574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1896574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
189749b5e25fSSatish Balay   }
18986c6c5352SBarry Smith   a->nz = a->i[mbs];
189949b5e25fSSatish Balay 
190049b5e25fSSatish Balay   /* read in nonzero values */
190187828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
190249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
190349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
190449b5e25fSSatish Balay 
190549b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
190649b5e25fSSatish Balay   nzcount = 0; jcount = 0;
190749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
190849b5e25fSSatish Balay     nzcountb = nzcount;
190949b5e25fSSatish Balay     nmask    = 0;
191049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
191149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
191249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19132d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
191403630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19152d703238SHong Zhang       }
19162d703238SHong Zhang     }
19172d703238SHong Zhang     /* sort the masked values */
19182d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19192d703238SHong Zhang 
19202d703238SHong Zhang     /* set "j" values into matrix */
19212d703238SHong Zhang     maskcount = 1;
19222d703238SHong Zhang     for (j=0; j<nmask; j++) {
192349b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
192449b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
192549b5e25fSSatish Balay     }
1926574b2666SHong Zhang 
192749b5e25fSSatish Balay     /* set "a" values into matrix */
192849b5e25fSSatish Balay     ishift = bs2*a->i[i];
192949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
193049b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
193149b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1932574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1933574b2666SHong Zhang         if (tmp >= i){
193449b5e25fSSatish Balay           block     = mask[tmp] - 1;
193549b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
193649b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1937574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1938574b2666SHong Zhang         }
1939574b2666SHong Zhang         nzcountb++;
194049b5e25fSSatish Balay       }
194149b5e25fSSatish Balay     }
194249b5e25fSSatish Balay     /* zero out the mask elements we set */
194349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
194449b5e25fSSatish Balay   }
19456c6c5352SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
194649b5e25fSSatish Balay 
194749b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1948574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
194949b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
195049b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
195149b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
195249b5e25fSSatish Balay 
19539abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19549abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
195549b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19569abb65ffSKris Buschelman   *A = B;
195749b5e25fSSatish Balay   PetscFunctionReturn(0);
195849b5e25fSSatish Balay }
1959574b2666SHong Zhang 
1960d06b337dSHong Zhang #undef __FUNCT__
1961d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
196213f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1963d06b337dSHong Zhang {
1964d06b337dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
1965d06b337dSHong Zhang   MatScalar      *aa=a->a,*v,*v1;
1966d06b337dSHong Zhang   PetscScalar    *x,*b,*t,sum,d;
19676849ba73SBarry Smith   PetscErrorCode ierr;
1968899cda47SBarry Smith   PetscInt       m=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j;
1969521d7252SBarry Smith   PetscInt       nz,nz1,*vj,*vj1,i;
1970d06b337dSHong Zhang 
1971d06b337dSHong Zhang   PetscFunctionBegin;
197251f519a2SBarry Smith   if (flag & SOR_EISENSTAT) SETERRQ(PETSC_ERR_SUP,"No support yet for Eisenstat");
197351f519a2SBarry Smith 
1974b965ef7fSBarry Smith   its = its*lits;
197577431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1976d06b337dSHong Zhang 
1977d06b337dSHong Zhang   if (bs > 1)
1978d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1979d06b337dSHong Zhang 
19801ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1981d06b337dSHong Zhang   if (xx != bb) {
19821ebc52fbSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1983d06b337dSHong Zhang   } else {
1984d06b337dSHong Zhang     b = x;
1985d06b337dSHong Zhang   }
1986d06b337dSHong Zhang 
1987e2ee2a47SBarry Smith   if (!a->relax_work) {
1988e2ee2a47SBarry Smith     ierr = PetscMalloc(m*sizeof(PetscScalar),&a->relax_work);CHKERRQ(ierr);
1989e2ee2a47SBarry Smith   }
1990e2ee2a47SBarry Smith   t = a->relax_work;
1991d06b337dSHong Zhang 
1992d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19932798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1994290bbb0aSBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
1995d06b337dSHong Zhang 
1996d06b337dSHong Zhang       for (i=0; i<m; i++){
199744706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1998d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1999d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2000d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2001e6b907acSBarry Smith         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2002d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2003d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
2004d06b337dSHong Zhang       }
2005d06b337dSHong Zhang     }
2006d06b337dSHong Zhang 
20072798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
200895d750ceSBarry Smith       if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){
200995d750ceSBarry Smith         t = b;
2010d06b337dSHong Zhang       }
201195d750ceSBarry Smith 
2012d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2013d06b337dSHong Zhang         d  = *(aa + ai[i]);
2014d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2015d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2016d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2017e6b907acSBarry Smith         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2018d06b337dSHong Zhang         sum = t[i];
2019d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
2020d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2021d06b337dSHong Zhang       }
202295d750ceSBarry Smith       t = a->relax_work;
2023d06b337dSHong Zhang     }
2024d06b337dSHong Zhang     its--;
2025d06b337dSHong Zhang   }
2026d06b337dSHong Zhang 
2027d06b337dSHong Zhang   while (its--) {
202844706875SHong Zhang     /*
202944706875SHong Zhang        forward sweep:
203044706875SHong Zhang        for i=0,...,m-1:
203144706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
203244706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
203344706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2034d06b337dSHong Zhang 
203544706875SHong Zhang     */
20362798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2037290bbb0aSBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
2038d06b337dSHong Zhang 
2039d06b337dSHong Zhang       for (i=0; i<m; i++){
204044706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2041d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2042d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2043d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2044d06b337dSHong Zhang         sum = t[i];
2045e6b907acSBarry Smith         ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr);
2046d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2047d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2048d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
2049d06b337dSHong Zhang       }
2050d06b337dSHong Zhang     }
2051d06b337dSHong Zhang 
20522798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
205344706875SHong Zhang       /*
205444706875SHong Zhang        backward sweep:
205544706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
205644706875SHong Zhang        for i=m-1,...,0:
205744706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
205844706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
205944706875SHong Zhang       */
206095d750ceSBarry Smith       /* if there was a forward sweep done above then I thing the next two for loops are not needed */
2061290bbb0aSBarry Smith       ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr);
2062d06b337dSHong Zhang 
2063d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2064d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2065d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2066d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2067efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2068e6b907acSBarry Smith         while (nz--) t[*vj++] -= x[i]*(*v++);
2069d06b337dSHong Zhang       }
2070d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2071d06b337dSHong Zhang         d  = *(aa + ai[i]);
2072d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2073d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2074d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2075e6b907acSBarry Smith         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2076d06b337dSHong Zhang         sum = t[i];
2077d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
2078d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2079d06b337dSHong Zhang       }
2080d06b337dSHong Zhang     }
2081d06b337dSHong Zhang   }
2082d06b337dSHong Zhang 
20831ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2084d06b337dSHong Zhang   if (bb != xx) {
20851ebc52fbSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2086d06b337dSHong Zhang   }
2087d06b337dSHong Zhang   PetscFunctionReturn(0);
2088d06b337dSHong Zhang }
2089d06b337dSHong Zhang 
2090c75a6043SHong Zhang #undef __FUNCT__
2091c75a6043SHong Zhang #define __FUNCT__ "MatCreateSeqSBAIJWithArrays"
2092c75a6043SHong Zhang /*@
2093c75a6043SHong Zhang      MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2094c75a6043SHong Zhang               (upper triangular entries in CSR format) provided by the user.
2095c75a6043SHong Zhang 
2096c75a6043SHong Zhang      Collective on MPI_Comm
2097c75a6043SHong Zhang 
2098c75a6043SHong Zhang    Input Parameters:
2099c75a6043SHong Zhang +  comm - must be an MPI communicator of size 1
2100c75a6043SHong Zhang .  bs - size of block
2101c75a6043SHong Zhang .  m - number of rows
2102c75a6043SHong Zhang .  n - number of columns
2103c75a6043SHong Zhang .  i - row indices
2104c75a6043SHong Zhang .  j - column indices
2105c75a6043SHong Zhang -  a - matrix values
2106c75a6043SHong Zhang 
2107c75a6043SHong Zhang    Output Parameter:
2108c75a6043SHong Zhang .  mat - the matrix
2109c75a6043SHong Zhang 
2110c75a6043SHong Zhang    Level: intermediate
2111c75a6043SHong Zhang 
2112c75a6043SHong Zhang    Notes:
2113c75a6043SHong Zhang        The i, j, and a arrays are not copied by this routine, the user must free these arrays
2114c75a6043SHong Zhang     once the matrix is destroyed
2115c75a6043SHong Zhang 
2116c75a6043SHong Zhang        You cannot set new nonzero locations into this matrix, that will generate an error.
2117c75a6043SHong Zhang 
2118c75a6043SHong Zhang        The i and j indices are 0 based
2119c75a6043SHong Zhang 
2120c75a6043SHong Zhang .seealso: MatCreate(), MatCreateMPISBAIJ(), MatCreateSeqSBAIJ()
2121c75a6043SHong Zhang 
2122c75a6043SHong Zhang @*/
2123c75a6043SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat)
2124c75a6043SHong Zhang {
2125c75a6043SHong Zhang   PetscErrorCode ierr;
2126c75a6043SHong Zhang   PetscInt       ii;
2127c75a6043SHong Zhang   Mat_SeqSBAIJ   *sbaij;
2128c75a6043SHong Zhang 
2129c75a6043SHong Zhang   PetscFunctionBegin;
2130c75a6043SHong Zhang   if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
2131c75a6043SHong Zhang   if (i[0]) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2132c75a6043SHong Zhang 
2133c75a6043SHong Zhang   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
2134c75a6043SHong Zhang   ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr);
2135c75a6043SHong Zhang   ierr = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr);
2136c75a6043SHong Zhang   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr);
2137c75a6043SHong Zhang   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2138c75a6043SHong Zhang   ierr = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr);
2139c75a6043SHong Zhang 
2140c75a6043SHong Zhang   sbaij->i = i;
2141c75a6043SHong Zhang   sbaij->j = j;
2142c75a6043SHong Zhang   sbaij->a = a;
2143c75a6043SHong Zhang   sbaij->singlemalloc = PETSC_FALSE;
2144c75a6043SHong Zhang   sbaij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2145e6b907acSBarry Smith   sbaij->free_a       = PETSC_FALSE;
2146e6b907acSBarry Smith   sbaij->free_ij      = PETSC_FALSE;
2147c75a6043SHong Zhang 
2148c75a6043SHong Zhang   for (ii=0; ii<m; ii++) {
2149c75a6043SHong Zhang     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2150c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG)
2151c75a6043SHong Zhang     if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
2152c75a6043SHong Zhang #endif
2153c75a6043SHong Zhang   }
2154c75a6043SHong Zhang #if defined(PETSC_USE_DEBUG)
2155c75a6043SHong Zhang   for (ii=0; ii<sbaij->i[m]; ii++) {
2156c75a6043SHong Zhang     if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2157c75a6043SHong Zhang     if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
2158c75a6043SHong Zhang   }
2159c75a6043SHong Zhang #endif
2160c75a6043SHong Zhang 
2161c75a6043SHong Zhang   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2162c75a6043SHong Zhang   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2163c75a6043SHong Zhang   PetscFunctionReturn(0);
2164c75a6043SHong Zhang }
2165d06b337dSHong Zhang 
2166d06b337dSHong Zhang 
2167d06b337dSHong Zhang 
216849b5e25fSSatish Balay 
216949b5e25fSSatish Balay 
2170